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...
 
virtual ~GlobalLinSysIterative ()
 
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...
 
int m_maxiter
 maximum iterations 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 50 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 67 of file GlobalLinSysIterative.cpp.

70 : GlobalLinSys(pKey, pExpList, pLocToGloMap),
74{
75 m_tolerance = pLocToGloMap->GetIterativeTolerance();
76 m_isAbsoluteTolerance = pLocToGloMap->IsAbsoluteTolerance();
77 m_maxiter = pLocToGloMap->GetMaxIterations();
78 m_linSysIterSolver = pLocToGloMap->GetLinSysIterSolver();
79
81 m_expList.lock()->GetComm()->GetRowComm();
82 m_root = (vComm->GetRank()) ? false : true;
83
84 m_numSuccessiveRHS = pLocToGloMap->GetSuccessiveRHS();
88
89 // Check for advection matrix and switch to GMRES, if not already used
90 m_matrixType = StdRegions::MatrixTypeMap[pKey.GetMatrixType()];
92 m_matrixType.find("AdvectionDiffusionReaction") != string::npos;
93
95 !m_linSysIterSolver.compare("ConjugateGradient"))
96 {
97 m_linSysIterSolver = "GMRES";
99 false,
100 "Detected ConjugateGradient solver and a "
101 "Advection-Diffusion-Reaction matrix. "
102 "Switchted to a GMRES solver for this non-symmetric matrix type. "
103 "Change LinSysIterSolver to GMRES in the session file to suppress "
104 "this warning.");
105 }
106
107 if (m_isAconjugate && m_linSysIterSolver.compare("GMRES") == 0 &&
108 m_linSysIterSolver.compare("GMRESLoc") == 0)
109 {
110 WARNINGL0(false, "To use A-conjugate projection, the matrix "
111 "should be symmetric positive definite.");
112 }
113
115 this);
117 this);
119 this);
120}
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:222
void DefineAssembleLoc(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:122
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:108
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition: NekSys.h:115
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:124
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:57
static PreconditionerSharedPtr NullPreconditionerSharedPtr
static const NekDouble kNekUnsetDouble
const char *const MatrixTypeMap[]
Definition: StdRegions.hpp:145
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_maxiter, m_NekSysOp, m_numSuccessiveRHS, m_root, m_tolerance, m_useProjection, Nektar::StdRegions::MatrixTypeMap, and WARNINGL0.

◆ ~GlobalLinSysIterative()

Nektar::MultiRegions::GlobalLinSysIterative::~GlobalLinSysIterative ( )
virtual

Definition at line 122 of file GlobalLinSysIterative.cpp.

123{
124}

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

156 {
157 m_precon->DoAssembleLoc(pInput, pOutput, ZeroDir);
158 }

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

62 {
63 v_DoMatrixMultiply(pInput, pOutput);
64 }
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 145 of file GlobalLinSysIterative.h.

148 {
149 boost::ignore_unused(controlFlag);
150
151 v_DoMatrixMultiply(pInput, pOutput);
152 }

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

141 {
142 m_precon->DoPreconditioner(pInput, pOutput, isLocal);
143 }

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

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

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

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

Definition at line 463 of file GlobalLinSysIterative.cpp.

464{
466 {
467 m_rhs_magnitude = 1.0;
468 return;
469 }
470
471 NekDouble vExchange(0.0);
472 if (m_map.size() > 0)
473 {
474 vExchange =
475 Vmath::Dot2(pIn.GetDimension(), &pIn[0], &pIn[0], &m_map[0]);
476 }
477
478 m_expList.lock()->GetComm()->GetRowComm()->AllReduce(
480
481 // To ensure that very different rhs values are not being
482 // used in subsequent solvers such as the velocit solve in
483 // INC NS. If this works we then need to work out a better
484 // way to control this.
485 NekDouble new_rhs_mag = (vExchange > 1e-6) ? vExchange : 1.0;
486
488 {
489 m_rhs_magnitude = new_rhs_mag;
490 }
491 else
492 {
494 (1.0 - m_rhs_mag_sm) * new_rhs_mag);
495 }
496}
unsigned int GetDimension() const
Returns the number of dimensions for the point.
Definition: NekVector.cpp:201

References Vmath::Dot2(), Nektar::NekVector< DataType >::GetDimension(), Nektar::NekConstants::kNekUnsetDouble, Nektar::MultiRegions::GlobalLinSys::m_expList, m_isAbsoluteTolerance, m_map, m_rhs_mag_sm, m_rhs_magnitude, and Nektar::LibUtilities::ReduceSum.

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

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

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

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

◆ m_coeffMatrix

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

Definition at line 100 of file GlobalLinSysIterative.h.

Referenced by UpdateKnownSolutions().

◆ m_coeffMatrixFactor

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

Definition at line 101 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ m_ipivot

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

Definition at line 102 of file GlobalLinSysIterative.h.

Referenced by DoProjection(), and UpdateKnownSolutions().

◆ m_isAbsoluteTolerance

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

Definition at line 108 of file GlobalLinSysIterative.h.

Referenced by GlobalLinSysIterative(), and Set_Rhs_Magnitude().

◆ m_isAconjugate

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

◆ m_isNonSymmetricLinSys

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

Definition at line 106 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 105 of file GlobalLinSysIterative.h.

Referenced by GlobalLinSysIterative().

◆ m_maxiter

int Nektar::MultiRegions::GlobalLinSysIterative::m_maxiter
protected

maximum iterations

Definition at line 71 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 84 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 98 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 80 of file GlobalLinSysIterative.h.

Referenced by Set_Rhs_Magnitude(), and 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 92 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 86 of file GlobalLinSysIterative.h.

◆ m_useProjection

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