Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GlobalLinSysIterativeStaticCond.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: GlobalLinSysIterativeStaticCond.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Implementation to linear solver using single-
33 // or multi-level static condensation
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
43 
44 namespace Nektar
45 {
46  namespace MultiRegions
47  {
48  /**
49  * @class GlobalLinSysIterativeStaticCond
50  *
51  * Solves a linear system iteratively using single- or multi-level
52  * static condensation.
53  */
54 
55  /**
56  * Registers the class with the Factory.
57  */
60  "IterativeStaticCond",
62  "Iterative static condensation.");
63 
66  "IterativeMultiLevelStaticCond",
68  "Iterative multi-level static condensation.");
69 
70 
73  "LocalMatrixStorageStrategy",
74  "Sparse");
77  "LocalMatrixStorageStrategy",
78  "Contiguous",
81  "LocalMatrixStorageStrategy",
82  "Non-contiguous",
85  "LocalMatrixStorageStrategy",
86  "Sparse",
88  };
89 
90  /**
91  * For a matrix system of the form @f[
92  * \left[ \begin{array}{cc}
93  * \boldsymbol{A} & \boldsymbol{B}\\
94  * \boldsymbol{C} & \boldsymbol{D}
95  * \end{array} \right]
96  * \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2}
97  * \end{array}\right]
98  * = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2}
99  * \end{array}\right],
100  * @f]
101  * where @f$\boldsymbol{D}@f$ and
102  * @f$(\boldsymbol{A-BD^{-1}C})@f$ are invertible, store and assemble
103  * a static condensation system, according to a given local to global
104  * mapping. #m_linSys is constructed by AssembleSchurComplement().
105  * @param mKey Associated matrix key.
106  * @param pLocMatSys LocalMatrixSystem
107  * @param locToGloMap Local to global mapping.
108  */
110  const GlobalLinSysKey &pKey,
111  const boost::weak_ptr<ExpList> &pExpList,
112  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
113  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
114  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
115  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
116  {
119  "This constructor is only valid when using static "
120  "condensation");
122  == pLocToGloMap->GetGlobalSysSolnType(),
123  "The local to global map is not set up for the requested "
124  "solution type");
125  }
126 
127 
128  /**
129  *
130  */
132  const GlobalLinSysKey &pKey,
133  const boost::weak_ptr<ExpList> &pExpList,
134  const DNekScalBlkMatSharedPtr pSchurCompl,
135  const DNekScalBlkMatSharedPtr pBinvD,
136  const DNekScalBlkMatSharedPtr pC,
137  const DNekScalBlkMatSharedPtr pInvD,
138  const boost::shared_ptr<AssemblyMap> &pLocToGloMap,
139  const PreconditionerSharedPtr pPrecon)
140  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
141  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
142  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
143  {
144  m_schurCompl = pSchurCompl;
145  m_S1Blk = pSchurCompl;
146  m_BinvD = pBinvD;
147  m_C = pC;
148  m_invD = pInvD;
149  m_precon = pPrecon;
150  }
151 
152 
154  {
156 
157  // Allocate memory for top-level structure
159 
160  // Setup Block Matrix systems
161  int n, n_exp = m_expList.lock()->GetNumElmts();
162 
163  MatrixStorage blkmatStorage = eDIAGONAL;
164  const Array<OneD,const unsigned int>& nbdry_size
165  = m_locToGloMap->GetNumLocalBndCoeffsPerPatch();
166 
168  ::AllocateSharedPtr(nbdry_size, nbdry_size, blkmatStorage);
169 
170  // Preserve original matrix in m_S1Blk
171  for (n = 0; n < n_exp; ++n)
172  {
173  DNekScalMatSharedPtr mat = m_schurCompl->GetBlock(n, n);
174  m_S1Blk->SetBlock(n, n, mat);
175  }
176 
177  // Build preconditioner
178  m_precon->BuildPreconditioner();
179 
180  // Do transform of Schur complement matrix
181  for (n = 0; n < n_exp; ++n)
182  {
183  if (m_linSysKey.GetMatrixType() !=
185  {
186  DNekScalMatSharedPtr mat = m_S1Blk->GetBlock(n, n);
187  DNekScalMatSharedPtr t = m_precon->TransformedSchurCompl(
188  m_expList.lock()->GetOffset_Elmt_Id(n), mat);
189  m_schurCompl->SetBlock(n, n, t);
190  }
191  }
192 
193  // Construct this level
195  }
196 
197  /**
198  *
199  */
201  {
202 
203  }
204 
206  v_GetStaticCondBlock(unsigned int n)
207  {
208  DNekScalBlkMatSharedPtr schurComplBlock;
209  int scLevel = m_locToGloMap->GetStaticCondLevel();
210  DNekScalBlkMatSharedPtr sc = scLevel == 0 ? m_S1Blk : m_schurCompl;
211  DNekScalMatSharedPtr localMat = sc->GetBlock(n,n);
212  unsigned int nbdry = localMat->GetRows();
213  unsigned int nblks = 1;
214  unsigned int esize[1] = {nbdry};
215 
216  schurComplBlock = MemoryManager<DNekScalBlkMat>
217  ::AllocateSharedPtr(nblks, nblks, esize, esize);
218  schurComplBlock->SetBlock(0, 0, localMat);
219 
220  return schurComplBlock;
221  }
222 
223  /**
224  * Assemble the schur complement matrix from the block matrices stored
225  * in #m_blkMatrices and the given local to global mapping information.
226  * @param locToGloMap Local to global mapping information.
227  */
229  const AssemblyMapSharedPtr pLocToGloMap)
230  {
231  int i,j,n,cnt,gid1,gid2;
232  NekDouble sign1,sign2;
233 
234  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
235  DoGlobalMatOp(m_linSysKey.GetMatrixType());
236 
237  // Set up unique map
238  v_UniqueMap();
239 
240  // Build precon again if we in multi-level static condensation (a
241  // bit of a hack)
244  {
246  m_precon->BuildPreconditioner();
247  }
248 
249  if (!doGlobalOp)
250  {
252  return;
253  }
254 
255  int nBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
256  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
257  unsigned int rows = nBndDofs - NumDirBCs;
258  unsigned int cols = nBndDofs - NumDirBCs;
259 
260  // COO sparse storage to assist in assembly
261  COOMatType gmat_coo;
262 
263  // Get the matrix storage structure
264  // (whether to store only one triangular part, if symmetric)
265  MatrixStorage matStorage = eFULL;
266 
267  // assemble globally
268  DNekScalMatSharedPtr loc_mat;
269  int loc_lda;
270  for(n = cnt = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
271  {
272  loc_mat = m_schurCompl->GetBlock(n,n);
273  loc_lda = loc_mat->GetRows();
274 
275  // Set up Matrix;
276  for(i = 0; i < loc_lda; ++i)
277  {
278  gid1 = pLocToGloMap->GetLocalToGlobalBndMap (cnt + i)
279  - NumDirBCs;
280  sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
281 
282  if(gid1 >= 0)
283  {
284  for(j = 0; j < loc_lda; ++j)
285  {
286  gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt+j)
287  - NumDirBCs;
288  sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt+j);
289 
290  if (gid2 >= 0)
291  {
292  gmat_coo[std::make_pair(gid1,gid2)] +=
293  sign1*sign2*(*loc_mat)(i,j);
294  }
295  }
296  }
297  }
298  cnt += loc_lda;
299  }
300 
302  sparseStorage (1);
303 
304  BCOMatType partMat;
305  convertCooToBco(rows, cols, 1, gmat_coo, partMat);
306 
307  sparseStorage[0] =
309  AllocateSharedPtr(rows, cols, 1, partMat, matStorage );
310 
311  // Create block diagonal matrix
313  AllocateSharedPtr(sparseStorage);
314  }
315 
316 
317  /**
318  * Populates sparse block-diagonal schur complement matrix from
319  * the block matrices stored in #m_blkMatrices.
320  */
322  {
323  LocalMatrixStorageStrategy storageStrategy =
324  m_expList.lock()->GetSession()->
325  GetSolverInfoAsEnum<LocalMatrixStorageStrategy>(
326  "LocalMatrixStorageStrategy");
327 
328  switch(storageStrategy)
329  {
332  {
333  size_t storageSize = 0;
334  int nBlk = m_schurCompl->GetNumberOfBlockRows();
335 
336  m_scale = Array<OneD, NekDouble> (nBlk, 1.0);
337  m_rows = Array<OneD, unsigned int> (nBlk, 0U);
338 
339  // Determine storage requirements for dense blocks.
340  for (int i = 0; i < nBlk; ++i)
341  {
342  m_rows[i] = m_schurCompl->GetBlock(i,i)->GetRows();
343  m_scale[i] = m_schurCompl->GetBlock(i,i)->Scale();
344  storageSize += m_rows[i] * m_rows[i];
345  }
346 
347  // Assemble dense storage blocks.
348  DNekScalMatSharedPtr loc_mat;
349  m_denseBlocks.resize(nBlk);
350  double *ptr = 0;
351 
352  if (MultiRegions::eContiguous == storageStrategy)
353  {
354  m_storage.resize (storageSize);
355  ptr = &m_storage[0];
356  }
357 
358  for (unsigned int n = 0; n < nBlk; ++n)
359  {
360  loc_mat = m_schurCompl->GetBlock(n,n);
361 
362  if (MultiRegions::eContiguous == storageStrategy)
363  {
364  int loc_lda = loc_mat->GetRows();
365  int blockSize = loc_lda * loc_lda;
366  m_denseBlocks[n] = ptr;
367  for(int i = 0; i < loc_lda; ++i)
368  {
369  for(int j = 0; j < loc_lda; ++j)
370  {
371  ptr[j*loc_lda+i] = (*loc_mat)(i,j);
372  }
373  }
374  ptr += blockSize;
376  m_expList.lock()->GetOffset_Elmt_Id(n));
377  }
378  else
379  {
380  m_denseBlocks[n] = loc_mat->GetRawPtr();
381  }
382  }
383  break;
384  }
386  {
387  DNekScalMatSharedPtr loc_mat;
388  int loc_lda;
389  int blockSize = 0;
390 
391  // First run through to split the set of local matrices into
392  // partitions of fixed block size, and count number of local
393  // matrices that belong to each partition.
394  std::vector<std::pair<int,int> > partitions;
395  for(int n = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
396  {
397  loc_mat = m_schurCompl->GetBlock(n,n);
398  loc_lda = loc_mat->GetRows();
399 
400  ASSERTL1(loc_lda >= 0,
401  boost::lexical_cast<std::string>(n) + "-th "
402  "matrix block in Schur complement has "
403  "rank 0!");
404 
405  if (blockSize == loc_lda)
406  {
407  partitions[partitions.size()-1].first++;
408  }
409  else
410  {
411  blockSize = loc_lda;
412  partitions.push_back(make_pair(1,loc_lda));
413  }
414  }
415 
416  MatrixStorage matStorage = eFULL;
417 
418  // Create a vector of sparse storage holders
420  sparseStorage (partitions.size());
421 
422  for (int part = 0, n = 0; part < partitions.size(); ++part)
423  {
424  BCOMatType partMat;
425 
426  for(int k = 0; k < partitions[part].first; ++k, ++n)
427  {
428  loc_mat = m_schurCompl->GetBlock(n,n);
429  loc_lda = loc_mat->GetRows();
430 
431  ASSERTL1(loc_lda == partitions[part].second,
432  boost::lexical_cast<std::string>(n) + "-th"
433  " matrix block in Schur complement has "
434  "unexpected rank");
435 
436  NekDouble scale = loc_mat->Scale();
437  if(fabs(scale-1.0) > NekConstants::kNekZeroTol)
438  {
439  Array<OneD, NekDouble> matarray(loc_lda*loc_lda);
440  Vmath::Smul(loc_lda*loc_lda,scale,
441  loc_mat->GetRawPtr(),1,&matarray[0],1);
442  partMat[make_pair(k,k)] = BCOEntryType(matarray);
443  }
444  else // scale factor is 1.0
445  {
446  partMat[make_pair(k,k)] = BCOEntryType(
447  loc_lda*loc_lda, loc_mat->GetRawPtr());
448  }
449 
451  m_expList.lock()->GetOffset_Elmt_Id(n));
452  }
453 
454  sparseStorage[part] =
457  partitions[part].first, partitions[part].first,
458  partitions[part].second, partMat, matStorage );
459  }
460 
461  // Create block diagonal matrix
463  AllocateSharedPtr(sparseStorage);
464 
465  break;
466  }
467  default:
468  ErrorUtil::NekError("Solver info property \
469  LocalMatrixStorageStrategy takes values \
470  Contiguous, Non-contiguous and Sparse");
471  }
472  }
473 
474  /**
475  *
476  */
478  const Array<OneD, NekDouble>& pInput,
479  Array<OneD, NekDouble>& pOutput)
480  {
481  int nLocal = m_locToGloMap->GetNumLocalBndCoeffs();
482  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
483  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
484  DoGlobalMatOp(m_linSysKey.GetMatrixType());
485 
486  if(doGlobalOp)
487  {
488  // Do matrix multiply globally
489  Array<OneD, NekDouble> in = pInput + nDir;
490  Array<OneD, NekDouble> out = pOutput + nDir;
491 
492  m_sparseSchurCompl->Multiply(in,out);
493  m_locToGloMap->UniversalAssembleBnd(pOutput, nDir);
494  }
495  else if (m_sparseSchurCompl)
496  {
497  // Do matrix multiply locally using block-diagonal sparse matrix
498  Array<OneD, NekDouble> tmp = m_wsp + nLocal;
499 
500  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
501  m_sparseSchurCompl->Multiply(m_wsp,tmp);
502  m_locToGloMap->AssembleBnd(tmp, pOutput);
503  }
504  else
505  {
506  // Do matrix multiply locally, using direct BLAS calls
507  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
508  int i, cnt;
509  Array<OneD, NekDouble> tmpout = m_wsp + nLocal;
510  for (i = cnt = 0; i < m_denseBlocks.size(); cnt += m_rows[i], ++i)
511  {
512  const int rows = m_rows[i];
513  Blas::Dgemv('N', rows, rows,
514  m_scale[i], m_denseBlocks[i], rows,
515  m_wsp.get()+cnt, 1,
516  0.0, tmpout.get()+cnt, 1);
517  }
518  m_locToGloMap->AssembleBnd(tmpout, pOutput);
519  }
520  }
521 
523  {
524  m_map = m_locToGloMap->GetGlobalToUniversalBndMapUnique();
525  }
526 
528  int scLevel,
529  NekVector<NekDouble> &F_GlobBnd)
530  {
531  if (scLevel == 0)
532  {
533  // When matrices are supplied to the constructor at the top
534  // level, the preconditioner is never set up.
535  if (!m_precon)
536  {
538  m_precon->BuildPreconditioner();
539  }
540 
541  Set_Rhs_Magnitude(F_GlobBnd);
542 
543  return m_S1Blk;
544  }
545  else
546  {
547  // for multilevel iterative solver always use rhs
548  // vector value with no weighting
550 
551  return m_schurCompl;
552  }
553  }
554 
556  Array<OneD, NekDouble>& pInOut,
557  int offset)
558  {
559  m_precon->DoTransformToLowEnergy(pInOut, offset);
560  }
561 
563  Array<OneD, NekDouble>& pInOut)
564  {
565  m_precon->DoTransformFromLowEnergy(pInOut);
566  }
567 
569  const GlobalLinSysKey &mkey,
570  const boost::weak_ptr<ExpList> &pExpList,
571  const DNekScalBlkMatSharedPtr pSchurCompl,
572  const DNekScalBlkMatSharedPtr pBinvD,
573  const DNekScalBlkMatSharedPtr pC,
574  const DNekScalBlkMatSharedPtr pInvD,
575  const boost::shared_ptr<AssemblyMap> &l2gMap)
576  {
578  GlobalLinSysIterativeStaticCond>::AllocateSharedPtr(
579  mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
580  m_precon);
581  sys->Initialise(l2gMap);
582  return sys;
583  }
584  }
585 }
Array< OneD, NekDouble > m_wsp
Workspace array for matrix multiplication.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
virtual DNekScalBlkMatSharedPtr v_PreSolve(int scLevel, NekVector< NekDouble > &F_GlobBnd)
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
std::map< CoordType, NekDouble > COOMatType
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
void Set_Rhs_Magnitude(const NekVector< NekDouble > &pIn)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
void v_AssembleSchurComplement(const boost::shared_ptr< AssemblyMap > locToGloMap)
Assemble the Schur complement matrix.
std::vector< double > m_storage
Dense storage for block Schur complement matrix.
DNekScalBlkMatSharedPtr m_C
Block matrix.
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Perform a Shur-complement matrix multiply operation.
boost::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:62
void Initialise(const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Definition: GlobalLinSys.h:214
DNekScalBlkMatSharedPtr m_invD
Block matrix.
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
Array< OneD, NekDouble > BCOEntryType
GlobalLinSysIterativeStaticCond(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
static const NekDouble kNekZeroTol
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
boost::shared_ptr< GlobalLinSysIterativeStaticCond > GlobalLinSysIterativeStaticCondSharedPtr
virtual GlobalLinSysStaticCondSharedPtr v_Recurse(const GlobalLinSysKey &mkey, const boost::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const boost::shared_ptr< AssemblyMap > &locToGloMap)
Array< OneD, NekDouble > m_scale
Scaling factors for local matrices.
virtual void v_BasisInvTransform(Array< OneD, NekDouble > &pInOut)
std::map< CoordType, BCOEntryType > BCOMatType
Array< OneD, int > m_map
Global to universal unique map.
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
Array< OneD, unsigned int > m_rows
Ranks of local matrices.
static GlobalLinSysSharedPtr create(const GlobalLinSysKey &pLinSysKey, const boost::weak_ptr< ExpList > &pExpList, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
double NekDouble
std::vector< const double * > m_denseBlocks
Vector of pointers to local matrix data.
void PrepareLocalSchurComplement()
Prepares local representation of Schur complement stored as a sparse block-diagonal matrix...
void SetupTopLevel(const boost::shared_ptr< AssemblyMap > &locToGloMap)
Set up the storage for the Schur complement or the top level of the multi-level Schur complement...
static const NekDouble kNekUnsetDouble
Describe a linear system.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
virtual void v_BasisTransform(Array< OneD, NekDouble > &pInOut, int offset)
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:127
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
A global linear system.
Definition: GlobalLinSys.h:74
DNekSmvBsrDiagBlkMatSharedPtr m_sparseSchurCompl
Sparse representation of Schur complement matrix at this level.
virtual void v_DropStaticCondBlock(unsigned int n)
Releases the static condensation block matrices from NekManager of n-th expansion using the matrix ke...
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
void convertCooToBco(const unsigned int blkRows, const unsigned int blkColumns, const unsigned int blkDim, const COOMatType &cooMat, BCOMatType &bcoMat)
Definition: SparseUtils.cpp:48
boost::shared_ptr< AssemblyMap > m_locToGloMap
Local to global map.
GlobalLinSysFactory & GetGlobalLinSysFactory()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
virtual DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n)
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...
DNekScalBlkMatSharedPtr m_BinvD
Block matrix.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map...
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215
boost::shared_ptr< GlobalLinSysStaticCond > GlobalLinSysStaticCondSharedPtr
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129