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
35#include <boost/core/ignore_unused.hpp>
36
42
45
46namespace Nektar
47{
48namespace MultiRegions
49{
50std::string GlobalLinSys::lookupIds[12] = {
52 "GlobalSysSoln", "DirectFull", MultiRegions::eDirectFullMatrix),
54 "GlobalSysSoln", "DirectStaticCond", MultiRegions::eDirectStaticCond),
56 "GlobalSysSoln", "DirectMultiLevelStaticCond",
59 "GlobalSysSoln", "IterativeFull", MultiRegions::eIterativeFull),
61 "GlobalSysSoln", "IterativeStaticCond",
64 "GlobalSysSoln", "IterativeMultiLevelStaticCond",
67 "GlobalSysSoln", "XxtFull", MultiRegions::eXxtFullMatrix),
69 "GlobalSysSoln", "XxtStaticCond", MultiRegions::eXxtStaticCond),
71 "GlobalSysSoln", "XxtMultiLevelStaticCond",
74 "GlobalSysSoln", "PETScFull", MultiRegions::ePETScFullMatrix),
76 "GlobalSysSoln", "PETScStaticCond", MultiRegions::ePETScStaticCond),
78 "GlobalSysSoln", "PETScMultiLevelStaticCond",
80
81#ifdef NEKTAR_USE_SCOTCH
82std::string GlobalLinSys::def =
84 "GlobalSysSoln", "DirectMultiLevelStaticCond");
85#else
86std::string GlobalLinSys::def =
88 "DirectStaticCond");
89#endif
90
91/**
92 *
93 */
95{
96 static GlobalLinSysFactory instance;
97 return instance;
98}
99
100/**
101 * @class GlobalLinSys
102 *
103 * Consider the linear system
104 * \f$\boldsymbol{M\hat{u}}_g=\boldsymbol{\hat{f}}\f$.
105 * Distinguishing between the boundary and interior components of
106 * \f$\boldsymbol{\hat{u}}_g\f$ and \f$\boldsymbol{\hat{f}}\f$ using
107 * \f$\boldsymbol{\hat{u}}_b\f$,\f$\boldsymbol{\hat{u}}_i\f$ and
108 * \f$\boldsymbol{\hat{f}}_b\f$,\f$\boldsymbol{\hat{f}}_i\f$
109 * respectively, this system can be split into its constituent parts as
110 * \f[\left[\begin{array}{cc}
111 * \boldsymbol{M}_b&\boldsymbol{M}_{c1}\\
112 * \boldsymbol{M}_{c2}&\boldsymbol{M}_i\\
113 * \end{array}\right]
114 * \left[\begin{array}{c}
115 * \boldsymbol{\hat{u}_b}\\
116 * \boldsymbol{\hat{u}_i}\\
117 * \end{array}\right]=
118 * \left[\begin{array}{c}
119 * \boldsymbol{\hat{f}_b}\\
120 * \boldsymbol{\hat{f}_i}\\
121 * \end{array}\right]\f]
122 * where \f$\boldsymbol{M}_b\f$ represents the components of
123 * \f$\boldsymbol{M}\f$ resulting from boundary-boundary mode
124 * interactions,
125 * \f$\boldsymbol{M}_{c1}\f$ and \f$\boldsymbol{M}_{c2}\f$ represent the
126 * components resulting from coupling between the boundary-interior
127 * modes, and \f$\boldsymbol{M}_i\f$ represents the components of
128 * \f$\boldsymbol{M}\f$ resulting from interior-interior mode
129 * interactions.
130 *
131 * The solution of the linear system can now be determined in two steps:
132 * \f{eqnarray*}
133 * \mathrm{step 1:}&\quad&(\boldsymbol{M}_b-\boldsymbol{M}_{c1}
134 * \boldsymbol{M}_i^{-1}\boldsymbol{M}_{c2}) \boldsymbol{\hat{u}_b} =
135 * \boldsymbol{\hat{f}}_b - \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}
136 * \boldsymbol{\hat{f}}_i,\nonumber \\
137 * \mathrm{step 2:}&\quad&\boldsymbol{\hat{u}_i}=\boldsymbol{M}_i^{-1}
138 * \left( \boldsymbol{\hat{f}}_i
139 * - \boldsymbol{M}_{c2}\boldsymbol{\hat{u}_b}
140 * \right). \nonumber \\ \f}
141 * As the inverse of \f$\boldsymbol{M}_i^{-1}\f$ is
142 * \f[ \boldsymbol{M}_i^{-1} = \left [\underline{\boldsymbol{M}^e_i}
143 * \right ]^{-1} = \underline{[\boldsymbol{M}^e_i]}^{-1} \f]
144 * and the following operations can be evaluated as,
145 * \f{eqnarray*}
146 * \boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}\boldsymbol{\hat{f}}_i &
147 * =& \boldsymbol{\mathcal{A}}_b^T \underline{\boldsymbol{M}^e_{c1}}
148 * \underline{[\boldsymbol{M}^e_i]}^{-1} \boldsymbol{\hat{f}}_i \\
149 * \boldsymbol{M}_{c2} \boldsymbol{\hat{u}_b} &=&
150 * \underline{\boldsymbol{M}^e_{c2}} \boldsymbol{\mathcal{A}}_b
151 * \boldsymbol{\hat{u}_b}.\f}
152 * where \f$\boldsymbol{\mathcal{A}}_b \f$ is the permutation matrix
153 * which scatters from global to local degrees of freedom, only the
154 * following four matrices should be constructed:
155 * - \f$\underline{[\boldsymbol{M}^e_i]}^{-1}\f$
156 * - \f$\underline{\boldsymbol{M}^e_{c1}}
157 * \underline{[\boldsymbol{M}^e_i]}^{-1}\f$
158 * - \f$\underline{\boldsymbol{M}^e_{c2}}\f$
159 * - The Schur complement: \f$\boldsymbol{M}_{\mathrm{Schur}}=
160 * \quad\boldsymbol{M}_b-\boldsymbol{M}_{c1}\boldsymbol{M}_i^{-1}
161 * \boldsymbol{M}_{c2}\f$
162 *
163 * The first three matrices are just a concatenation of the
164 * corresponding local matrices and they can be created as such. They
165 * also allow for an elemental evaluation of the operations concerned.
166 *
167 * The global Schur complement however should be assembled from the
168 * concatenation of the local elemental Schur complements, that is,
169 * \f[ \boldsymbol{M}_{\mathrm{Schur}}=\boldsymbol{M}_b
170 * - \boldsymbol{M}_{c1}
171 * \boldsymbol{M}_i^{-1} \boldsymbol{M}_{c2} =
172 * \boldsymbol{\mathcal{A}}_b^T \left [\underline{\boldsymbol{M}^e_b -
173 * \boldsymbol{M}^e_{c1} [\boldsymbol{M}^e_i]^{-1}
174 * (\boldsymbol{M}^e_{c2})} \right ] \boldsymbol{\mathcal{A}}_b \f]
175 * and it is the only matrix operation that need to be evaluated on a
176 * global level when using static condensation.
177 * However, due to the size and sparsity of the matrix
178 * \f$\boldsymbol{\mathcal{A}}_b\f$, it is more efficient to assemble
179 * the global Schur matrix using the mapping array bmap\f$[e][i]\f$
180 * contained in the input argument \a locToGloMap. The global Schur
181 * complement is then constructed as:
182 * \f[\boldsymbol{M}_{\mathrm{Schur}}\left[\mathrm{\a bmap}[e][i]\right]
183 * \left[\mathrm{\a bmap}[e][j]\right]=\mathrm{\a bsign}[e][i]\cdot
184 * \mathrm{\a bsign}[e][j]
185 * \cdot\boldsymbol{M}^e_{\mathrm{Schur}}[i][j]\f]
186 * All four matrices are stored in the \a GlobalLinSys returned by this
187 * function.
188 */
189
190/**
191 * Given a block matrix, construct a global matrix system according to
192 * a local to global mapping. #m_linSys is constructed by
193 * AssembleFullMatrix().
194 * @param pkey Associated linear system key.
195 * @param locToGloMap Local to global mapping.
196 */
198 const std::weak_ptr<ExpList> &pExpList,
199 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
200 : m_linSysKey(pKey), m_expList(pExpList),
201 m_robinBCInfo(m_expList.lock()->GetRobinBCInfo()),
202 m_verbose(
203 m_expList.lock()->GetSession()->DefinesCmdLineArgument("verbose"))
204{
205 boost::ignore_unused(pLocToGloMap);
206}
207
208/**
209 * @brief Create a preconditioner object from the parameters defined in
210 * the supplied assembly map.
211 *
212 * @param asmMap Assembly map used to construct the global system.
213 */
215{
216 std::string PreconType = asmMap->GetPreconType();
217 return GetPreconFactory().CreateInstance(PreconType, GetSharedThisPtr(),
218 asmMap);
219}
220
221/**
222 * @brief Get the number of blocks in this system.
223 *
224 * At the top level this corresponds to the number of elements in the
225 * expansion list.
226 */
228{
229 return m_expList.lock()->GetExpSize();
230}
231
232/**
233 Assemble the matrix key for each block n
234**/
235
237{
238
239 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
240
241 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
242
243 // need to be initialised with zero size for non variable
244 // coefficient case
245 StdRegions::VarCoeffMap vVarCoeffMap;
246
248
249 // setup variable factors
250 if (m_linSysKey.GetNVarFactors() > 0)
251 {
252 if (m_linSysKey.GetVarFactors().count(
254 {
255 vConstFactorMap[StdRegions::eFactorSVVDiffCoeff] =
257
260 "VarCoeffSVVCuroffRatio is set but "
261 " not FactorSVVCutoffRatio");
262
263 vConstFactorMap[StdRegions::eFactorSVVCutoffRatio] =
265 }
266
267 if (m_linSysKey.GetVarFactors().count(
269 {
273 }
274
275 if (m_linSysKey.GetVarFactors().count(
277 {
278 vConstFactorMap[StdRegions::eFactorSVVDGKerDiffCoeff] =
281 }
282 }
283
284 // retrieve variable coefficients
285 if (m_linSysKey.GetNVarCoeffs() > 0)
286 {
288 expList->GetPhys_Offset(n),
289 vExp->GetTotPoints());
290 }
291
293 vExp->DetShapeType(), *vExp, vConstFactorMap,
294 vVarCoeffMap);
295 return matkey;
296}
297
298/**
299 * @brief Retrieves the block matrix from n-th expansion using the
300 * matrix key provided by the #m_linSysKey.
301 *
302 * @param n Number of the expansion.
303 * @return Block matrix for the specified expansion.
304 */
306{
307 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
308 DNekScalMatSharedPtr loc_mat;
309 loc_mat = vExp->GetLocMatrix(GetBlockMatrixKey(n));
310
311 // apply robin boundary conditions to the matrix.
312 if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
313 {
315
316 // declare local matrix from scaled matrix.
317 int rows = loc_mat->GetRows();
318 int cols = loc_mat->GetColumns();
319 const NekDouble *dat = loc_mat->GetRawPtr();
320 DNekMatSharedPtr new_mat =
322 Blas::Dscal(rows * cols, loc_mat->Scale(), new_mat->GetRawPtr(), 1);
323
324 // add local matrix contribution
325 for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
326 {
327 vExp->AddRobinMassMatrix(rBC->m_robinID,
328 rBC->m_robinPrimitiveCoeffs, new_mat);
329 }
330
331 // redeclare loc_mat to point to new_mat plus the scalar.
332 loc_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
333 }
334
335 // finally return the matrix.
336 return loc_mat;
337}
338
339/**
340 * @brief Retrieves a the static condensation block matrices from n-th
341 * expansion using the matrix key provided by the #m_linSysKey.
342 *
343 * @param n Number of the expansion
344 * @return 2x2 Block matrix holding the static condensation
345 * matrices for the n-th expansion.
346 */
348{
349
350 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
352 loc_mat = vExp->GetLocStaticCondMatrix(GetBlockMatrixKey(n));
353
354 if (m_robinBCInfo.count(n) != 0) // add robin mass matrix
355 {
356 DNekScalMatSharedPtr tmp_mat;
358
359 tmp_mat = loc_mat->GetBlock(0, 0);
360
361 // declare local matrix from scaled matrix.
362 int rows = tmp_mat->GetRows();
363 int cols = tmp_mat->GetColumns();
364 const NekDouble *dat = tmp_mat->GetRawPtr();
365 DNekMatSharedPtr new_mat =
367 Blas::Dscal(rows * cols, tmp_mat->Scale(), new_mat->GetRawPtr(), 1);
368
369 // add local matrix contribution
370 for (rBC = m_robinBCInfo.find(n)->second; rBC; rBC = rBC->next)
371 {
372 vExp->AddRobinMassMatrix(rBC->m_robinID,
373 rBC->m_robinPrimitiveCoeffs, new_mat);
374 }
375
376 // redeclare loc_mat to point to new_mat plus the scalar.
377 tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, new_mat);
378 DNekScalBlkMatSharedPtr new_loc_mat;
379 unsigned int exp_size[] = {tmp_mat->GetRows(),
380 loc_mat->GetBlock(1, 1)->GetRows()};
381 unsigned int nblks = 2;
383 nblks, nblks, exp_size, exp_size);
384
385 new_loc_mat->SetBlock(0, 0, tmp_mat);
386 new_loc_mat->SetBlock(0, 1, loc_mat->GetBlock(0, 1));
387 new_loc_mat->SetBlock(1, 0, loc_mat->GetBlock(1, 0));
388 new_loc_mat->SetBlock(1, 1, loc_mat->GetBlock(1, 1));
389 loc_mat = new_loc_mat;
390 }
391
392 return loc_mat;
393}
394
395/**
396 * @brief Releases the static condensation block matrices from NekManager
397 * of n-th expansion using the matrix key provided by the #m_linSysKey.
398 *
399 * @param n Number of the expansion
400 */
402{
403 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
404 vExp->DropLocStaticCondMatrix(GetBlockMatrixKey(n));
405}
406
407/**
408 * @brief Releases the local block matrix from NekManager
409 * of n-th expansion using the matrix key provided by the #m_linSysKey.
410 *
411 * @param n Number of the expansion
412 */
413void GlobalLinSys::v_DropBlock(unsigned int n)
414{
415 LocalRegions::ExpansionSharedPtr vExp = m_expList.lock()->GetExp(n);
416 vExp->DropLocMatrix(GetBlockMatrixKey(n));
417}
418
420{
421 NEKERROR(ErrorUtil::efatal, "Method does not exist");
422}
423
425 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
426{
427 boost::ignore_unused(pLocToGloMap);
428 NEKERROR(ErrorUtil::efatal, "Method does not exist");
429}
430} // namespace MultiRegions
431} // namespace Nektar
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Provides a generic Factory class.
Definition: NekFactory.hpp:105
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
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:124
static std::string lookupIds[]
Definition: GlobalLinSys.h:158
LocalRegions::MatrixKey GetBlockMatrixKey(unsigned int n)
std::shared_ptr< GlobalLinSys > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: GlobalLinSys.h:102
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:126
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:122
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:151
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:60
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
VarCoeffMap RestrictCoeffMap(const VarCoeffMap &m, size_t offset, size_t cnt)
Definition: StdRegions.hpp:355
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:352
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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