Nektar++
GlobalLinSys.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GlobalLinSys.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: GlobalLinSys definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
40
43
45{
46std::string GlobalLinSys::lookupIds[12] = {
48 "GlobalSysSoln", "DirectFull", MultiRegions::eDirectFullMatrix),
50 "GlobalSysSoln", "DirectStaticCond", MultiRegions::eDirectStaticCond),
52 "GlobalSysSoln", "DirectMultiLevelStaticCond",
55 "GlobalSysSoln", "IterativeFull", MultiRegions::eIterativeFull),
57 "GlobalSysSoln", "IterativeStaticCond",
60 "GlobalSysSoln", "IterativeMultiLevelStaticCond",
63 "GlobalSysSoln", "XxtFull", MultiRegions::eXxtFullMatrix),
65 "GlobalSysSoln", "XxtStaticCond", MultiRegions::eXxtStaticCond),
67 "GlobalSysSoln", "XxtMultiLevelStaticCond",
70 "GlobalSysSoln", "PETScFull", MultiRegions::ePETScFullMatrix),
72 "GlobalSysSoln", "PETScStaticCond", MultiRegions::ePETScStaticCond),
74 "GlobalSysSoln", "PETScMultiLevelStaticCond",
76
77#ifdef NEKTAR_USE_SCOTCH
78std::string GlobalLinSys::def =
80 "GlobalSysSoln", "DirectMultiLevelStaticCond");
81#else
82std::string GlobalLinSys::def =
84 "DirectStaticCond");
85#endif
86
87/**
88 *
89 */
91{
92 static GlobalLinSysFactory instance;
93 return instance;
94}
95
96/**
97 * @class GlobalLinSys
98 *
99 * Consider the linear system
100 * \f$\boldsymbol{M\hat{u}}_g=\boldsymbol{\hat{f}}\f$.
101 * Distinguishing between the boundary and interior components of
102 * \f$\boldsymbol{\hat{u}}_g\f$ and \f$\boldsymbol{\hat{f}}\f$ using
103 * \f$\boldsymbol{\hat{u}}_b\f$,\f$\boldsymbol{\hat{u}}_i\f$ and
104 * \f$\boldsymbol{\hat{f}}_b\f$,\f$\boldsymbol{\hat{f}}_i\f$
105 * respectively, this system can be split into its constituent parts as
106 * \f[\left[\begin{array}{cc}
107 * \boldsymbol{M}_b&\boldsymbol{M}_{c1}\\
108 * \boldsymbol{M}_{c2}&\boldsymbol{M}_i\\
109 * \end{array}\right]
110 * \left[\begin{array}{c}
111 * \boldsymbol{\hat{u}_b}\\
112 * \boldsymbol{\hat{u}_i}\\
113 * \end{array}\right]=
114 * \left[\begin{array}{c}
115 * \boldsymbol{\hat{f}_b}\\
116 * \boldsymbol{\hat{f}_i}\\
117 * \end{array}\right]\f]
118 * where \f$\boldsymbol{M}_b\f$ represents the components of
119 * \f$\boldsymbol{M}\f$ resulting from boundary-boundary mode
120 * interactions,
121 * \f$\boldsymbol{M}_{c1}\f$ and \f$\boldsymbol{M}_{c2}\f$ represent the
122 * components resulting from coupling between the boundary-interior
123 * modes, and \f$\boldsymbol{M}_i\f$ represents the components of
124 * \f$\boldsymbol{M}\f$ resulting from interior-interior mode
125 * interactions.
126 *
127 * The solution of the linear system can now be determined in two steps:
128 * \f{eqnarray*}
129 * \mathrm{step 1:}&\quad&(\boldsymbol{M}_b-\boldsymbol{M}_{c1}
130 * \boldsymbol{M}_i^{-1}\boldsymbol{M}_{c2}) \boldsymbol{\hat{u}_b} =
131 * \boldsymbol{\hat{f}}_b - \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}
132 * \boldsymbol{\hat{f}}_i,\nonumber \\
133 * \mathrm{step 2:}&\quad&\boldsymbol{\hat{u}_i}=\boldsymbol{M}_i^{-1}
134 * \left( \boldsymbol{\hat{f}}_i
135 * - \boldsymbol{M}_{c2}\boldsymbol{\hat{u}_b}
136 * \right). \nonumber \\ \f}
137 * As the inverse of \f$\boldsymbol{M}_i^{-1}\f$ is
138 * \f[ \boldsymbol{M}_i^{-1} = \left [\underline{\boldsymbol{M}^e_i}
139 * \right ]^{-1} = \underline{[\boldsymbol{M}^e_i]}^{-1} \f]
140 * and the following operations can be evaluated as,
141 * \f{eqnarray*}
142 * \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}\boldsymbol{\hat{f}}_i &
143 * =& \boldsymbol{\mathcal{A}}_b^T \underline{\boldsymbol{M}^e_{c1}}
144 * \underline{[\boldsymbol{M}^e_i]}^{-1} \boldsymbol{\hat{f}}_i \\
145 * \boldsymbol{M}_{c2} \boldsymbol{\hat{u}_b} &=&
146 * \underline{\boldsymbol{M}^e_{c2}} \boldsymbol{\mathcal{A}}_b
147 * \boldsymbol{\hat{u}_b}.\f}
148 * where \f$\boldsymbol{\mathcal{A}}_b \f$ is the permutation matrix
149 * which scatters from global to local degrees of freedom, only the
150 * following four matrices should be constructed:
151 * - \f$\underline{[\boldsymbol{M}^e_i]}^{-1}\f$
152 * - \f$\underline{\boldsymbol{M}^e_{c1}}
153 * \underline{[\boldsymbol{M}^e_i]}^{-1}\f$
154 * - \f$\underline{\boldsymbol{M}^e_{c2}}\f$
155 * - The Schur complement: \f$\boldsymbol{M}_{\mathrm{Schur}}=
156 * \quad\boldsymbol{M}_b-\boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}
157 * \boldsymbol{M}_{c2}\f$
158 *
159 * The first three matrices are just a concatenation of the
160 * corresponding local matrices and they can be created as such. They
161 * also allow for an elemental evaluation of the operations concerned.
162 *
163 * The global Schur complement however should be assembled from the
164 * concatenation of the local elemental Schur complements, that is,
165 * \f[ \boldsymbol{M}_{\mathrm{Schur}}=\boldsymbol{M}_b
166 * - \boldsymbol{M}_{c1}
167 * \boldsymbol{M}_i^{-1} \boldsymbol{M}_{c2} =
168 * \boldsymbol{\mathcal{A}}_b^T \left [\underline{\boldsymbol{M}^e_b -
169 * \boldsymbol{M}^e_{c1} [\boldsymbol{M}^e_i]^{-1}
170 * (\boldsymbol{M}^e_{c2})} \right ] \boldsymbol{\mathcal{A}}_b \f]
171 * and it is the only matrix operation that need to be evaluated on a
172 * global level when using static condensation.
173 * However, due to the size and sparsity of the matrix
174 * \f$\boldsymbol{\mathcal{A}}_b\f$, it is more efficient to assemble
175 * the global Schur matrix using the mapping array bmap\f$[e][i]\f$
176 * contained in the input argument \a locToGloMap. The global Schur
177 * complement is then constructed as:
178 * \f[\boldsymbol{M}_{\mathrm{Schur}}\left[\mathrm{\a bmap}[e][i]\right]
179 * \left[\mathrm{\a bmap}[e][j]\right]=\mathrm{\a bsign}[e][i]\cdot
180 * \mathrm{\a bsign}[e][j]
181 * \cdot\boldsymbol{M}^e_{\mathrm{Schur}}[i][j]\f]
182 * All four matrices are stored in the \a GlobalLinSys returned by this
183 * function.
184 */
185
186/**
187 * Given a block matrix, construct a global matrix system according to
188 * a local to global mapping. #m_linSys is constructed by
189 * AssembleFullMatrix().
190 * @param pkey Associated linear system key.
191 * @param locToGloMap Local to global mapping.
192 */
194 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExpList,
195 [[maybe_unused]] const std::shared_ptr<AssemblyMap> &pLocToGloMap)
196 : m_linSysKey(pKey), m_expList(pExpList),
197 m_robinBCInfo(m_expList.lock()->GetRobinBCInfo()),
198 m_verbose(
199 m_expList.lock()->GetSession()->DefinesCmdLineArgument("verbose"))
200{
201}
202
203/**
204 * @brief Create a preconditioner object from the parameters defined in
205 * the supplied assembly map.
206 *
207 * @param asmMap Assembly map used to construct the global system.
208 */
210{
211 std::string PreconType = asmMap->GetPreconType();
212 return GetPreconFactory().CreateInstance(PreconType, GetSharedThisPtr(),
213 asmMap);
214}
215
216/**
217 * @brief Get the number of blocks in this system.
218 *
219 * At the top level this corresponds to the number of elements in the
220 * expansion list.
221 */
223{
224 return m_expList.lock()->GetExpSize();
225}
226
227/**
228 Assemble the matrix key for each block n
229**/
230
232{
233
234 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
235
236 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
237
238 // need to be initialised with zero size for non variable
239 // coefficient case
240 StdRegions::VarCoeffMap vVarCoeffMap;
241
243
244 // setup variable factors
245 if (m_linSysKey.GetNVarFactors() > 0)
246 {
247 if (m_linSysKey.GetVarFactors().count(
249 {
250 vConstFactorMap[StdRegions::eFactorSVVDiffCoeff] =
252
255 "VarCoeffSVVCuroffRatio is set but "
256 " not FactorSVVCutoffRatio");
257
258 vConstFactorMap[StdRegions::eFactorSVVCutoffRatio] =
260 }
261
262 if (m_linSysKey.GetVarFactors().count(
264 {
268 }
269
270 if (m_linSysKey.GetVarFactors().count(
272 {
273 vConstFactorMap[StdRegions::eFactorSVVDGKerDiffCoeff] =
276 }
277 }
278
279 // retrieve variable coefficients
280 if (m_linSysKey.GetNVarCoeffs() > 0)
281 {
283 expList->GetPhys_Offset(n),
284 vExp->GetTotPoints());
285 }
286
288 vExp->DetShapeType(), *vExp, vConstFactorMap,
289 vVarCoeffMap);
290 return matkey;
291}
292
293/**
294 * @brief Retrieves the block matrix from n-th expansion using the
295 * matrix key provided by the #m_linSysKey.
296 *
297 * @param n Number of the expansion.
298 * @return Block matrix for the specified expansion.
299 */
301{
302 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
303 DNekScalMatSharedPtr loc_mat;
304 loc_mat = vExp->GetLocMatrix(GetBlockMatrixKey(n));
305
306 // apply robin boundary conditions to the matrix.
307 if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
308 {
310
311 // declare local matrix from scaled matrix.
312 int rows = loc_mat->GetRows();
313 int cols = loc_mat->GetColumns();
314 const NekDouble *dat = loc_mat->GetRawPtr();
315 DNekMatSharedPtr new_mat =
317 Blas::Dscal(rows * cols, loc_mat->Scale(), new_mat->GetRawPtr(), 1);
318
319 // add local matrix contribution
320 for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
321 {
322 vExp->AddRobinMassMatrix(rBC->m_robinID,
323 rBC->m_robinPrimitiveCoeffs, new_mat);
324 }
325
326 // redeclare loc_mat to point to new_mat plus the scalar.
327 loc_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
328 }
329
330 // finally return the matrix.
331 return loc_mat;
332}
333
334/**
335 * @brief Retrieves a the static condensation block matrices from n-th
336 * expansion using the matrix key provided by the #m_linSysKey.
337 *
338 * @param n Number of the expansion
339 * @return 2x2 Block matrix holding the static condensation
340 * matrices for the n-th expansion.
341 */
343{
344
345 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
347 loc_mat = vExp->GetLocStaticCondMatrix(GetBlockMatrixKey(n));
348
349 if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
350 {
351 DNekScalMatSharedPtr tmp_mat;
353
354 tmp_mat = loc_mat->GetBlock(0, 0);
355
356 // declare local matrix from scaled matrix.
357 int rows = tmp_mat->GetRows();
358 int cols = tmp_mat->GetColumns();
359 const NekDouble *dat = tmp_mat->GetRawPtr();
360 DNekMatSharedPtr new_mat =
362 Blas::Dscal(rows * cols, tmp_mat->Scale(), new_mat->GetRawPtr(), 1);
363
364 // add local matrix contribution
365 for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
366 {
367 vExp->AddRobinMassMatrix(rBC->m_robinID,
368 rBC->m_robinPrimitiveCoeffs, new_mat);
369 }
370
371 // redeclare loc_mat to point to new_mat plus the scalar.
372 tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
373 DNekScalBlkMatSharedPtr new_loc_mat;
374 unsigned int exp_size[] = {tmp_mat->GetRows(),
375 loc_mat->GetBlock(1, 1)->GetRows()};
376 unsigned int nblks = 2;
378 nblks, nblks, exp_size, exp_size);
379
380 new_loc_mat->SetBlock(0, 0, tmp_mat);
381 new_loc_mat->SetBlock(0, 1, loc_mat->GetBlock(0, 1));
382 new_loc_mat->SetBlock(1, 0, loc_mat->GetBlock(1, 0));
383 new_loc_mat->SetBlock(1, 1, loc_mat->GetBlock(1, 1));
384 loc_mat = new_loc_mat;
385 }
386
387 return loc_mat;
388}
389
390/**
391 * @brief Releases the static condensation block matrices from NekManager
392 * of n-th expansion using the matrix key provided by the #m_linSysKey.
393 *
394 * @param n Number of the expansion
395 */
397{
398 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
399 vExp->DropLocStaticCondMatrix(GetBlockMatrixKey(n));
400}
401
402/**
403 * @brief Releases the local block matrix from NekManager
404 * of n-th expansion using the matrix key provided by the #m_linSysKey.
405 *
406 * @param n Number of the expansion
407 */
408void GlobalLinSys::v_DropBlock(unsigned int n)
409{
410 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
411 vExp->DropLocMatrix(GetBlockMatrixKey(n));
412}
413
415{
416 NEKERROR(ErrorUtil::efatal, "Method does not exist");
417}
418
420 [[maybe_unused]] const std::shared_ptr<AssemblyMap> &pLocToGloMap)
421{
422 NEKERROR(ErrorUtil::efatal, "Method does not exist");
423}
424} // namespace Nektar::MultiRegions
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Provides a generic Factory class.
Definition: NekFactory.hpp:104
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual int v_GetNumBlocks()
Get the number of blocks in this system.
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:122
static std::string lookupIds[]
Definition: GlobalLinSys.h:156
LocalRegions::MatrixKey GetBlockMatrixKey(unsigned int n)
std::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:100
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
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.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:120
virtual DNekScalMatSharedPtr v_GetBlock(unsigned int n)
Retrieves the block matrix from n-th expansion using the matrix key provided by the m_linSysKey.
virtual void v_DropBlock(unsigned int n)
Releases the local block matrix from NekManager of n-th expansion using the matrix key provided by th...
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map.
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n)
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...
virtual void v_Initialise(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
virtual void v_DropStaticCondBlock(unsigned int n)
Releases the static condensation block matrices from NekManager of n-th expansion using the matrix ke...
const Array< OneD, const NekDouble > & GetVarFactors(const StdRegions::ConstFactorType &coeff) const
const StdRegions::ConstFactorMap & GetConstFactors() const
Returns all the constants.
const StdRegions::VarCoeffMap & GetVarCoeffs() const
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:58
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
VarCoeffMap RestrictCoeffMap(const VarCoeffMap &m, size_t offset, size_t cnt)
Definition: StdRegions.hpp:349
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble