Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
Nektar::MultiRegions::GlobalLinSysIterative Class Referenceabstract

A global linear system. More...

#include <GlobalLinSysIterative.h>

Inheritance diagram for Nektar::MultiRegions::GlobalLinSysIterative:
[legend]

Public Member Functions

 GlobalLinSysIterative (const GlobalLinSysKey &pKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve. More...
 
 ~GlobalLinSysIterative () override
 
void DoMatrixMultiply (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
- Public Member Functions inherited from Nektar::MultiRegions::GlobalLinSys
 GlobalLinSys (const GlobalLinSysKey &pKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 Constructor for full direct matrix solve. More...
 
virtual ~GlobalLinSys ()
 
const GlobalLinSysKeyGetKey (void) const
 Returns the key associated with the system. More...
 
const std::weak_ptr< ExpList > & GetLocMat (void) const
 
void InitObject ()
 
void Initialise (const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 
void Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve the linear system for given input and output vectors using a specified local to global map. More...
 
std::shared_ptr< GlobalLinSysGetSharedThisPtr ()
 Returns a shared pointer to the current object. More...
 
int GetNumBlocks ()
 
DNekScalMatSharedPtr GetBlock (unsigned int n)
 
void DropBlock (unsigned int n)
 
DNekScalBlkMatSharedPtr GetStaticCondBlock (unsigned int n)
 
void DropStaticCondBlock (unsigned int n)
 
void SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir=0)
 Solve the linear system for given input and output vectors. More...
 

Protected Member Functions

void DoProjection (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir, const NekDouble tol, const bool isAconjugate)
 projection technique More...
 
void Set_Rhs_Magnitude (const NekVector< NekDouble > &pIn)
 
virtual void v_UniqueMap ()=0
 
virtual void v_DoMatrixMultiply (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)=0
 
- Protected Member Functions inherited from Nektar::MultiRegions::GlobalLinSys
virtual void v_Solve (const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)=0
 Solve a linear system based on mapping. More...
 
virtual void v_SolveLinearSystem (const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir)=0
 Solve a basic matrix system. More...
 
virtual void v_InitObject ()
 
virtual void v_Initialise (const std::shared_ptr< AssemblyMap > &pLocToGloMap)
 
virtual int v_GetNumBlocks ()
 Get the number of blocks in this system. More...
 
virtual DNekScalMatSharedPtr v_GetBlock (unsigned int n)
 Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual void v_DropBlock (unsigned int n)
 Releases the local block matrix from NekManager of n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock (unsigned int n)
 Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided by the m_linSysKey. More...
 
virtual void v_DropStaticCondBlock (unsigned int n)
 Releases the static condensation block matrices from NekManager of n-th expansion using the matrix key provided by the m_linSysKey. More...
 
PreconditionerSharedPtr CreatePrecon (AssemblyMapSharedPtr asmMap)
 Create a preconditioner object from the parameters defined in the supplied assembly map. More...
 

Protected Attributes

Array< OneD, int > m_map
 Global to universal unique map. More...
 
NekDouble m_tolerance
 Tolerance of iterative solver. More...
 
NekDouble m_rhs_magnitude
 dot product of rhs to normalise stopping criterion More...
 
NekDouble m_rhs_mag_sm
 cnt to how many times rhs_magnitude is called More...
 
PreconditionerSharedPtr m_precon
 
std::string m_precontype
 
int m_totalIterations
 
bool m_useProjection
 Whether to apply projection technique. More...
 
bool m_root
 Root if parallel. More...
 
std::string m_linSysIterSolver
 Iterative solver: Conjugate Gradient, GMRES. More...
 
std::vector< Array< OneD, NekDouble > > m_prevLinSol
 Storage for solutions to previous linear problems. More...
 
std::vector< Array< OneD, NekDouble > > m_prevBasis
 
DNekMatSharedPtr m_coeffMatrix
 
Array< OneD, NekDoublem_coeffMatrixFactor
 
Array< OneD, int > m_ipivot
 
int m_numSuccessiveRHS
 
bool m_isAconjugate
 
std::string m_matrixType
 
bool m_isNonSymmetricLinSys
 
int m_numPrevSols
 
bool m_isAbsoluteTolerance
 
LibUtilities::NekSysOperators m_NekSysOp
 
LibUtilities::NekLinSysIterSharedPtr m_linsol
 
- Protected Attributes inherited from Nektar::MultiRegions::GlobalLinSys
const GlobalLinSysKey m_linSysKey
 Key associated with this linear system. More...
 
const std::weak_ptr< ExpListm_expList
 Local Matrix System. More...
 
const std::map< int, RobinBCInfoSharedPtrm_robinBCInfo
 Robin boundary info. More...
 
bool m_verbose
 

Static Protected Attributes

static std::string IteratSolverlookupIds []
 
static std::string IteratSolverdef
 

Private Member Functions

void UpdateKnownSolutions (const int pGlobalBndDofs, const Array< OneD, const NekDouble > &pSolution, const int pNumDirBndDofs, const bool isAconjugate)
 
int ResetKnownSolutionsToLatestOne ()
 
void DoPreconditionerFlag (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
 
void DoMatrixMultiplyFlag (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &controlFlag)
 
void DoAssembleLocFlag (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &ZeroDir)
 

Detailed Description

A global linear system.

Solves a linear system using iterative methods.

Definition at line 48 of file GlobalLinSysIterative.h.

Constructor & Destructor Documentation

◆ GlobalLinSysIterative()

Nektar::MultiRegions::GlobalLinSysIterative::GlobalLinSysIterative ( const GlobalLinSysKey pKey,
const std::weak_ptr< ExpList > &  pExpList,
const std::shared_ptr< AssemblyMap > &  pLocToGloMap 
)

Constructor for full direct matrix solve.

Constructor for full iterative matrix solve.

Definition at line 65 of file GlobalLinSysIterative.cpp.

68 : GlobalLinSys(pKey, pExpList, pLocToGloMap),
72{
73 m_tolerance = pLocToGloMap->GetIterativeTolerance();
74 m_isAbsoluteTolerance = pLocToGloMap->IsAbsoluteTolerance();
75 m_linSysIterSolver = pLocToGloMap->GetLinSysIterSolver();
76
78 m_expList.lock()->GetComm()->GetRowComm();
79 m_root = (vComm->GetRank()) ? false : true;
80
81 m_numSuccessiveRHS = pLocToGloMap->GetSuccessiveRHS();
85
86 // Check for advection matrix and switch to GMRES, if not already used
87 m_matrixType = StdRegions::MatrixTypeMap[pKey.GetMatrixType()];
89 m_matrixType.find("AdvectionDiffusionReaction") != string::npos;
90
92 !m_linSysIterSolver.compare("ConjugateGradient"))
93 {
94 m_linSysIterSolver = "GMRES";
96 false,
97 "Detected ConjugateGradient solver and a "
98 "Advection-Diffusion-Reaction matrix. "
99 "Switchted to a GMRES solver for this non-symmetric matrix type. "
100 "Change LinSysIterSolver to GMRES in the session file to suppress "
101 "this warning.");
102 }
103
104 if (m_isAconjugate && m_linSysIterSolver.compare("GMRES") == 0 &&
105 m_linSysIterSolver.compare("GMRESLoc") == 0)
106 {
107 WARNINGL0(false, "To use A-conjugate projection, the matrix "
108 "should be symmetric positive definite.");
109 }
110
112 this);
114 this);
116 this);
117}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:120
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:106
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:113
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:122
GlobalLinSys(const GlobalLinSysKey &pKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
void DoAssembleLocFlag(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &ZeroDir)
void DoPreconditionerFlag(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
void DoMatrixMultiplyFlag(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &controlFlag)
bool m_useProjection
Whether to apply projection technique.
NekDouble m_tolerance
Tolerance of iterative solver.
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
NekDouble m_rhs_mag_sm
cnt to how many times rhs_magnitude is called
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
static PreconditionerSharedPtr NullPreconditionerSharedPtr
static const NekDouble kNekUnsetDouble
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:139
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298

References tinysimd::abs(), Nektar::LibUtilities::NekSysOperators::DefineAssembleLoc(), Nektar::LibUtilities::NekSysOperators::DefineNekSysLhsEval(), Nektar::LibUtilities::NekSysOperators::DefineNekSysPrecon(), DoAssembleLocFlag(), DoMatrixMultiplyFlag(), DoPreconditionerFlag(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::MultiRegions::GlobalLinSys::m_expList, m_isAbsoluteTolerance, m_isAconjugate, m_isNonSymmetricLinSys, m_linSysIterSolver, m_matrixType, m_NekSysOp, m_numSuccessiveRHS, m_root, m_tolerance, m_useProjection, Nektar::StdRegions::MatrixTypeMap, and WARNINGL0.

◆ ~GlobalLinSysIterative()

Nektar::MultiRegions::GlobalLinSysIterative::~GlobalLinSysIterative ( )
override

Definition at line 119 of file GlobalLinSysIterative.cpp.

120{
121}

Member Function Documentation

◆ DoAssembleLocFlag()

void Nektar::MultiRegions::GlobalLinSysIterative::DoAssembleLocFlag ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  ZeroDir 
)
inlineprivate

Definition at line 147 of file GlobalLinSysIterative.h.

149 {
150 m_precon->DoAssembleLoc(pInput, pOutput, ZeroDir);
151 }

References m_precon.

Referenced by GlobalLinSysIterative().

◆ DoMatrixMultiply()

void Nektar::MultiRegions::GlobalLinSysIterative::DoMatrixMultiply ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
inline

Definition at line 58 of file GlobalLinSysIterative.h.

60 {
61 v_DoMatrixMultiply(pInput, pOutput);
62 }
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)=0

References v_DoMatrixMultiply().

◆ DoMatrixMultiplyFlag()

void Nektar::MultiRegions::GlobalLinSysIterative::DoMatrixMultiplyFlag ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  controlFlag 
)
inlineprivate

Definition at line 140 of file GlobalLinSysIterative.h.

143 {
144 v_DoMatrixMultiply(pInput, pOutput);
145 }

References v_DoMatrixMultiply().

Referenced by DoProjection(), GlobalLinSysIterative(), and UpdateKnownSolutions().

◆ DoPreconditionerFlag()

void Nektar::MultiRegions::GlobalLinSysIterative::DoPreconditionerFlag ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  isLocal = false 
)
inlineprivate

Definition at line 133 of file GlobalLinSysIterative.h.

136 {
137 m_precon->DoPreconditioner(pInput, pOutput, isLocal);
138 }

References m_precon.

Referenced by GlobalLinSysIterative().

◆ DoProjection()

void Nektar::MultiRegions::GlobalLinSysIterative::DoProjection ( const int  nGlobal,
const Array< OneD, const NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const int  nDir,
const NekDouble  tol,
const bool  isAconjugate 
)
protected

projection technique

This method implements projection techniques in order to speed up successive linear solves with right-hand sides arising from time-dependent discretisations. (P.F.Fischer, Comput. Methods Appl. Mech. Engrg. 163, 1998)

Definition at line 129 of file GlobalLinSysIterative.cpp.

133{
134 int numIterations = 0;
135 if (0 == m_numPrevSols)
136 {
137 // no previous solutions found
138 numIterations =
139 m_linsol->SolveSystem(nGlobal, pInput, pOutput, nDir, tol);
140 }
141 else
142 {
143 // Get the communicator for performing data exchanges
145 m_expList.lock()->GetComm()->GetRowComm();
146
147 // Get vector sizes
148 int nNonDir = nGlobal - nDir;
149
150 // check the input vector (rhs) is not zero
151 Array<OneD, NekDouble> tmp;
152
153 NekDouble rhsNorm =
154 Vmath::Dot2(nNonDir, pInput + nDir, pInput + nDir, m_map + nDir);
155
156 vComm->AllReduce(rhsNorm, Nektar::LibUtilities::ReduceSum);
157
158 if (rhsNorm < tol * tol * m_rhs_magnitude)
159 {
160 Vmath::Zero(nNonDir, tmp = pOutput + nDir, 1);
161 if (m_verbose && m_root)
162 {
163 cout << "No iterations made"
164 << " using tolerance of " << tol
165 << " (error = " << sqrt(rhsNorm / m_rhs_magnitude)
166 << ", rhs_mag = " << sqrt(m_rhs_magnitude) << ")" << endl;
167 }
168 return;
169 }
170
171 // Create NekVector wrappers for linear algebra operations
172 NekVector<NekDouble> b(nNonDir, pInput + nDir, eWrapper);
173 NekVector<NekDouble> x(nNonDir, tmp = pOutput + nDir, eWrapper);
174 // Allocate array storage
175 Array<OneD, NekDouble> px_s(nGlobal, 0.0);
176 Array<OneD, NekDouble> pb_s(nGlobal, 0.0);
177 Array<OneD, NekDouble> tmpAx_s(nGlobal, 0.0);
178 Array<OneD, NekDouble> tmpx_s(nGlobal, 0.0);
179
180 NekVector<NekDouble> pb(nNonDir, tmp = pb_s + nDir, eWrapper);
181 NekVector<NekDouble> px(nNonDir, tmp = px_s + nDir, eWrapper);
182 NekVector<NekDouble> tmpAx(nNonDir, tmp = tmpAx_s + nDir, eWrapper);
183 NekVector<NekDouble> tmpx(nNonDir, tmp = tmpx_s + nDir, eWrapper);
184
185 // notation follows the paper cited:
186 // \alpha_i = \tilda{x_i}^T b^n
187 // projected x, px = \sum \alpha_i \tilda{x_i}
188
189 Array<OneD, NekDouble> alpha(m_prevBasis.size(), 0.0);
190 Array<OneD, NekDouble> alphaback(m_prevBasis.size(), 0.0);
191 for (int i = 0; i < m_prevBasis.size(); i++)
192 {
193 alpha[i] = Vmath::Dot2(nNonDir, m_prevBasis[i], pInput + nDir,
194 m_map + nDir);
195 }
196 vComm->AllReduce(alpha, Nektar::LibUtilities::ReduceSum);
197 int n = m_prevBasis.size(), info = -1;
198 Vmath::Vcopy(m_prevBasis.size(), alpha, 1, alphaback, 1);
199 Lapack::Dsptrs('U', n, 1, m_coeffMatrixFactor.get(), m_ipivot.get(),
200 alpha.get(), n, info);
201 if (info != 0)
202 {
203 // Dsptrs fails, only keep the latest solution
204 int latest = ResetKnownSolutionsToLatestOne();
205 alpha[0] = alphaback[latest];
206 }
207 for (int i = 0; i < m_prevBasis.size(); ++i)
208 {
209 NekVector<NekDouble> xi(nNonDir, m_prevLinSol[i], eWrapper);
210 px += alpha[i] * xi;
211 }
212
213 // pb = b^n - A px
214 Vmath::Vcopy(nNonDir, pInput.get() + nDir, 1, pb_s.get() + nDir, 1);
215
216 DoMatrixMultiplyFlag(px_s, tmpAx_s, false);
217
218 pb -= tmpAx;
219
220 if (m_verbose)
221 {
222 if (m_root)
223 {
224 cout << "SuccessiveRHS: " << m_prevBasis.size()
225 << "-bases projection reduces L2-norm of RHS from "
226 << std::sqrt(rhsNorm) << " to ";
227 }
228 NekDouble tmprhsNorm =
229 Vmath::Dot2(nNonDir, pb_s + nDir, pb_s + nDir, m_map + nDir);
230 vComm->AllReduce(tmprhsNorm, Nektar::LibUtilities::ReduceSum);
231 if (m_root)
232 {
233 cout << std::sqrt(tmprhsNorm) << endl;
234 }
235 }
236
237 // solve the system with projected rhs
238 numIterations = m_linsol->SolveSystem(nGlobal, pb_s, tmpx_s, nDir, tol);
239
240 // remainder solution + projection of previous solutions
241 x = tmpx + px;
242 }
243 // save the auxiliary solution to prev. known solutions
244 if (numIterations)
245 {
246 UpdateKnownSolutions(nGlobal, pOutput, nDir, isAconjugate);
247 }
248}
std::vector< Array< OneD, NekDouble > > m_prevBasis
LibUtilities::NekLinSysIterSharedPtr m_linsol
void UpdateKnownSolutions(const int pGlobalBndDofs, const Array< OneD, const NekDouble > &pSolution, const int pNumDirBndDofs, const bool isAconjugate)
std::vector< Array< OneD, NekDouble > > m_prevLinSol
Storage for solutions to previous linear problems.
Array< OneD, int > m_map
Global to universal unique map.
static void Dsptrs(const char &uplo, const int &n, const int &nrhs, const double *ap, const int *ipiv, double *b, const int &ldb, int &info)
Solve a real packed-symmetric matrix problem using Bunch-Kaufman pivoting.
Definition: Lapack.hpp:154
double NekDouble
T Dot2(int n, const T *w, const T *x, const int *y)
dot product
Definition: Vmath.hpp:790
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 DoMatrixMultiplyFlag(), Vmath::Dot2(), Lapack::Dsptrs(), Nektar::eWrapper, m_coeffMatrixFactor, Nektar::MultiRegions::GlobalLinSys::m_expList, m_ipivot, m_linsol, m_map, m_numPrevSols, m_prevBasis, m_prevLinSol, m_rhs_magnitude, m_root, Nektar::MultiRegions::GlobalLinSys::m_verbose, Nektar::LibUtilities::ReduceSum, ResetKnownSolutionsToLatestOne(), tinysimd::sqrt(), UpdateKnownSolutions(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by Nektar::MultiRegions::GlobalLinSysIterativeFull::v_SolveLinearSystem(), and Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_SolveLinearSystem().

◆ ResetKnownSolutionsToLatestOne()

int Nektar::MultiRegions::GlobalLinSysIterative::ResetKnownSolutionsToLatestOne ( )
private

Definition at line 250 of file GlobalLinSysIterative.cpp.

251{
252 if (m_numPrevSols == 0)
253 {
254 return -1;
255 }
257 Array<OneD, NekDouble> b = m_prevBasis[latest];
258 Array<OneD, NekDouble> x = m_prevLinSol[latest];
259 m_prevBasis.clear();
260 m_prevLinSol.clear();
261 m_prevBasis.push_back(b);
262 m_prevLinSol.push_back(x);
263 m_numPrevSols = 1;
264 return latest;
265}

References m_numPrevSols, m_numSuccessiveRHS, m_prevBasis, and m_prevLinSol.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ Set_Rhs_Magnitude()

void Nektar::MultiRegions::GlobalLinSysIterative::Set_Rhs_Magnitude ( const NekVector< NekDouble > &  pIn)
protected

◆ UpdateKnownSolutions()

void Nektar::MultiRegions::GlobalLinSysIterative::UpdateKnownSolutions ( const int  nGlobal,
const Array< OneD, const NekDouble > &  newX,
const int  nDir,
const bool  isAconjugate 
)
private

Updates the storage of previously known solutions. Performs normalisation of input vector wrt A-norm.

Definition at line 271 of file GlobalLinSysIterative.cpp.

274{
275 // Get vector sizes
276 int nNonDir = nGlobal - nDir;
277 int insertLocation = m_numPrevSols % m_numSuccessiveRHS;
278 int fullbuffer = (m_prevBasis.size() == m_numSuccessiveRHS);
279
280 // Get the communicator for performing data exchanges
282 m_expList.lock()->GetComm()->GetRowComm();
283
284 Array<OneD, NekDouble> tmpAx_s(nGlobal, 0.0);
285 Array<OneD, NekDouble> y_s(m_prevBasis.size() - fullbuffer, 0.0);
286 Array<OneD, NekDouble> invMy_s(y_s.size(), 0.0);
287 Array<OneD, int> ipivot(m_numSuccessiveRHS);
288 Array<OneD, NekDouble> tmp, newBasis;
289
290 DoMatrixMultiplyFlag(newX, tmpAx_s, false);
291
292 if (isAconjugate)
293 {
294 newBasis = newX + nDir;
295 }
296 else
297 {
298 newBasis = tmpAx_s + nDir;
299 }
300
301 // Check the solution is non-zero
302 NekDouble solNorm =
303 Vmath::Dot2(nNonDir, newBasis, tmpAx_s + nDir, m_map + nDir);
304 vComm->AllReduce(solNorm, Nektar::LibUtilities::ReduceSum);
305
307 {
308 return;
309 }
310
311 // normalisation of A x
312 Vmath::Smul(nNonDir, 1.0 / sqrt(solNorm), tmpAx_s + nDir, 1,
313 tmp = tmpAx_s + nDir, 1);
314
315 for (int i = 0; i < m_prevBasis.size(); ++i)
316 {
317 if (i == insertLocation)
318 {
319 continue;
320 }
321 int skip = i > insertLocation;
322 y_s[i - skip] =
323 Vmath::Dot2(nNonDir, m_prevBasis[i], tmpAx_s + nDir, m_map + nDir);
324 }
325 vComm->AllReduce(y_s, Nektar::LibUtilities::ReduceSum);
326
327 // check if linearly dependent
328 DNekMatSharedPtr tilCoeffMatrix;
329 if (fullbuffer && m_numSuccessiveRHS > 1)
330 {
333 for (int i = 0; i < m_numSuccessiveRHS; ++i)
334 {
335 if (i == insertLocation)
336 {
337 continue;
338 }
339 int iskip = i > insertLocation;
340 for (int j = i; j < m_numSuccessiveRHS; ++j)
341 {
342 if (j == insertLocation)
343 {
344 continue;
345 }
346 int jskip = j > insertLocation;
347 tilCoeffMatrix->SetValue(i - iskip, j - jskip,
348 m_coeffMatrix->GetValue(i, j));
349 }
350 }
351 }
352 else if (!fullbuffer && m_prevBasis.size())
353 {
355 m_prevBasis.size(), m_prevBasis.size(), 0.0, eSYMMETRIC);
356 Vmath::Vcopy(tilCoeffMatrix->GetStorageSize(), m_coeffMatrix->GetPtr(),
357 1, tilCoeffMatrix->GetPtr(), 1);
358 }
359
360 int n, info1 = 0, info2 = 0, info3 = 0;
361 if (y_s.size())
362 {
363 n = tilCoeffMatrix->GetRows();
364 Array<OneD, NekDouble> tilCoeffMatrixFactor(
365 tilCoeffMatrix->GetStorageSize());
366 Vmath::Vcopy(tilCoeffMatrix->GetStorageSize(), tilCoeffMatrix->GetPtr(),
367 1, tilCoeffMatrixFactor, 1);
368 Lapack::Dsptrf('U', n, tilCoeffMatrixFactor.get(), ipivot.get(), info1);
369 if (info1 == 0)
370 {
371 Vmath::Vcopy(n, y_s, 1, invMy_s, 1);
372 Lapack::Dsptrs('U', n, 1, tilCoeffMatrixFactor.get(), ipivot.get(),
373 invMy_s.get(), n, info2);
374 }
375 }
376 if (info1 || info2)
377 {
378 int latest = ResetKnownSolutionsToLatestOne();
379 y_s[0] = y_s[latest - (latest > insertLocation)];
380 invMy_s[0] = y_s[0];
381 insertLocation = m_numPrevSols % m_numSuccessiveRHS;
382 fullbuffer = (m_prevBasis.size() == m_numSuccessiveRHS);
383 }
384 NekDouble residual = 1.;
385 NekDouble epsilon = 10. * NekConstants::kNekZeroTol;
386 for (int i = 0; i < m_prevBasis.size() - fullbuffer; i++)
387 {
388 residual -= y_s[i] * invMy_s[i];
389 }
390 if (m_verbose && m_root)
391 {
392 cout << "SuccessiveRHS: residual " << residual;
393 }
394 if (residual < epsilon)
395 {
396 if (m_verbose && m_root)
397 {
398 cout << " < " << epsilon << ", reject" << endl;
399 }
400 return;
401 }
402
403 // calculate new coefficient matrix and its factor
404 DNekMatSharedPtr newCoeffMatrix;
405 if (fullbuffer)
406 {
409 Vmath::Vcopy(m_coeffMatrix->GetStorageSize(), m_coeffMatrix->GetPtr(),
410 1, newCoeffMatrix->GetPtr(), 1);
411 newCoeffMatrix->SetValue(insertLocation, insertLocation, 1.);
412 for (int i = 0; i < m_numSuccessiveRHS; ++i)
413 {
414 if (i == insertLocation)
415 {
416 continue;
417 }
418 int iskip = i > insertLocation;
419 newCoeffMatrix->SetValue(insertLocation, i, y_s[i - iskip]);
420 }
421 }
422 else
423 {
425 m_prevBasis.size() + 1, m_prevBasis.size() + 1, 0.0, eSYMMETRIC);
426 newCoeffMatrix->SetValue(insertLocation, insertLocation, 1.);
427 for (int i = 0; i < m_prevBasis.size(); ++i)
428 {
429 newCoeffMatrix->SetValue(insertLocation, i, y_s[i]);
430 for (int j = i; j < m_prevBasis.size(); ++j)
431 {
432 newCoeffMatrix->SetValue(i, j, m_coeffMatrix->GetValue(i, j));
433 }
434 }
435 }
436 n = newCoeffMatrix->GetRows();
437 Array<OneD, NekDouble> coeffMatrixFactor(newCoeffMatrix->GetStorageSize());
438 Vmath::Vcopy(newCoeffMatrix->GetStorageSize(), newCoeffMatrix->GetPtr(), 1,
439 coeffMatrixFactor, 1);
440 Lapack::Dsptrf('U', n, coeffMatrixFactor.get(), ipivot.get(), info3);
441 if (info3)
442 {
443 if (m_verbose && m_root)
444 {
445 cout << " >= " << epsilon << ", reject (Dsptrf fails)" << endl;
446 }
447 return;
448 }
449 if (m_verbose && m_root)
450 {
451 cout << " >= " << epsilon << ", accept" << endl;
452 }
453
454 // if success, update basis, rhs, coefficient matrix, and its factor
455 if (m_prevBasis.size() < m_numSuccessiveRHS)
456 {
457 m_prevBasis.push_back(tmp = tmpAx_s + nDir);
458 if (isAconjugate)
459 {
460 m_prevLinSol.push_back(tmp);
461 }
462 else
463 {
464 Array<OneD, NekDouble> solution(nNonDir, 0.0);
465 m_prevLinSol.push_back(solution);
466 }
467 }
468 Vmath::Smul(nNonDir, 1. / sqrt(solNorm), tmp = newX + nDir, 1,
469 m_prevLinSol[insertLocation], 1);
470 if (!isAconjugate)
471 {
472 m_prevBasis[insertLocation] = tmpAx_s + nDir;
473 }
474 m_coeffMatrix = newCoeffMatrix;
475 m_coeffMatrixFactor = coeffMatrixFactor;
476 m_ipivot = ipivot;
478}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static void Dsptrf(const char &uplo, const int &n, double *ap, int *ipiv, int &info)
factor a real packed-symmetric matrix using Bunch-Kaufman pivoting.
Definition: Lapack.hpp:128
static const NekDouble kNekSparseNonZeroTol
static const NekDouble kNekZeroTol
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), DoMatrixMultiplyFlag(), Vmath::Dot2(), Lapack::Dsptrf(), Lapack::Dsptrs(), Nektar::eSYMMETRIC, Nektar::NekConstants::kNekSparseNonZeroTol, Nektar::NekConstants::kNekZeroTol, m_coeffMatrix, m_coeffMatrixFactor, Nektar::MultiRegions::GlobalLinSys::m_expList, m_ipivot, m_map, m_numPrevSols, m_numSuccessiveRHS, m_prevBasis, m_prevLinSol, m_root, Nektar::MultiRegions::GlobalLinSys::m_verbose, Nektar::LibUtilities::ReduceSum, ResetKnownSolutionsToLatestOne(), Vmath::Smul(), tinysimd::sqrt(), and Vmath::Vcopy().

Referenced by DoProjection().

◆ v_DoMatrixMultiply()

virtual void Nektar::MultiRegions::GlobalLinSysIterative::v_DoMatrixMultiply ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
protectedpure virtual

◆ v_UniqueMap()

virtual void Nektar::MultiRegions::GlobalLinSysIterative::v_UniqueMap ( )
protectedpure virtual

Member Data Documentation

◆ IteratSolverdef

std::string Nektar::MultiRegions::GlobalLinSysIterative::IteratSolverdef
staticprotected
Initial value:
=
"ConjugateGradient")
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.

Definition at line 110 of file GlobalLinSysIterative.h.

◆ IteratSolverlookupIds

std::string Nektar::MultiRegions::GlobalLinSysIterative::IteratSolverlookupIds
staticprotected
Initial value:
= {
"LinSysIterSolver", "ConjugateGradient",
"LinSysIterSolver", "ConjugateGradientLoc",
"LinSysIterSolver", "GMRESLoc", MultiRegions::eGMRESLoc),
}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
@ eGMRESLoc
GMRES in Local storage.
@ eConjugateGradient
Conjugate Gradient.

Definition at line 109 of file GlobalLinSysIterative.h.

◆ m_coeffMatrix

DNekMatSharedPtr Nektar::MultiRegions::GlobalLinSysIterative::m_coeffMatrix
protected

Definition at line 95 of file GlobalLinSysIterative.h.

Referenced by UpdateKnownSolutions().

◆ m_coeffMatrixFactor

Array<OneD, NekDouble> Nektar::MultiRegions::GlobalLinSysIterative::m_coeffMatrixFactor
protected

Definition at line 96 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ m_ipivot

Array<OneD, int> Nektar::MultiRegions::GlobalLinSysIterative::m_ipivot
protected

Definition at line 97 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ m_isAbsoluteTolerance

bool Nektar::MultiRegions::GlobalLinSysIterative::m_isAbsoluteTolerance
protected

◆ m_isAconjugate

bool Nektar::MultiRegions::GlobalLinSysIterative::m_isAconjugate
protected

◆ m_isNonSymmetricLinSys

bool Nektar::MultiRegions::GlobalLinSysIterative::m_isNonSymmetricLinSys
protected

Definition at line 101 of file GlobalLinSysIterative.h.

Referenced by GlobalLinSysIterative().

◆ m_linsol

LibUtilities::NekLinSysIterSharedPtr Nektar::MultiRegions::GlobalLinSysIterative::m_linsol
protected

◆ m_linSysIterSolver

std::string Nektar::MultiRegions::GlobalLinSysIterative::m_linSysIterSolver
protected

◆ m_map

Array<OneD, int> Nektar::MultiRegions::GlobalLinSysIterative::m_map
protected

◆ m_matrixType

std::string Nektar::MultiRegions::GlobalLinSysIterative::m_matrixType
protected

Definition at line 100 of file GlobalLinSysIterative.h.

Referenced by GlobalLinSysIterative().

◆ m_NekSysOp

LibUtilities::NekSysOperators Nektar::MultiRegions::GlobalLinSysIterative::m_NekSysOp
protected

◆ m_numPrevSols

int Nektar::MultiRegions::GlobalLinSysIterative::m_numPrevSols
protected

◆ m_numSuccessiveRHS

int Nektar::MultiRegions::GlobalLinSysIterative::m_numSuccessiveRHS
protected

◆ m_precon

PreconditionerSharedPtr Nektar::MultiRegions::GlobalLinSysIterative::m_precon
protected

◆ m_precontype

std::string Nektar::MultiRegions::GlobalLinSysIterative::m_precontype
protected

Definition at line 79 of file GlobalLinSysIterative.h.

◆ m_prevBasis

std::vector<Array<OneD, NekDouble> > Nektar::MultiRegions::GlobalLinSysIterative::m_prevBasis
protected

◆ m_prevLinSol

std::vector<Array<OneD, NekDouble> > Nektar::MultiRegions::GlobalLinSysIterative::m_prevLinSol
protected

Storage for solutions to previous linear problems.

Definition at line 93 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), ResetKnownSolutionsToLatestOne(), and UpdateKnownSolutions().

◆ m_rhs_mag_sm

NekDouble Nektar::MultiRegions::GlobalLinSysIterative::m_rhs_mag_sm
protected

cnt to how many times rhs_magnitude is called

Definition at line 75 of file GlobalLinSysIterative.h.

Referenced by Nektar::MultiRegions::GlobalLinSysIterativeStaticCond::v_PreSolve().

◆ m_rhs_magnitude

NekDouble Nektar::MultiRegions::GlobalLinSysIterative::m_rhs_magnitude
protected

◆ m_root

bool Nektar::MultiRegions::GlobalLinSysIterative::m_root
protected

Root if parallel.

Definition at line 87 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), GlobalLinSysIterative(), and UpdateKnownSolutions().

◆ m_tolerance

NekDouble Nektar::MultiRegions::GlobalLinSysIterative::m_tolerance
protected

◆ m_totalIterations

int Nektar::MultiRegions::GlobalLinSysIterative::m_totalIterations
protected

Definition at line 81 of file GlobalLinSysIterative.h.

◆ m_useProjection

bool Nektar::MultiRegions::GlobalLinSysIterative::m_useProjection
protected