Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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 using namespace std;
45 
46 namespace Nektar
47 {
48  namespace MultiRegions
49  {
50  /**
51  * @class GlobalLinSysIterativeStaticCond
52  *
53  * Solves a linear system iteratively using single- or multi-level
54  * static condensation.
55  */
56 
57  /**
58  * Registers the class with the Factory.
59  */
60  string GlobalLinSysIterativeStaticCond::className
62  "IterativeStaticCond",
63  GlobalLinSysIterativeStaticCond::create,
64  "Iterative static condensation.");
65 
66  string GlobalLinSysIterativeStaticCond::className2
68  "IterativeMultiLevelStaticCond",
69  GlobalLinSysIterativeStaticCond::create,
70  "Iterative multi-level static condensation.");
71 
72 
73  std::string GlobalLinSysIterativeStaticCond::storagedef =
74  LibUtilities::SessionReader::RegisterDefaultSolverInfo(
75  "LocalMatrixStorageStrategy",
76  "Sparse");
77  std::string GlobalLinSysIterativeStaticCond::storagelookupIds[3] = {
78  LibUtilities::SessionReader::RegisterEnumValue(
79  "LocalMatrixStorageStrategy",
80  "Contiguous",
82  LibUtilities::SessionReader::RegisterEnumValue(
83  "LocalMatrixStorageStrategy",
84  "Non-contiguous",
86  LibUtilities::SessionReader::RegisterEnumValue(
87  "LocalMatrixStorageStrategy",
88  "Sparse",
90  };
91 
92  /**
93  * For a matrix system of the form @f[
94  * \left[ \begin{array}{cc}
95  * \boldsymbol{A} & \boldsymbol{B}\\
96  * \boldsymbol{C} & \boldsymbol{D}
97  * \end{array} \right]
98  * \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2}
99  * \end{array}\right]
100  * = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2}
101  * \end{array}\right],
102  * @f]
103  * where @f$\boldsymbol{D}@f$ and
104  * @f$(\boldsymbol{A-BD^{-1}C})@f$ are invertible, store and assemble
105  * a static condensation system, according to a given local to global
106  * mapping. #m_linSys is constructed by AssembleSchurComplement().
107  * @param mKey Associated matrix key.
108  * @param pLocMatSys LocalMatrixSystem
109  * @param locToGloMap Local to global mapping.
110  */
111  GlobalLinSysIterativeStaticCond::GlobalLinSysIterativeStaticCond(
112  const GlobalLinSysKey &pKey,
113  const boost::weak_ptr<ExpList> &pExpList,
114  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
115  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
116  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
117  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
118  {
121  "This constructor is only valid when using static "
122  "condensation");
124  == pLocToGloMap->GetGlobalSysSolnType(),
125  "The local to global map is not set up for the requested "
126  "solution type");
127  }
128 
129 
130  /**
131  *
132  */
134  const GlobalLinSysKey &pKey,
135  const boost::weak_ptr<ExpList> &pExpList,
136  const DNekScalBlkMatSharedPtr pSchurCompl,
137  const DNekScalBlkMatSharedPtr pBinvD,
138  const DNekScalBlkMatSharedPtr pC,
139  const DNekScalBlkMatSharedPtr pInvD,
140  const boost::shared_ptr<AssemblyMap> &pLocToGloMap,
141  const PreconditionerSharedPtr pPrecon)
142  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
143  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
144  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
145  {
146  m_schurCompl = pSchurCompl;
147  m_S1Blk = pSchurCompl;
148  m_BinvD = pBinvD;
149  m_C = pC;
150  m_invD = pInvD;
151  m_precon = pPrecon;
152  }
153 
154 
156  {
158 
159  // Allocate memory for top-level structure
161 
162  // Setup Block Matrix systems
163  int n, n_exp = m_expList.lock()->GetNumElmts();
164 
165  MatrixStorage blkmatStorage = eDIAGONAL;
166  const Array<OneD,const unsigned int>& nbdry_size
167  = m_locToGloMap->GetNumLocalBndCoeffsPerPatch();
168 
170  ::AllocateSharedPtr(nbdry_size, nbdry_size, blkmatStorage);
171 
172  // Preserve original matrix in m_S1Blk
173  for (n = 0; n < n_exp; ++n)
174  {
175  DNekScalMatSharedPtr mat = m_schurCompl->GetBlock(n, n);
176  m_S1Blk->SetBlock(n, n, mat);
177  }
178 
179  // Build preconditioner
180  m_precon->BuildPreconditioner();
181 
182  // Do transform of Schur complement matrix
183  for (n = 0; n < n_exp; ++n)
184  {
185  if (m_linSysKey.GetMatrixType() !=
187  {
188  DNekScalMatSharedPtr mat = m_S1Blk->GetBlock(n, n);
189  DNekScalMatSharedPtr t = m_precon->TransformedSchurCompl(n, mat);
190  m_schurCompl->SetBlock(n, n, t);
191  }
192  }
193 
194  // Construct this level
196  }
197 
198  /**
199  *
200  */
202  {
203 
204  }
205 
207  v_GetStaticCondBlock(unsigned int n)
208  {
209  DNekScalBlkMatSharedPtr schurComplBlock;
210  int scLevel = m_locToGloMap->GetStaticCondLevel();
211  DNekScalBlkMatSharedPtr sc = scLevel == 0 ? m_S1Blk : m_schurCompl;
212  DNekScalMatSharedPtr localMat = sc->GetBlock(n,n);
213  unsigned int nbdry = localMat->GetRows();
214  unsigned int nblks = 1;
215  unsigned int esize[1] = {nbdry};
216 
217  schurComplBlock = MemoryManager<DNekScalBlkMat>
218  ::AllocateSharedPtr(nblks, nblks, esize, esize);
219  schurComplBlock->SetBlock(0, 0, localMat);
220 
221  return schurComplBlock;
222  }
223 
224  /**
225  * Assemble the schur complement matrix from the block matrices stored
226  * in #m_blkMatrices and the given local to global mapping information.
227  * @param locToGloMap Local to global mapping information.
228  */
230  const AssemblyMapSharedPtr pLocToGloMap)
231  {
232  int i,j,n,cnt,gid1,gid2;
233  NekDouble sign1,sign2;
234 
235  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
236  DoGlobalMatOp(m_linSysKey.GetMatrixType());
237 
238  // Set up unique map
239  v_UniqueMap();
240 
241  // Build precon again if we in multi-level static condensation (a
242  // bit of a hack)
245  {
247  m_precon->BuildPreconditioner();
248  }
249 
250  if (!doGlobalOp)
251  {
253  return;
254  }
255 
256  int nBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
257  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
258  unsigned int rows = nBndDofs - NumDirBCs;
259  unsigned int cols = nBndDofs - NumDirBCs;
260 
261  // COO sparse storage to assist in assembly
262  COOMatType gmat_coo;
263 
264  // Get the matrix storage structure
265  // (whether to store only one triangular part, if symmetric)
266  MatrixStorage matStorage = eFULL;
267 
268  // assemble globally
269  DNekScalMatSharedPtr loc_mat;
270  int loc_lda;
271  for(n = cnt = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
272  {
273  loc_mat = m_schurCompl->GetBlock(n,n);
274  loc_lda = loc_mat->GetRows();
275 
276  // Set up Matrix;
277  for(i = 0; i < loc_lda; ++i)
278  {
279  gid1 = pLocToGloMap->GetLocalToGlobalBndMap (cnt + i)
280  - NumDirBCs;
281  sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
282 
283  if(gid1 >= 0)
284  {
285  for(j = 0; j < loc_lda; ++j)
286  {
287  gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt+j)
288  - NumDirBCs;
289  sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt+j);
290 
291  if (gid2 >= 0)
292  {
293  gmat_coo[std::make_pair(gid1,gid2)] +=
294  sign1*sign2*(*loc_mat)(i,j);
295  }
296  }
297  }
298  }
299  cnt += loc_lda;
300  }
301 
303  sparseStorage (1);
304 
305  BCOMatType partMat;
306  convertCooToBco(rows, cols, 1, gmat_coo, partMat);
307 
308  sparseStorage[0] =
310  AllocateSharedPtr(rows, cols, 1, partMat, matStorage );
311 
312  // Create block diagonal matrix
314  AllocateSharedPtr(sparseStorage);
315  }
316 
317 
318  /**
319  * Populates sparse block-diagonal schur complement matrix from
320  * the block matrices stored in #m_blkMatrices.
321  */
323  {
324  LocalMatrixStorageStrategy storageStrategy =
325  m_expList.lock()->GetSession()->
326  GetSolverInfoAsEnum<LocalMatrixStorageStrategy>(
327  "LocalMatrixStorageStrategy");
328 
329  switch(storageStrategy)
330  {
333  {
334  size_t storageSize = 0;
335  int nBlk = m_schurCompl->GetNumberOfBlockRows();
336 
337  m_scale = Array<OneD, NekDouble> (nBlk, 1.0);
338  m_rows = Array<OneD, unsigned int> (nBlk, 0U);
339 
340  // Determine storage requirements for dense blocks.
341  for (int i = 0; i < nBlk; ++i)
342  {
343  m_rows[i] = m_schurCompl->GetBlock(i,i)->GetRows();
344  m_scale[i] = m_schurCompl->GetBlock(i,i)->Scale();
345  storageSize += m_rows[i] * m_rows[i];
346  }
347 
348  // Assemble dense storage blocks.
349  DNekScalMatSharedPtr loc_mat;
350  m_denseBlocks.resize(nBlk);
351  double *ptr = 0;
352 
353  if (MultiRegions::eContiguous == storageStrategy)
354  {
355  m_storage.resize (storageSize);
356  ptr = &m_storage[0];
357  }
358 
359  for (unsigned int n = 0; n < nBlk; ++n)
360  {
361  loc_mat = m_schurCompl->GetBlock(n,n);
362 
363  if (MultiRegions::eContiguous == storageStrategy)
364  {
365  int loc_lda = loc_mat->GetRows();
366  int blockSize = loc_lda * loc_lda;
367  m_denseBlocks[n] = ptr;
368  for(int i = 0; i < loc_lda; ++i)
369  {
370  for(int j = 0; j < loc_lda; ++j)
371  {
372  ptr[j*loc_lda+i] = (*loc_mat)(i,j);
373  }
374  }
375  ptr += blockSize;
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  }
452 
453  sparseStorage[part] =
456  partitions[part].first, partitions[part].first,
457  partitions[part].second, partMat, matStorage );
458  }
459 
460  // Create block diagonal matrix
462  AllocateSharedPtr(sparseStorage);
463 
464  break;
465  }
466  default:
467  ErrorUtil::NekError("Solver info property \
468  LocalMatrixStorageStrategy takes values \
469  Contiguous, Non-contiguous and Sparse");
470  }
471  }
472 
473  /**
474  *
475  */
477  const Array<OneD, NekDouble>& pInput,
478  Array<OneD, NekDouble>& pOutput)
479  {
480  int nLocal = m_locToGloMap->GetNumLocalBndCoeffs();
481  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
482  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
483  DoGlobalMatOp(m_linSysKey.GetMatrixType());
484 
485  if(doGlobalOp)
486  {
487  // Do matrix multiply globally
488  Array<OneD, NekDouble> in = pInput + nDir;
489  Array<OneD, NekDouble> out = pOutput + nDir;
490 
491  m_sparseSchurCompl->Multiply(in,out);
492  m_locToGloMap->UniversalAssembleBnd(pOutput, nDir);
493  }
494  else if (m_sparseSchurCompl)
495  {
496  // Do matrix multiply locally using block-diagonal sparse matrix
497  Array<OneD, NekDouble> tmp = m_wsp + nLocal;
498 
499  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
500  m_sparseSchurCompl->Multiply(m_wsp,tmp);
501  m_locToGloMap->AssembleBnd(tmp, pOutput);
502  }
503  else
504  {
505  // Do matrix multiply locally, using direct BLAS calls
506  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
507  int i, cnt;
508  Array<OneD, NekDouble> tmpout = m_wsp + nLocal;
509  for (i = cnt = 0; i < m_denseBlocks.size(); cnt += m_rows[i], ++i)
510  {
511  const int rows = m_rows[i];
512  Blas::Dgemv('N', rows, rows,
513  m_scale[i], m_denseBlocks[i], rows,
514  m_wsp.get()+cnt, 1,
515  0.0, tmpout.get()+cnt, 1);
516  }
517  m_locToGloMap->AssembleBnd(tmpout, pOutput);
518  }
519  }
520 
522  {
523  m_map = m_locToGloMap->GetGlobalToUniversalBndMapUnique();
524  }
525 
527  int scLevel,
528  NekVector<NekDouble> &F_GlobBnd)
529  {
530  if (scLevel == 0)
531  {
532  // When matrices are supplied to the constructor at the top
533  // level, the preconditioner is never set up.
534  if (!m_precon)
535  {
537  m_precon->BuildPreconditioner();
538  }
539 
540  Set_Rhs_Magnitude(F_GlobBnd);
541 
542  return m_S1Blk;
543  }
544  else
545  {
546  // for multilevel iterative solver always use rhs
547  // vector value with no weighting
549 
550  return m_schurCompl;
551  }
552  }
553 
555  Array<OneD, NekDouble>& pInOut,
556  int offset)
557  {
558  m_precon->DoTransformToLowEnergy(pInOut, offset);
559  }
560 
562  Array<OneD, NekDouble>& pInOut)
563  {
564  m_precon->DoTransformFromLowEnergy(pInOut);
565  }
566 
568  const GlobalLinSysKey &mkey,
569  const boost::weak_ptr<ExpList> &pExpList,
570  const DNekScalBlkMatSharedPtr pSchurCompl,
571  const DNekScalBlkMatSharedPtr pBinvD,
572  const DNekScalBlkMatSharedPtr pC,
573  const DNekScalBlkMatSharedPtr pInvD,
574  const boost::shared_ptr<AssemblyMap> &l2gMap)
575  {
577  GlobalLinSysIterativeStaticCond>::AllocateSharedPtr(
578  mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
579  m_precon);
580  sys->Initialise(l2gMap);
581  return sys;
582  }
583  }
584 }
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 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.
STL namespace.
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:216
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:213
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.
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...
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:228
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