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 bool isAconjugate)
 projection technique More...
 
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_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_isAbsoluteTolerance = pLocToGloMap->IsAbsoluteTolerance();
74 m_linSysIterSolver = pLocToGloMap->GetLinSysIterSolver();
75
77 m_expList.lock()->GetComm()->GetRowComm();
78 m_root = (vComm->GetRank()) ? false : true;
79
80 m_numSuccessiveRHS = pLocToGloMap->GetSuccessiveRHS();
84
85 // Check for advection matrix and switch to GMRES, if not already used
86 m_matrixType = StdRegions::MatrixTypeMap[pKey.GetMatrixType()];
88 m_matrixType.find("AdvectionDiffusionReaction") != string::npos;
89
91 !m_linSysIterSolver.compare("ConjugateGradient"))
92 {
93 m_linSysIterSolver = "GMRES";
95 false,
96 "Detected ConjugateGradient solver and a "
97 "Advection-Diffusion-Reaction matrix. "
98 "Switchted to a GMRES solver for this non-symmetric matrix type. "
99 "Change LinSysIterSolver to GMRES in the session file to suppress "
100 "this warning.");
101 }
102
103 if (m_isAconjugate && m_linSysIterSolver.compare("GMRES") == 0 &&
104 m_linSysIterSolver.compare("GMRESLoc") == 0)
105 {
106 WARNINGL0(false, "To use A-conjugate projection, the matrix "
107 "should be symmetric positive definite.");
108 }
109
111 this);
113 this);
115 this);
116}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:121
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:107
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:114
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_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_useProjection, Nektar::StdRegions::MatrixTypeMap, and WARNINGL0.

◆ ~GlobalLinSysIterative()

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

Definition at line 118 of file GlobalLinSysIterative.cpp.

119{
120}

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 142 of file GlobalLinSysIterative.h.

144 {
145 m_precon->DoAssembleLoc(pInput, pOutput, ZeroDir);
146 }

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 135 of file GlobalLinSysIterative.h.

138 {
139 v_DoMatrixMultiply(pInput, pOutput);
140 }

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 128 of file GlobalLinSysIterative.h.

131 {
132 m_precon->DoPreconditioner(pInput, pOutput, isLocal);
133 }

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 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 128 of file GlobalLinSysIterative.cpp.

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

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

References m_numPrevSols, m_numSuccessiveRHS, m_prevBasis, and m_prevLinSol.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ 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 269 of file GlobalLinSysIterative.cpp.

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

◆ m_coeffMatrix

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

Definition at line 92 of file GlobalLinSysIterative.h.

Referenced by UpdateKnownSolutions().

◆ m_coeffMatrixFactor

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

Definition at line 93 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ m_ipivot

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

Definition at line 94 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 98 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 97 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 76 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 90 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 72 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 84 of file GlobalLinSysIterative.h.

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

◆ m_totalIterations

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

Definition at line 78 of file GlobalLinSysIterative.h.

◆ m_useProjection

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