Nektar++
PreconditionerBlock.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File PreconditionerBlock.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: Block Preconditioner definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
41 #include <LocalRegions/MatrixKey.h>
42 #include <LocalRegions/SegExp.h>
43 #include <math.h>
44 
45 using namespace std;
46 
47 namespace Nektar
48 {
49  namespace MultiRegions
50  {
51  /**
52  * Registers the class with the Factory.
53  */
54  string PreconditionerBlock::className
56  "Block",
57  PreconditionerBlock::create,
58  "Block Diagonal Preconditioning");
59  /**
60  * @class Block Preconditioner
61  *
62  * This class implements Block preconditioning for the conjugate
63  * gradient matrix solver.
64  */
65 
66  PreconditionerBlock::PreconditionerBlock(
67  const std::shared_ptr<GlobalLinSys> &plinsys,
68  const AssemblyMapSharedPtr &pLocToGloMap)
69  : Preconditioner(plinsys, pLocToGloMap)
70  {
71  }
72 
74  {
75  GlobalSysSolnType solvertype =
76  m_locToGloMap.lock()->GetGlobalSysSolnType();
78  solvertype == MultiRegions::ePETScStaticCond,
79  "Solver type not valid");
80  }
81 
83  {
84  GlobalLinSysKey key = m_linsys.lock()->GetKey();
85 
86  // Different setup for HDG and CG.
88  {
90  }
91  else
92  {
94  }
95  }
96 
97  /**
98  * \brief Construct a block preconditioner from \f$\mathbf{S}_{1}\f$ for
99  * the continuous Galerkin system.
100  *
101  * The preconditioner is defined as:
102  *
103  * \f[\mathbf{M}^{-1}=\left[\begin{array}{ccc}
104  * \mathrm{Diag}[(\mathbf{S_{1}})_{vv}] & & \\
105  * & (\mathbf{S}_{1})_{eb} & \\
106  * & & (\mathbf{S}_{1})_{fb} \end{array}\right] \f]
107  *
108  * where \f$\mathbf{S}_{1}\f$ is the local Schur complement matrix for
109  * each element and the subscript denotes the portion of the Schur
110  * complement associated with the vertex, edge and face blocks
111  * respectively.
112  */
114  {
115  ExpListSharedPtr expList = m_linsys.lock()->GetLocMat().lock();
117  DNekScalBlkMatSharedPtr loc_mat;
118  DNekScalMatSharedPtr bnd_mat;
119 
120  int i, j, k, n, cnt, gId;
121  int meshVertId, meshEdgeId, meshFaceId;
122 
123  auto asmMap = m_locToGloMap.lock();
124 
125  const int nExp = expList->GetExpSize();
126  const int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
127 
128  // Grab periodic geometry information.
129  PeriodicMap periodicVerts, periodicEdges, periodicFaces;
130  expList->GetPeriodicEntities(
131  periodicVerts, periodicEdges, periodicFaces);
132 
133  // The vectors below are of size 3 to have separate storage for
134  // vertices, edges and faces.
135 
136  // Maps from geometry ID to the matrix representing the extracted
137  // portion of S_1. For example idMats[2] folds the S_1 face blocks.
138  vector<map<int, vector<NekDouble> > > idMats(3);
139 
140  // Maps from the global ID, as obtained from AssemblyMapCG's
141  // localToGlobalMap, to the geometry ID.
142  vector<map<int, int> > gidMeshIds(3);
143 
144  // Maps from the global ID to the number of degrees of freedom for
145  // this geometry object.
146  vector<map<int, int> > gidDofs(3);
147 
148  // Array containing maximum information needed for the universal
149  // numbering later. For i = 0,1,2 for each geometry dimension:
150  //
151  // maxVertIds[2*i] = maximum geometry ID at dimension i
152  // maxVertIds[2*i+1] = maximum number of degrees of freedom for all
153  // elements of dimension i.
154  Array<OneD, int> maxVertIds(6, -1);
155 
156  // Figure out mapping from each elemental contribution to offset in
157  // (vert,edge,face) triples.
158  for (cnt = n = 0; n < nExp; ++n)
159  {
160  exp = expList->GetExp(n);
161 
162  // Grab reference to local Schur complement matrix.
163  DNekScalMatSharedPtr schurMat =
164  m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0,0);
165 
166  // Process vertices to extract relevant portion of the Schur
167  // complement matrix.
168  for (i = 0; i < exp->GetNverts(); ++i)
169  {
170  meshVertId = exp->GetGeom()->GetVid(i);
171  int locId = exp->GetVertexMap(i);
172 
173  // Get the global ID of this vertex.
174  gId = asmMap->GetLocalToGlobalMap(
175  cnt + locId) - nDirBnd;
176 
177  // Ignore all Dirichlet vertices.
178  if (gId < 0)
179  {
180  continue;
181  }
182 
183  gidDofs[0][gId] = 1;
184 
185  // Extract vertex value from Schur complement matrix.
186  NekDouble vertVal = (*schurMat)(locId,locId);
187 
188  // See if we have processed this vertex from another
189  // element.
190  auto gIt = idMats[0].find(gId);
191 
192  if (gIt == idMats[0].end())
193  {
194  // If not then put our 'matrix' inside idMats.
195  idMats[0][gId] = vector<NekDouble>(1, vertVal);
196  }
197  else
198  {
199  // Otherwise combine with the value that is already
200  // there (i.e. do assembly on this degree of freedom).
201  gIt->second[0] += vertVal;
202  }
203 
204  // Now check to see if the vertex is periodic. If it is,
205  // then we change meshVertId to be the minimum of all the
206  // other periodic vertices, so that we don't end up
207  // duplicating the matrix in our final block matrix.
208  auto pIt = periodicVerts.find(meshVertId);
209  if (pIt != periodicVerts.end())
210  {
211  for (j = 0; j < pIt->second.size(); ++j)
212  {
213  meshVertId = min(meshVertId, pIt->second[j].id);
214  }
215  }
216 
217  // Finally record the other information we need into the
218  // other matrices.
219  gidMeshIds[0][gId] = meshVertId;
220  maxVertIds[0] = max(maxVertIds[0], meshVertId);
221  maxVertIds[1] = 1;
222  }
223 
224  // Process edges. This logic is mostly the same as the previous
225  // block.
226  for (i = 0; i < exp->GetNedges(); ++i)
227  {
228  meshEdgeId = exp->GetGeom()->GetEid(i);
229 
230  Array<OneD, unsigned int> bmap, bmap2;
232  StdRegions::Orientation edgeOrient =
233  exp->GetGeom()->GetEorient(i);
234 
235  // Check if this edge is periodic. We may need to flip
236  // orientation if it is.
237  auto pIt = periodicEdges.find(meshEdgeId);
238  if (pIt != periodicEdges.end())
239  {
240  pair<int, StdRegions::Orientation> idOrient =
242  meshEdgeId, edgeOrient, pIt->second);
243  meshEdgeId = idOrient.first;
244  edgeOrient = idOrient.second;
245  }
246 
247  // Grab edge interior map, and the edge inverse boundary
248  // map, so that we can extract this edge from the Schur
249  // complement matrix.
250  exp->GetEdgeInteriorMap(i, edgeOrient, bmap, sign);
251  bmap2 = exp->GetEdgeInverseBoundaryMap(i);
252 
253  // Allocate temporary storage for the extracted edge matrix.
254  const int nEdgeCoeffs = bmap.num_elements();
255  vector<NekDouble> tmpStore(nEdgeCoeffs*nEdgeCoeffs);
256 
257  gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
258 
259  for (j = 0; j < nEdgeCoeffs; ++j)
260  {
261  // We record the minimum ID from the edge for our
262  // maps. This follows the logic that the assembly map
263  // ordering will always give us a contiguous ordering of
264  // global degrees of freedom for edge interior
265  // coefficients.
266  gId = min(gId,
267  asmMap->GetLocalToGlobalMap(
268  cnt + bmap[j])
269  - nDirBnd);
270 
271  // Ignore Dirichlet edges.
272  if (gId < 0)
273  {
274  continue;
275  }
276 
277  const NekDouble sign1 = sign[j];
278 
279  // Extract this edge, along with sign array for assembly
280  // later.
281  for (k = 0; k < nEdgeCoeffs; ++k)
282  {
283  tmpStore[k+j*nEdgeCoeffs] =
284  sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
285  }
286  }
287 
288  if (gId < 0)
289  {
290  continue;
291  }
292 
293  gidDofs[1][gId] = nEdgeCoeffs;
294 
295  // Assemble this edge matrix with another one, if it exists.
296  auto gIt = idMats[1].find(gId);
297  if (gIt == idMats[1].end())
298  {
299  idMats[1][gId] = tmpStore;
300  }
301  else
302  {
303  ASSERTL1(tmpStore.size() == gIt->second.size(),
304  "Number of modes mismatch");
305  Vmath::Vadd(nEdgeCoeffs*nEdgeCoeffs, &gIt->second[0], 1,
306  &tmpStore[0], 1, &gIt->second[0], 1);
307  }
308 
309  gidMeshIds[1][gId] = meshEdgeId;
310  maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
311  maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
312  }
313 
314  // Process faces. This logic is mostly the same as the previous
315  // block.
316  for (i = 0; i < exp->GetNfaces(); ++i)
317  {
318  meshFaceId = exp->GetGeom()->GetFid(i);
319 
320  Array<OneD, unsigned int> bmap, bmap2;
322  StdRegions::Orientation faceOrient =
323  exp->GetGeom()->GetForient(i);
324 
325  // Check if this face is periodic. We may need to flip
326  // orientation if it is.
327  auto pIt = periodicFaces.find(meshFaceId);
328  if (pIt != periodicFaces.end())
329  {
330  meshFaceId = min(meshFaceId, pIt->second[0].id);
331  faceOrient = DeterminePeriodicFaceOrient(
332  faceOrient, pIt->second[0].orient);
333  }
334 
335  exp->GetFaceInteriorMap(i, faceOrient, bmap, sign);
336  bmap2 = exp->GetFaceInverseBoundaryMap(i);
337 
338  // Allocate temporary storage for the extracted face matrix.
339  const int nFaceCoeffs = bmap.num_elements();
340  vector<NekDouble> tmpStore(nFaceCoeffs*nFaceCoeffs);
341 
342  gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
343 
344  for (j = 0; j < nFaceCoeffs; ++j)
345  {
346  gId = min(gId,
347  asmMap->GetLocalToGlobalMap(cnt + bmap[j])
348  - nDirBnd);
349 
350  // Ignore Dirichlet faces.
351  if (gId < 0)
352  {
353  continue;
354  }
355 
356  const NekDouble sign1 = sign[j];
357 
358  // Extract this face, along with sign array for assembly
359  // later.
360  for (k = 0; k < nFaceCoeffs; ++k)
361  {
362  tmpStore[k+j*nFaceCoeffs] =
363  sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
364  }
365  }
366 
367  if (gId < 0)
368  {
369  continue;
370  }
371 
372  gidDofs[2][gId] = nFaceCoeffs;
373 
374  // Assemble this face matrix with another one, if it exists.
375  auto gIt = idMats[2].find(gId);
376  if (gIt == idMats[2].end())
377  {
378  idMats[2][gId] = tmpStore;
379  }
380  else
381  {
382  ASSERTL1(tmpStore.size() == gIt->second.size(),
383  "Number of modes mismatch");
384  Vmath::Vadd(nFaceCoeffs*nFaceCoeffs, &gIt->second[0], 1,
385  &tmpStore[0], 1, &gIt->second[0], 1);
386  }
387 
388  gidMeshIds[2][gId] = meshFaceId;
389  maxVertIds[4] = max(maxVertIds[4], meshFaceId);
390  maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
391  }
392 
393  cnt += exp->GetNcoeffs();
394  }
395 
396  // Perform a reduction to find maximum vertex, edge and face
397  // geometry IDs.
398  m_comm = expList->GetSession()->GetComm()->GetRowComm();
399  m_comm->AllReduce(maxVertIds, LibUtilities::ReduceMax);
400 
401  // Concatenate all matrices into contiguous storage and figure out
402  // universal ID numbering.
403  vector<NekDouble> storageBuf;
404  vector<long> globalToUniversal;
405 
406  for (i = 0, cnt = 1; i < 3; ++i)
407  {
408  const int maxDofs = maxVertIds[2*i+1];
409 
410  // Note that iterating over the map uses the property that keys
411  // are accessed in order of ascending order, putting everything
412  // in the right order for the global system.
413  for (auto &gIt : idMats[i])
414  {
415  // Copy matrix into storage.
416  storageBuf.insert(storageBuf.end(),
417  gIt.second.begin(), gIt.second.end());
418 
419  // Get mesh ID from global ID number.
420  ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
421  "Unable to find global ID " +
422  boost::lexical_cast<string>(gIt.first) +
423  " inside map");
424  meshVertId = gidMeshIds[i][gIt.first];
425 
426  for (j = 0; j < gIt.second.size(); ++j)
427  {
428  globalToUniversal.push_back(
429  cnt + meshVertId*maxDofs*maxDofs + j);
430  }
431 
432  // Free up the temporary storage.
433  gIt.second.clear();
434  }
435 
436  cnt += (maxVertIds[2*i]+1)*maxDofs*maxDofs;
437  }
438 
439  ASSERTL1(storageBuf.size() == globalToUniversal.size(),
440  "Storage buffer and global to universal map size does "
441  "not match");
442 
443  Array<OneD, NekDouble> storageData(
444  storageBuf.size(), &storageBuf[0]);
445  Array<OneD, long> globalToUniversalMap(
446  globalToUniversal.size(), &globalToUniversal[0]);
447 
448  // Use GS to assemble data between processors.
449  Gs::gs_data *tmpGs = Gs::Init(
450  globalToUniversalMap, m_comm,
451  expList->GetSession()->DefinesCmdLineArgument("verbose"));
452  Gs::Gather(storageData, Gs::gs_add, tmpGs);
453 
454  // Figure out what storage we need in the block matrix.
456  1 + idMats[1].size() + idMats[2].size());
457 
458  // Vertex block is a diagonal matrix.
459  n_blks[0] = idMats[0].size();
460 
461  // Now extract number of rows in each edge and face block from the
462  // gidDofs map.
463  cnt = 1;
464  for (i = 1; i < 3; ++i)
465  {
466  for (auto &gIt : idMats[i])
467  {
468  ASSERTL1(gidDofs[i].count(gIt.first) > 0,
469  "Unable to find number of degrees of freedom for "
470  "global ID " + boost::lexical_cast<string>(
471  gIt.first));
472 
473  n_blks[cnt++] = gidDofs[i][gIt.first];
474  }
475  }
476 
477  // Allocate storage for the block matrix.
479  ::AllocateSharedPtr(n_blks, n_blks, eDIAGONAL);
480 
481  // We deal with the vertex matrix separately since all vertices can
482  // be condensed into a single, block-diagonal matrix.
484  ::AllocateSharedPtr(n_blks[0], n_blks[0], 0.0, eDIAGONAL);
485 
486  // Fill the vertex matrix with the inverse of each vertex value.
487  cnt = 0;
488  for (auto gIt = idMats[0].begin(); gIt != idMats[0].end();
489  ++gIt, ++cnt)
490  {
491  (*vertMat)(cnt, cnt) = 1.0/storageData[cnt];
492  }
493 
494  // Put the vertex matrix in the block matrix.
495  m_blkMat->SetBlock(0,0,vertMat);
496 
497  // Finally, grab the matrices from the block storage, invert them
498  // and place them in the correct position inside the block matrix.
499  int cnt2 = 1;
500  for (i = 1; i < 3; ++i)
501  {
502  for (auto &gIt : idMats[i])
503  {
504  int nDofs = gidDofs[i][gIt.first];
505 
507  ::AllocateSharedPtr(nDofs, nDofs);
508 
509  for (j = 0; j < nDofs; ++j)
510  {
511  for (k = 0; k < nDofs; ++k)
512  {
513  (*tmp)(j,k) = storageData[k+j*nDofs + cnt];
514  }
515  }
516 
517  cnt += nDofs*nDofs;
518 
519  tmp->Invert();
520  m_blkMat->SetBlock(cnt2, cnt2, tmp);
521  ++cnt2;
522  }
523  }
524  }
525 
526  /**
527  * @brief Construct a block preconditioner for the hybridized
528  * discontinuous Galerkin system.
529  *
530  * This system is built in a similar fashion to its continuous variant
531  * found in PreconditionerBlock::BlockPreconditionerCG. In this setting
532  * however, the matrix is constructed as:
533  *
534  * \f[ M^{-1} = \mathrm{Diag}[ (\mathbf{S_{1}})_{f}^{-1} ] \f]
535  *
536  * where each matrix is the Schur complement system restricted to a
537  * single face of the trace system.
538  */
540  {
541  std::shared_ptr<MultiRegions::ExpList>
542  expList = ((m_linsys.lock())->GetLocMat()).lock();
543  std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
545  DNekScalBlkMatSharedPtr loc_mat;
546  DNekScalMatSharedPtr bnd_mat;
547 
548  AssemblyMapDGSharedPtr asmMap = std::dynamic_pointer_cast<
549  AssemblyMapDG>(m_locToGloMap.lock());
550 
551  int i, j, k, n, cnt, cnt2;
552 
553  // Figure out number of Dirichlet trace elements
554  int nTrace = expList->GetTrace()->GetExpSize();
555  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
556 
557  for (cnt = n = 0; n < nTrace; ++n)
558  {
559  if (cnt >= nDir)
560  {
561  break;
562  }
563 
564  cnt += trace->GetExp(n)->GetNcoeffs();
565  }
566 
567  nDir = n;
568 
569  // Allocate storage for block matrix. Need number of unique faces in
570  // trace space.
571  int maxTraceSize = -1;
572  map<int, int> arrayOffsets;
573 
574  for (cnt = 0, n = nDir; n < nTrace; ++n)
575  {
576  int nCoeffs = trace->GetExp(n)->GetNcoeffs();
577  int nCoeffs2 = nCoeffs * nCoeffs;
578 
579  arrayOffsets[n] = cnt;
580  cnt += nCoeffs2;
581 
582  if (maxTraceSize < nCoeffs2)
583  {
584  maxTraceSize = nCoeffs2;
585  }
586  }
587 
588  // Find maximum trace size.
589  m_comm = expList->GetSession()->GetComm()->GetRowComm();
590  m_comm->AllReduce(maxTraceSize, LibUtilities::ReduceMax);
591 
592  // Zero matrix storage.
593  Array<OneD, NekDouble> tmpStore(cnt, 0.0);
594 
595  // Assemble block matrices for each trace element.
596  for (cnt = n = 0; n < expList->GetExpSize(); ++n)
597  {
598  int elmt = n;
599  locExpansion = expList->GetExp(elmt);
600 
602  asmMap->GetElmtToTrace()[elmt];
603 
604  // Block matrix (lambda)
605  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
606  bnd_mat = loc_mat->GetBlock(0,0);
607 
608  int nFacets = locExpansion->GetNumBases() == 2 ?
609  locExpansion->GetNedges() : locExpansion->GetNfaces();
610 
611  for (cnt2 = i = 0; i < nFacets; ++i)
612  {
613  int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
614  int elmtId = elmtToTraceMap[i]->GetElmtId ();
615 
616  // Ignore Dirichlet edges.
617  if (elmtId < nDir)
618  {
619  cnt += nCoeffs;
620  cnt2 += nCoeffs;
621  continue;
622  }
623 
624  NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
625 
626  for (j = 0; j < nCoeffs; ++j)
627  {
628  NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(
629  cnt + j);
630  for (k = 0; k < nCoeffs; ++k)
631  {
632  NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(
633  cnt + k);
634  off[j*nCoeffs + k] +=
635  (*bnd_mat)(cnt2+j, cnt2+k) * sign1 * sign2;
636  }
637  }
638 
639  cnt += nCoeffs;
640  cnt2 += nCoeffs;
641  }
642  }
643 
644  // Set up IDs for universal numbering.
645  Array<OneD, long> uniIds(tmpStore.num_elements());
646  for (cnt = 0, n = nDir; n < nTrace; ++n)
647  {
648  LocalRegions::ExpansionSharedPtr traceExp = trace->GetExp(n);
649  int nCoeffs = traceExp->GetNcoeffs();
650  int geomId = traceExp->GetGeom()->GetGlobalID();
651 
652  for (i = 0; i < nCoeffs*nCoeffs; ++i)
653  {
654  uniIds[cnt++] = geomId * maxTraceSize + i + 1;
655  }
656  }
657 
658  // Assemble matrices across partitions.
659  Gs::gs_data *gsh = Gs::Init(
660  uniIds, m_comm,
661  expList->GetSession()->DefinesCmdLineArgument("verbose"));
662  Gs::Gather(tmpStore, Gs::gs_add, gsh);
663 
664  // Set up diagonal block matrix
665  Array<OneD, unsigned int> n_blks(nTrace - nDir);
666  for (n = 0; n < nTrace - nDir; ++n)
667  {
668  n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
669  }
670 
672  ::AllocateSharedPtr(n_blks, n_blks, eDIAGONAL);
673 
674  for (n = 0; n < nTrace - nDir; ++n)
675  {
676  int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
678  ::AllocateSharedPtr(nCoeffs, nCoeffs);
679  NekDouble *off = &tmpStore[arrayOffsets[n+nDir]];
680 
681  for (i = 0; i < nCoeffs; ++i)
682  {
683  for (j = 0; j < nCoeffs; ++j)
684  {
685  (*tmp)(i,j) = off[i*nCoeffs + j];
686  }
687  }
688 
689  tmp->Invert();
690  m_blkMat->SetBlock(n, n, tmp);
691  }
692  }
693 
694  /**
695  * @brief Apply preconditioner to \p pInput and store the result in \p
696  * pOutput.
697  */
699  const Array<OneD, NekDouble>& pInput,
700  Array<OneD, NekDouble>& pOutput)
701  {
702  int nDir = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
703  int nGlobal = m_locToGloMap.lock()->GetNumGlobalBndCoeffs();
704  int nNonDir = nGlobal-nDir;
705  DNekBlkMat &M = (*m_blkMat);
706  NekVector<NekDouble> r(nNonDir,pInput,eWrapper);
707  NekVector<NekDouble> z(nNonDir,pOutput,eWrapper);
708  z = M * r;
709  }
710  }
711 }
LibUtilities::CommSharedPtr m_comm
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply preconditioner to pInput and store the result in pOutput.
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:245
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
const std::weak_ptr< GlobalLinSys > m_linsys
PreconFactory & GetPreconFactory()
STL namespace.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:167
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
StdRegions::Orientation DeterminePeriodicFaceOrient(StdRegions::Orientation faceOrient, StdRegions::Orientation perFaceOrient)
Determine relative orientation between two faces.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
double NekDouble
void BlockPreconditionerHDG(void)
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
Describe a linear system.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
void BlockPreconditionerCG(void)
Construct a block preconditioner from for the continuous Galerkin system.
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:48
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
pair< int, StdRegions::Orientation > DeterminePeriodicEdgeOrientId(int meshEdgeId, StdRegions::Orientation edgeOrient, const vector< PeriodicEntity > &periodicEdges)
Determine orientation of an edge to its periodic equivalents, as well as the ID of the representative...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
std::weak_ptr< AssemblyMap > m_locToGloMap
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:302