Nektar++
GlobalLinSysStaticCond.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GlobalLinSysStaticCond.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: Implementation to linear solver using single-
32// or multi-level static condensation
33//
34///////////////////////////////////////////////////////////////////////////////
35
41
42using namespace std;
43
44namespace Nektar
45{
46namespace MultiRegions
47{
48/**
49 * @class GlobalLinSysStaticCond
50 *
51 * Solves a linear system using single- or multi-level static
52 * condensation.
53 */
54
55/**
56 * For a matrix system of the form @f[
57 * \left[ \begin{array}{cc}
58 * \boldsymbol{A} & \boldsymbol{B}\\
59 * \boldsymbol{C} & \boldsymbol{D}
60 * \end{array} \right]
61 * \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2}
62 * \end{array}\right]
63 * = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2}
64 * \end{array}\right],
65 * @f]
66 * where @f$\boldsymbol{D}@f$ and
67 * @f$(\boldsymbol{A-BD^{-1}C})@f$ are invertible, store and assemble
68 * a static condensation system, according to a given local to global
69 * mapping. #m_linSys is constructed by AssembleSchurComplement().
70 * @param mKey Associated matrix key.
71 * @param pLocMatSys LocalMatrixSystem
72 * @param locToGloMap Local to global mapping.
73 */
75 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExpList,
76 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
77 : GlobalLinSys(pKey, pExpList, pLocToGloMap), m_locToGloMap(pLocToGloMap)
78{
79}
80
82{
83 // Allocate memory for top-level structure
85
86 // Construct this level
88}
89
90/**
91 *
92 */
94{
95}
96
97/**
98 *
99 */
101 const Array<OneD, const NekDouble> &pLocInput,
102 Array<OneD, NekDouble> &pLocOutput,
103 const AssemblyMapSharedPtr &pLocToGloMap,
104 const Array<OneD, const NekDouble> &dirForcing)
105{
106 boost::ignore_unused(dirForcing);
107 ASSERTL1(dirForcing.size() == 0,
108 "GlobalLinSysStaticCond: Not setup for dirForcing");
109
110 bool atLastLevel = pLocToGloMap->AtLastLevel();
111 int scLevel = pLocToGloMap->GetStaticCondLevel();
112
113 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
114 int nLocBndDofs = pLocToGloMap->GetNumLocalBndCoeffs();
115 int nGlobBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
116 int nDirBndDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
117 int nIntDofs = nGlobDofs - nGlobBndDofs;
118
119 if ((nGlobDofs - nDirBndDofs) == 0)
120 {
121 return; // nothing to solve;
122 }
123
124 Array<OneD, NekDouble> F_bnd, F_bnd1, F_int, V_bnd;
125 Array<OneD, NekDouble> tmp, unused;
126
127 F_bnd = m_wsp;
128 // unsued = Temporary storage required for DoMatrixMultiply, not used here
129 unused = m_wsp + 1 * nLocBndDofs;
130 F_bnd1 = m_wsp + 2 * nLocBndDofs;
131 V_bnd = m_wsp + 3 * nLocBndDofs;
132 F_int = m_wsp + 4 * nLocBndDofs;
133
134 pLocToGloMap->LocalToLocalBnd(pLocOutput, V_bnd);
135
136 NekVector<NekDouble> F_Int(nIntDofs, F_int, eWrapper);
137 NekVector<NekDouble> V_Bnd(nLocBndDofs, V_bnd, eWrapper);
138
139 if (nIntDofs)
140 {
141 m_locToGloMap.lock()->LocalToLocalInt(pLocInput, F_int);
142 }
143
144 // Boundary system solution
145 if (nGlobBndDofs - nDirBndDofs)
146 {
147 pLocToGloMap->LocalToLocalBnd(pLocInput, F_bnd);
148
149 // set up normalisation factor for right hand side on first SC level
150 v_PreSolve(scLevel, F_bnd);
151
152 // Gather boundary expansison into locbnd
153 NekVector<NekDouble> F_Bnd(nLocBndDofs, F_bnd1, eWrapper);
154
155 // construct boundary forcing
156 if (nIntDofs)
157 {
158 DNekScalBlkMat &BinvD = *m_BinvD;
159
160 F_Bnd = BinvD * F_Int;
161
162 Vmath::Vsub(nLocBndDofs, F_bnd, 1, F_bnd1, 1, F_bnd, 1);
163 }
164
165 if (atLastLevel)
166 {
167 // Transform to new basis if required
168 v_BasisFwdTransform(F_bnd);
169
170 DNekScalBlkMat &SchurCompl = *m_schurCompl;
171
172 v_CoeffsFwdTransform(V_bnd, V_bnd);
173
174 // subtract dirichlet boundary forcing
175 F_Bnd = SchurCompl * V_Bnd;
176
177 Vmath::Vsub(nLocBndDofs, F_bnd, 1, F_bnd1, 1, F_bnd, 1);
178
180
181 // Solve for difference from initial solution given inout;
182 SolveLinearSystem(nGlobBndDofs, F_bnd, F_bnd1, pLocToGloMap,
183 nDirBndDofs);
184
185 // Add back initial conditions onto difference
186 Vmath::Vadd(nLocBndDofs, V_bnd, 1, F_bnd1, 1, V_bnd, 1);
187
188 // Transform back to original basis
190
191 // put final bnd solution back in output array
192 m_locToGloMap.lock()->LocalBndToLocal(V_bnd, pLocOutput);
193 }
194 else // Process next level of recursion for multi level SC
195 {
196 AssemblyMapSharedPtr nextLevLocToGloMap =
197 pLocToGloMap->GetNextLevelLocalToGlobalMap();
198
199 // partially assemble F for next level and
200 // reshuffle V_bnd
201 nextLevLocToGloMap->PatchAssemble(F_bnd, F_bnd);
202 nextLevLocToGloMap->PatchLocalToGlobal(V_bnd, V_bnd);
203
204 m_recursiveSchurCompl->Solve(F_bnd, V_bnd, nextLevLocToGloMap);
205
206 // unpack V_bnd
207 nextLevLocToGloMap->PatchGlobalToLocal(V_bnd, V_bnd);
208
209 // place V_bnd in output array
210 m_locToGloMap.lock()->LocalBndToLocal(V_bnd, pLocOutput);
211 }
212 }
213
214 // solve interior system
215 if (nIntDofs)
216 {
217 Array<OneD, NekDouble> V_int(nIntDofs);
218 NekVector<NekDouble> V_Int(nIntDofs, V_int, eWrapper);
219
220 // get array of local solutions
221 DNekScalBlkMat &invD = *m_invD;
222 DNekScalBlkMat &C = *m_C;
223
224 F_Int = F_Int - C * V_Bnd;
225
226 Multiply(V_Int, invD, F_Int);
227
228 m_locToGloMap.lock()->LocalIntToLocal(V_int, pLocOutput);
229 }
230}
231
232/**
233 * If at the last level of recursion (or the only level in the case of
234 * single-level static condensation), assemble the Schur complement.
235 * For other levels, in the case of multi-level static condensation,
236 * the next level of the condensed system is computed.
237 * @param pLocToGloMap Local to global mapping.
238 */
240 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
241{
242 int nLocalBnd = m_locToGloMap.lock()->GetNumLocalBndCoeffs();
243 int nIntDofs = m_locToGloMap.lock()->GetNumLocalCoeffs() - nLocalBnd;
244 m_wsp = Array<OneD, NekDouble>(4 * nLocalBnd + nIntDofs, 0.0);
245
246 if (pLocToGloMap->AtLastLevel())
247 {
248 v_AssembleSchurComplement(pLocToGloMap);
249 }
250 else
251 {
253 pLocToGloMap->GetNextLevelLocalToGlobalMap());
254 }
255}
256
258{
259 return m_schurCompl->GetNumberOfBlockRows();
260}
261
262/**
263 * For the first level in multi-level static condensation, or the only
264 * level in the case of single-level static condensation, allocate the
265 * condensed matrices and populate them with the local matrices
266 * retrieved from the expansion list.
267 * @param
268 */
270 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
271{
272 int n;
273 int n_exp = m_expList.lock()->GetNumElmts();
274
275 const Array<OneD, const unsigned int> &nbdry_size =
276 pLocToGloMap->GetNumLocalBndCoeffsPerPatch();
277 const Array<OneD, const unsigned int> &nint_size =
278 pLocToGloMap->GetNumLocalIntCoeffsPerPatch();
279
280 // Setup Block Matrix systems
281 MatrixStorage blkmatStorage = eDIAGONAL;
283 nbdry_size, nbdry_size, blkmatStorage);
285 nbdry_size, nint_size, blkmatStorage);
287 nint_size, nbdry_size, blkmatStorage);
289 nint_size, nint_size, blkmatStorage);
290
291 for (n = 0; n < n_exp; ++n)
292 {
294 {
296 m_schurCompl->SetBlock(n, n, loc_mat);
297 }
298 else
299 {
300 DNekScalBlkMatSharedPtr loc_schur =
303 m_schurCompl->SetBlock(n, n, t = loc_schur->GetBlock(0, 0));
304 m_BinvD->SetBlock(n, n, t = loc_schur->GetBlock(0, 1));
305 m_C->SetBlock(n, n, t = loc_schur->GetBlock(1, 0));
306 m_invD->SetBlock(n, n, t = loc_schur->GetBlock(1, 1));
307 }
308 }
309}
310
311/**
312 *
313 */
315 const AssemblyMapSharedPtr &pLocToGloMap)
316{
317 int i, j, n, cnt;
318 DNekScalBlkMatSharedPtr blkMatrices[4];
319
320 // Create temporary matrices within an inner-local scope to ensure
321 // any references to the intermediate storage is lost before
322 // the recursive step, rather than at the end of the routine.
323 // This allows the schur complement matrix from this level to be
324 // disposed of in the next level after use without being retained
325 // due to lingering shared pointers.
326 {
327
328 const Array<OneD, const unsigned int> &nBndDofsPerPatch =
329 pLocToGloMap->GetNumLocalBndCoeffsPerPatch();
330 const Array<OneD, const unsigned int> &nIntDofsPerPatch =
331 pLocToGloMap->GetNumLocalIntCoeffsPerPatch();
332
333 // STEP 1:
334 // Based upon the schur complement of the the current level we
335 // will substructure this matrix in the form
336 // -- --
337 // | A B |
338 // | C D |
339 // -- --
340 // All matrices A,B,C and D are (diagonal) blockmatrices.
341 // However, as a start we will use an array of DNekMatrices as
342 // it is too hard to change the individual entries of a
343 // DNekScalBlkMatSharedPtr.
344
345 // In addition, we will also try to ensure that the memory of
346 // the blockmatrices will be contiguous. This will probably
347 // enhance the efficiency
348 // - Calculate the total number of entries in the blockmatrices
349 int nPatches = pLocToGloMap->GetNumPatches();
350 int nEntriesA = 0;
351 int nEntriesB = 0;
352 int nEntriesC = 0;
353 int nEntriesD = 0;
354
355 for (i = 0; i < nPatches; i++)
356 {
357 nEntriesA += nBndDofsPerPatch[i] * nBndDofsPerPatch[i];
358 nEntriesB += nBndDofsPerPatch[i] * nIntDofsPerPatch[i];
359 nEntriesC += nIntDofsPerPatch[i] * nBndDofsPerPatch[i];
360 nEntriesD += nIntDofsPerPatch[i] * nIntDofsPerPatch[i];
361 }
362
363 // Now create the DNekMatrices and link them to the memory
364 // allocated above
365 Array<OneD, DNekMatSharedPtr> substructuredMat[4] = {
366 Array<OneD, DNekMatSharedPtr>(nPatches), // Matrix A
367 Array<OneD, DNekMatSharedPtr>(nPatches), // Matrix B
368 Array<OneD, DNekMatSharedPtr>(nPatches), // Matrix C
369 Array<OneD, DNekMatSharedPtr>(nPatches)}; // Matrix D
370
371 // Initialise storage for the matrices. We do this separately
372 // for each matrix so the matrices may be independently
373 // deallocated when no longer required.
374 Array<OneD, NekDouble> storageA(nEntriesA, 0.0);
375 Array<OneD, NekDouble> storageB(nEntriesB, 0.0);
376 Array<OneD, NekDouble> storageC(nEntriesC, 0.0);
377 Array<OneD, NekDouble> storageD(nEntriesD, 0.0);
378
379 Array<OneD, NekDouble> tmparray;
380 PointerWrapper wType = eWrapper;
381 int cntA = 0;
382 int cntB = 0;
383 int cntC = 0;
384 int cntD = 0;
385
386 // Use symmetric storage for invD if possible
387 MatrixStorage storageTypeD = eFULL;
391 {
392 storageTypeD = eSYMMETRIC;
393 }
394
395 for (i = 0; i < nPatches; i++)
396 {
397 // Matrix A
398 tmparray = storageA + cntA;
399 substructuredMat[0][i] = MemoryManager<DNekMat>::AllocateSharedPtr(
400 nBndDofsPerPatch[i], nBndDofsPerPatch[i], tmparray, wType);
401 // Matrix B
402 tmparray = storageB + cntB;
403 substructuredMat[1][i] = MemoryManager<DNekMat>::AllocateSharedPtr(
404 nBndDofsPerPatch[i], nIntDofsPerPatch[i], tmparray, wType);
405 // Matrix C
406 tmparray = storageC + cntC;
407 substructuredMat[2][i] = MemoryManager<DNekMat>::AllocateSharedPtr(
408 nIntDofsPerPatch[i], nBndDofsPerPatch[i], tmparray, wType);
409 // Matrix D
410 tmparray = storageD + cntD;
411 substructuredMat[3][i] = MemoryManager<DNekMat>::AllocateSharedPtr(
412 nIntDofsPerPatch[i], nIntDofsPerPatch[i], tmparray, wType,
413 storageTypeD);
414
415 cntA += nBndDofsPerPatch[i] * nBndDofsPerPatch[i];
416 cntB += nBndDofsPerPatch[i] * nIntDofsPerPatch[i];
417 cntC += nIntDofsPerPatch[i] * nBndDofsPerPatch[i];
418 cntD += nIntDofsPerPatch[i] * nIntDofsPerPatch[i];
419 }
420
421 // Then, project SchurComplement onto
422 // the substructured matrices of the next level
424 DNekScalMatSharedPtr schurComplSubMat;
425 int schurComplSubMatnRows;
426 Array<OneD, const int> patchId, dofId;
429
430 for (n = cnt = 0; n < SchurCompl->GetNumberOfBlockRows(); ++n)
431 {
432 schurComplSubMat = SchurCompl->GetBlock(n, n);
433 schurComplSubMatnRows = schurComplSubMat->GetRows();
434
435 patchId =
436 pLocToGloMap->GetPatchMapFromPrevLevel()->GetPatchId() + cnt;
437 dofId = pLocToGloMap->GetPatchMapFromPrevLevel()->GetDofId() + cnt;
438 isBndDof =
439 pLocToGloMap->GetPatchMapFromPrevLevel()->IsBndDof() + cnt;
440 sign = pLocToGloMap->GetPatchMapFromPrevLevel()->GetSign() + cnt;
441
442 // Set up Matrix;
443 for (i = 0; i < schurComplSubMatnRows; ++i)
444 {
445 int pId = patchId[i];
446 Array<OneD, NekDouble> subMat0 =
447 substructuredMat[0][pId]->GetPtr();
448 Array<OneD, NekDouble> subMat1 =
449 substructuredMat[1][patchId[i]]->GetPtr();
450 Array<OneD, NekDouble> subMat2 =
451 substructuredMat[2][patchId[i]]->GetPtr();
452 DNekMatSharedPtr subMat3 = substructuredMat[3][patchId[i]];
453 int subMat0rows = substructuredMat[0][pId]->GetRows();
454 int subMat1rows = substructuredMat[1][pId]->GetRows();
455 int subMat2rows = substructuredMat[2][pId]->GetRows();
456
457 if (isBndDof[i])
458 {
459 for (j = 0; j < schurComplSubMatnRows; ++j)
460 {
461 ASSERTL0(patchId[i] == patchId[j],
462 "These values should be equal");
463
464 if (isBndDof[j])
465 {
466 subMat0[dofId[i] + dofId[j] * subMat0rows] +=
467 sign[i] * sign[j] * (*schurComplSubMat)(i, j);
468 }
469 else
470 {
471 subMat1[dofId[i] + dofId[j] * subMat1rows] +=
472 sign[i] * sign[j] * (*schurComplSubMat)(i, j);
473 }
474 }
475 }
476 else
477 {
478 for (j = 0; j < schurComplSubMatnRows; ++j)
479 {
480 ASSERTL0(patchId[i] == patchId[j],
481 "These values should be equal");
482
483 if (isBndDof[j])
484 {
485 subMat2[dofId[i] + dofId[j] * subMat2rows] +=
486 sign[i] * sign[j] * (*schurComplSubMat)(i, j);
487 }
488 else
489 {
490 if (storageTypeD == eSYMMETRIC)
491 {
492 if (dofId[i] <= dofId[j])
493 {
494 (*subMat3)(dofId[i], dofId[j]) +=
495 sign[i] * sign[j] *
496 (*schurComplSubMat)(i, j);
497 }
498 }
499 else
500 {
501 (*subMat3)(dofId[i], dofId[j]) +=
502 sign[i] * sign[j] *
503 (*schurComplSubMat)(i, j);
504 }
505 }
506 }
507 }
508 }
509 cnt += schurComplSubMatnRows;
510 }
511
512 // STEP 2: condense the system
513 // This can be done elementally (i.e. patch per patch)
514 for (i = 0; i < nPatches; i++)
515 {
516 if (nIntDofsPerPatch[i])
517 {
518 Array<OneD, NekDouble> subMat0 =
519 substructuredMat[0][i]->GetPtr();
520 Array<OneD, NekDouble> subMat1 =
521 substructuredMat[1][i]->GetPtr();
522 Array<OneD, NekDouble> subMat2 =
523 substructuredMat[2][i]->GetPtr();
524 int subMat0rows = substructuredMat[0][i]->GetRows();
525 int subMat1rows = substructuredMat[1][i]->GetRows();
526 int subMat2rows = substructuredMat[2][i]->GetRows();
527 int subMat2cols = substructuredMat[2][i]->GetColumns();
528
529 // 1. D -> InvD
530 substructuredMat[3][i]->Invert();
531 // 2. B -> BInvD
532 (*substructuredMat[1][i]) =
533 (*substructuredMat[1][i]) * (*substructuredMat[3][i]);
534 // 3. A -> A - BInvD*C (= schurcomplement)
535 // (*substructuredMat[0][i]) =(*substructuredMat[0][i])-
536 // (*substructuredMat[1][i])*(*substructuredMat[2][i]);
537 // Note: faster to use blas directly
538 if (subMat1rows && subMat2cols)
539 {
540 Blas::Dgemm('N', 'N', subMat1rows, subMat2cols, subMat2rows,
541 -1.0, &subMat1[0], subMat1rows, &subMat2[0],
542 subMat2rows, 1.0, &subMat0[0], subMat0rows);
543 }
544 }
545 }
546
547 // STEP 3: fill the block matrices. However, do note that we
548 // first have to convert them to a DNekScalMat in order to be
549 // compatible with the first level of static condensation
550 const Array<OneD, const unsigned int> &nbdry_size =
551 pLocToGloMap->GetNumLocalBndCoeffsPerPatch();
552 const Array<OneD, const unsigned int> &nint_size =
553 pLocToGloMap->GetNumLocalIntCoeffsPerPatch();
554 MatrixStorage blkmatStorage = eDIAGONAL;
555
557 nbdry_size, nbdry_size, blkmatStorage);
559 nbdry_size, nint_size, blkmatStorage);
561 nint_size, nbdry_size, blkmatStorage);
563 nint_size, nint_size, blkmatStorage);
564
565 DNekScalMatSharedPtr tmpscalmat;
566 for (i = 0; i < nPatches; i++)
567 {
568 for (j = 0; j < 4; j++)
569 {
571 1.0, substructuredMat[j][i]);
572 blkMatrices[j]->SetBlock(i, i, tmpscalmat);
573 }
574 }
575 }
576
577 // We've finished with the Schur complement matrix passed to this
578 // level, so return the memory to the system. The Schur complement
579 // matrix need only be retained at the last level. Save the other
580 // matrices at this level though.
581 m_schurCompl.reset();
582
584 v_Recurse(m_linSysKey, m_expList, blkMatrices[0], blkMatrices[1],
585 blkMatrices[2], blkMatrices[3], pLocToGloMap);
586}
587} // namespace MultiRegions
588} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A global linear system.
Definition: GlobalLinSys.h:72
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:124
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.
Definition: GlobalLinSys.h:192
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 DNekScalBlkMatSharedPtr v_GetStaticCondBlock(unsigned int n)
Retrieves a the static condensation block matrices from n-th expansion using the matrix key provided ...
void Initialise(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Definition: GlobalLinSys.h:205
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
virtual void v_CoeffsFwdTransform(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
std::weak_ptr< AssemblyMap > m_locToGloMap
Local to global map.
void SetupTopLevel(const std::shared_ptr< AssemblyMap > &locToGloMap)
Set up the storage for the Schur complement or the top level of the multi-level Schur complement.
Array< OneD, NekDouble > m_wsp
Workspace array for matrix multiplication.
virtual void v_AssembleSchurComplement(std::shared_ptr< AssemblyMap > pLoctoGloMap)
GlobalLinSysStaticCond(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
virtual void v_PreSolve(int scLevel, Array< OneD, NekDouble > &F_bnd)
virtual void v_BasisFwdTransform(Array< OneD, NekDouble > &pInOut)
GlobalLinSysStaticCondSharedPtr m_recursiveSchurCompl
Schur complement for Direct Static Condensation.
virtual void v_CoeffsBwdTransform(Array< OneD, NekDouble > &pInOut)
virtual void v_Initialise(const std::shared_ptr< AssemblyMap > &locToGloMap) override
Initialise this object.
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray) override
Solve the linear system for given input and output vectors using a specified local to global map.
DNekScalBlkMatSharedPtr m_BinvD
Block matrix.
virtual GlobalLinSysStaticCondSharedPtr v_Recurse(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const std::shared_ptr< AssemblyMap > &locToGloMap)=0
void ConstructNextLevelCondensedSystem(const std::shared_ptr< AssemblyMap > &locToGloMap)
DNekScalBlkMatSharedPtr m_C
Block matrix.
virtual int v_GetNumBlocks() override
Get the number of blocks in this system.
DNekScalBlkMatSharedPtr m_invD
Block matrix.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:385
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void Multiply(NekMatrix< ResultDataType, StandardMatrixTag > &result, const NekMatrix< LhsDataType, LhsMatrixType > &lhs, const ResultDataType &rhs)
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
PointerWrapper
Specifies if the pointer passed to a NekMatrix or NekVector is copied into an internal representation...
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414