Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Friends | List of all members
Nektar::LibUtilities::NekLinSysIterCGLoc Class Reference

#include <NekLinSysIterCGLoc.h>

Inheritance diagram for Nektar::LibUtilities::NekLinSysIterCGLoc:
[legend]

Public Member Functions

 NekLinSysIterCGLoc (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 Constructor for full direct matrix solve. More...
 
 ~NekLinSysIterCGLoc () override=default
 
- Public Member Functions inherited from Nektar::LibUtilities::NekLinSysIter
 NekLinSysIter (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
 ~NekLinSysIter () override=default
 
void SetUniversalUniqueMap (const Array< OneD, const int > &map)
 
void SetNekLinSysTolerance (const NekDouble in)
 
void SetNekLinSysMaxIterations (const unsigned int in)
 
int GetNekLinSysTolerance ()
 
bool IsLocal ()
 
- Public Member Functions inherited from Nektar::LibUtilities::NekSys
 NekSys (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
virtual ~NekSys ()=default
 
void InitObject ()
 
int SolveSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir=0)
 
const NekSysOperatorsGetSysOperators ()
 
void SetSysOperators (const NekSysOperators &in)
 
void SetFlagWarnings (bool in)
 
void SetRhsMagnitude (const NekDouble mag)
 

Static Public Member Functions

static NekLinSysIterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
- Static Public Member Functions inherited from Nektar::LibUtilities::NekLinSysIter
static NekLinSysIterSharedPtr CreateInstance (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 
- Static Public Member Functions inherited from Nektar::LibUtilities::NekSys
static NekSysSharedPtr CreateInstance (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

void v_InitObject () override
 
int v_SolveSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir) override
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekLinSysIter
void v_InitObject () override
 
void SetUniversalUniqueMap ()
 
void Set_Rhs_Magnitude (const Array< OneD, NekDouble > &pIn)
 
void ConvergenceCheck (const Array< OneD, const NekDouble > &Residual)
 
- Protected Member Functions inherited from Nektar::LibUtilities::NekSys
virtual void v_InitObject ()
 
virtual void v_SetSysOperators (const NekSysOperators &in)
 
virtual int v_SolveSystem (const int nGlobal, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int nDir)
 

Private Member Functions

void DoConjugateGradient (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Actual iterative solve. More...
 

Friends

class MemoryManager< NekLinSysIterCGLoc >
 Support creation through MemoryManager. More...
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::LibUtilities::NekLinSysIter
Array< OneD, int > m_map
 Global to universal unique map. More...
 
NekDouble m_NekLinSysTolerance
 
int m_NekLinSysMaxIterations
 
int m_totalIterations = 0
 
bool m_isLocal
 
- Protected Attributes inherited from Nektar::LibUtilities::NekSys
LibUtilities::CommSharedPtr m_rowComm
 
bool m_converged
 
bool m_root
 
bool m_verbose
 
bool m_FlagWarnings
 
int m_SysDimen
 
NekSysOperators m_operator
 
NekDouble m_rhs_magnitude = NekConstants::kNekUnsetDouble
 

Detailed Description

Solves a linear system using iterative methods.

Definition at line 48 of file NekLinSysIterCGLoc.h.

Constructor & Destructor Documentation

◆ NekLinSysIterCGLoc()

Nektar::LibUtilities::NekLinSysIterCGLoc::NekLinSysIterCGLoc ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vRowComm,
const int  nDimen,
const NekSysKey pKey 
)

Constructor for full direct matrix solve.

Definition at line 52 of file NekLinSysIterCGLoc.cpp.

56 : NekLinSysIter(pSession, vRowComm, nDimen, pKey)
57{
58 m_isLocal = true;
59}
NekLinSysIter(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)

References Nektar::LibUtilities::NekLinSysIter::m_isLocal.

◆ ~NekLinSysIterCGLoc()

Nektar::LibUtilities::NekLinSysIterCGLoc::~NekLinSysIterCGLoc ( )
overridedefault

Member Function Documentation

◆ create()

static NekLinSysIterSharedPtr Nektar::LibUtilities::NekLinSysIterCGLoc::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vRowComm,
const int  nDimen,
const NekSysKey pKey 
)
inlinestatic

Definition at line 54 of file NekLinSysIterCGLoc.h.

58 {
61 pSession, vRowComm, nDimen, pKey);
62 p->InitObject();
63 return p;
64 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< NekLinSysIterCGLoc > NekLinSysIterCGLocSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ DoConjugateGradient()

void Nektar::LibUtilities::NekLinSysIterCGLoc::DoConjugateGradient ( const int  nLocal,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
private

Actual iterative solve.

  Solve a global linear system using the conjugate gradient method.   We solve only for the non-Dirichlet modes. The operator is evaluated   using an auxiliary function m_operator.DoNekSysLhsEval defined by the   specific solver. Distributed math routines are used to support   parallel execution of the solver.     The implemented algorithm uses a reduced-communication reordering of   the standard PCG method (Demmel, Heath and Vorst, 1993)    

Parameters
pInputInput residual of all DOFs.  
pOutputSolution vector of all DOFs.  

Definition at line 91 of file NekLinSysIterCGLoc.cpp.

94{
95 // Allocate array storage
96 Array<OneD, NekDouble> w_A(nLocal, 0.0);
97 Array<OneD, NekDouble> s_A(nLocal, 0.0);
98 Array<OneD, NekDouble> p_A(nLocal, 0.0);
99 Array<OneD, NekDouble> r_A(nLocal, 0.0);
100 Array<OneD, NekDouble> q_A(nLocal, 0.0);
101 Array<OneD, NekDouble> wk(nLocal, 0.0);
102
103 NekDouble alpha;
105 NekDouble rho;
106 NekDouble rho_new;
107 NekDouble mu;
108 NekDouble eps;
109 Array<OneD, NekDouble> vExchange(3, 0.0);
110
111 // Copy initial residual from input
112 Vmath::Vcopy(nLocal, pInput, 1, r_A, 1);
113
114 // zero homogeneous out array ready for solution updates
115 // Should not be earlier in case input vector is same as
116 // output and above copy has been peformed
117 Vmath::Zero(nLocal, pOutput, 1);
118
119 // evaluate initial residual error for exit check
120 m_operator.DoAssembleLoc(r_A, wk, true);
121 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
122
123 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
124
125 eps = vExchange[2];
126
128 {
129 Set_Rhs_Magnitude(pInput);
130 }
131
133
134 // If input residual is less than tolerance skip solve.
136 {
137 if (m_verbose && m_root)
138 {
139 cout << "CG iterations made = " << m_totalIterations
140 << " using tolerance of " << m_NekLinSysTolerance
141 << " (error = " << sqrt(eps / m_rhs_magnitude)
142 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
143 }
144 return;
145 }
146
147 m_operator.DoNekSysPrecon(r_A, w_A, true);
148 m_operator.DoNekSysLhsEval(w_A, s_A);
149
150 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
151 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
152
153 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
154
155 rho = vExchange[0];
156 mu = vExchange[1];
157 beta = 0.0;
158 alpha = rho / mu;
160
161 // Continue until convergence
162 while (true)
163 {
165 {
166 if (m_root)
167 {
168 cout << "CG iterations made = " << m_totalIterations
169 << " using tolerance of " << m_NekLinSysTolerance
170 << " (error = " << sqrt(eps / m_rhs_magnitude)
171 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
172 }
174 "Exceeded maximum number of iterations");
175 }
176
177 // Compute new search direction p_k, q_k
178 Vmath::Svtvp(nLocal, beta, p_A, 1, w_A, 1, p_A, 1);
179 Vmath::Svtvp(nLocal, beta, q_A, 1, s_A, 1, q_A, 1);
180
181 // Update solution x_{k+1}
182 Vmath::Svtvp(nLocal, alpha, p_A, 1, pOutput, 1, pOutput, 1);
183
184 // Update residual vector r_{k+1}
185 Vmath::Svtvp(nLocal, -alpha, q_A, 1, r_A, 1, r_A, 1);
186
187 // Apply preconditioner
188 m_operator.DoNekSysPrecon(r_A, w_A, true);
189
190 // Perform the method-specific matrix-vector multiply operation.
191 m_operator.DoNekSysLhsEval(w_A, s_A);
192
193 // <r_{k+1}, w_{k+1}>
194 vExchange[0] = Vmath::Dot(nLocal, r_A, w_A);
195
196 // <s_{k+1}, w_{k+1}>
197 vExchange[1] = Vmath::Dot(nLocal, s_A, w_A);
198
199 // <r_{k+1}, r_{k+1}>
200 m_operator.DoAssembleLoc(r_A, wk, true);
201 vExchange[2] = Vmath::Dot(nLocal, wk, r_A);
202
203 // Perform inner-product exchanges
204 m_rowComm->AllReduce(vExchange, Nektar::LibUtilities::ReduceSum);
205
206 rho_new = vExchange[0];
207 mu = vExchange[1];
208 eps = vExchange[2];
209
211
212 // Test if norm is within tolerance
214 {
215 if (m_verbose && m_root)
216 {
217 cout << "CG iterations made = " << m_totalIterations
218 << " using tolerance of " << m_NekLinSysTolerance
219 << " (error = " << sqrt(eps / m_rhs_magnitude)
220 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
221 }
222 break;
223 }
224
225 // Compute search direction and solution coefficients
226 beta = rho_new / rho;
227 alpha = rho_new / (mu - rho_new * beta / alpha);
228 rho = rho_new;
229 }
230}
#define ROOTONLY_NEKERROR(type, msg)
Definition: ErrorUtil.hpp:205
void Set_Rhs_Magnitude(const Array< OneD, NekDouble > &pIn)
LibUtilities::CommSharedPtr m_rowComm
Definition: NekSys.h:291
NekSysOperators m_operator
Definition: NekSys.h:298
void DoAssembleLoc(InArrayType &xn, OutArrayType &xn1, const bool &flag=false) const
Definition: NekSys.h:162
void DoNekSysPrecon(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:149
void DoNekSysLhsEval(InArrayType &inarray, OutArrayType &outarray, const bool &flag=false) const
Definition: NekSys.h:142
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
static const NekDouble kNekUnsetDouble
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::beta, Nektar::LibUtilities::NekSysOperators::DoAssembleLoc(), Nektar::LibUtilities::NekSysOperators::DoNekSysLhsEval(), Nektar::LibUtilities::NekSysOperators::DoNekSysPrecon(), Vmath::Dot(), Nektar::ErrorUtil::efatal, Nektar::NekConstants::kNekUnsetDouble, Nektar::LibUtilities::NekLinSysIter::m_NekLinSysMaxIterations, Nektar::LibUtilities::NekLinSysIter::m_NekLinSysTolerance, Nektar::LibUtilities::NekSys::m_operator, Nektar::LibUtilities::NekSys::m_rhs_magnitude, Nektar::LibUtilities::NekSys::m_root, Nektar::LibUtilities::NekSys::m_rowComm, Nektar::LibUtilities::NekLinSysIter::m_totalIterations, Nektar::LibUtilities::NekSys::m_verbose, Nektar::LibUtilities::ReduceSum, ROOTONLY_NEKERROR, Nektar::LibUtilities::NekLinSysIter::Set_Rhs_Magnitude(), tinysimd::sqrt(), Vmath::Svtvp(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_SolveSystem().

◆ v_InitObject()

void Nektar::LibUtilities::NekLinSysIterCGLoc::v_InitObject ( )
overrideprotectedvirtual

◆ v_SolveSystem()

int Nektar::LibUtilities::NekLinSysIterCGLoc::v_SolveSystem ( const int  nGlobal,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const int  nDir 
)
overrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::NekSys.

Definition at line 69 of file NekLinSysIterCGLoc.cpp.

72{
73 DoConjugateGradient(nLocal, pInput, pOutput);
74
75 return m_totalIterations;
76}
void DoConjugateGradient(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Actual iterative solve.

References DoConjugateGradient(), and Nektar::LibUtilities::NekLinSysIter::m_totalIterations.

Friends And Related Function Documentation

◆ MemoryManager< NekLinSysIterCGLoc >

friend class MemoryManager< NekLinSysIterCGLoc >
friend

Support creation through MemoryManager.

Definition at line 46 of file NekLinSysIterCGLoc.h.

Member Data Documentation

◆ className

string Nektar::LibUtilities::NekLinSysIterCGLoc::className
static
Initial value:
=
"ConjugateGradientLoc", NekLinSysIterCGLoc::create,
"NekLinSysIterCG solver in Local Space.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static NekLinSysIterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vRowComm, const int nDimen, const NekSysKey &pKey)
NekLinSysIterFactory & GetNekLinSysIterFactory()

Definition at line 66 of file NekLinSysIterCGLoc.h.