Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Attributes | Private Member Functions | List of all members
Nektar::MultiRegions::PreconditionerBlock Class Reference

#include <PreconditionerBlock.h>

Inheritance diagram for Nektar::MultiRegions::PreconditionerBlock:
Inheritance graph
[legend]
Collaboration diagram for Nektar::MultiRegions::PreconditionerBlock:
Collaboration graph
[legend]

Public Member Functions

 PreconditionerBlock (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~PreconditionerBlock ()
 
- Public Member Functions inherited from Nektar::MultiRegions::Preconditioner
 Preconditioner (const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
 
virtual ~Preconditioner ()
 
void DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoPreconditionerWithNonVertOutput (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce=NullNekDouble1DArray)
 
void DoTransformToLowEnergy (Array< OneD, NekDouble > &pInOut, int offset)
 
void DoTransformToLowEnergy (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoTransformFromLowEnergy (Array< OneD, NekDouble > &pInOut)
 
void DoMultiplybyInverseTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void DoMultiplybyInverseTransposedTransformationMatrix (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 
void BuildPreconditioner ()
 
void InitObject ()
 
Array< OneD, NekDoubleAssembleStaticCondGlobalDiagonals ()
 Performs global assembly of diagonal entries to global Schur complement matrix. More...
 
const DNekScalBlkMatSharedPtrGetBlockTransformedSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockCMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockInvDMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockSchurCompl () const
 
const DNekScalBlkMatSharedPtrGetBlockTransformationMatrix () const
 
const DNekScalBlkMatSharedPtrGetBlockTransposedTransformationMatrix () const
 
DNekScalMatSharedPtr TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
 

Static Public Member Functions

static PreconditionerSharedPtr create (const boost::shared_ptr< GlobalLinSys > &plinsys, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Attributes

const boost::weak_ptr
< GlobalLinSys
m_linsys
 
PreconditionerType m_preconType
 
DNekBlkMatSharedPtr m_blkMat
 
boost::shared_ptr< AssemblyMapm_locToGloMap
 
- Protected Attributes inherited from Nektar::MultiRegions::Preconditioner
const boost::weak_ptr
< GlobalLinSys
m_linsys
 
PreconditionerType m_preconType
 
DNekMatSharedPtr m_preconditioner
 
boost::shared_ptr< AssemblyMapm_locToGloMap
 
LibUtilities::CommSharedPtr m_comm
 

Private Member Functions

void BlockPreconditionerCG (void)
 Construct a block preconditioner from $\mathbf{S}_{1}$ for the continuous Galerkin system. More...
 
void BlockPreconditionerHDG (void)
 
virtual void v_InitObject ()
 
virtual void v_DoPreconditioner (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
 Apply preconditioner to pInput and store the result in pOutput. More...
 
virtual void v_BuildPreconditioner ()
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::MultiRegions::Preconditioner
virtual DNekScalMatSharedPtr v_TransformedSchurCompl (int offset, const boost::shared_ptr< DNekScalMat > &loc_mat)
 Get block elemental transposed transformation matrix $\mathbf{R}^{T}$. More...
 

Detailed Description

Definition at line 52 of file PreconditionerBlock.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::PreconditionerBlock::PreconditionerBlock ( const boost::shared_ptr< GlobalLinSys > &  plinsys,
const AssemblyMapSharedPtr pLocToGloMap 
)

Definition at line 65 of file PreconditionerBlock.cpp.

68  : Preconditioner(plinsys, pLocToGloMap),
69  m_linsys(plinsys),
70  m_preconType(pLocToGloMap->GetPreconType()),
71  m_locToGloMap(pLocToGloMap)
72  {
73  }
Preconditioner(const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
boost::shared_ptr< AssemblyMap > m_locToGloMap
const boost::weak_ptr< GlobalLinSys > m_linsys
virtual Nektar::MultiRegions::PreconditionerBlock::~PreconditionerBlock ( )
inlinevirtual

Definition at line 75 of file PreconditionerBlock.h.

75 {}

Member Function Documentation

void Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerCG ( void  )
private

Construct a block preconditioner from $\mathbf{S}_{1}$ for the continuous Galerkin system.

The preconditioner is defined as:

\[\mathbf{M}^{-1}=\left[\begin{array}{ccc} Diag[(\mathbf{S_{1}})_{vv}] & & \\ & (\mathbf{S}_{1})_{eb} & \\ & & (\mathbf{S}_{1})_{fb} \end{array}\right] \]

where $\mathbf{S}_{1}$ is the local Schur complement matrix for each element and the subscript denotes the portion of the Schur complement associated with the vertex, edge and face blocks respectively.

Definition at line 113 of file PreconditionerBlock.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, Nektar::MultiRegions::DeterminePeriodicEdgeOrientId(), Nektar::MultiRegions::DeterminePeriodicFaceOrient(), Nektar::eDIAGONAL, Gs::Gather(), Gs::gs_add, Gs::Init(), Nektar::iterator, m_blkMat, Nektar::MultiRegions::Preconditioner::m_comm, m_linsys, m_locToGloMap, Nektar::LibUtilities::ReduceMax, sign, and Vmath::Vadd().

Referenced by v_BuildPreconditioner().

114  {
115  ExpListSharedPtr expList = m_linsys.lock()->GetLocMat().lock();
117  DNekScalBlkMatSharedPtr loc_mat;
118  DNekScalMatSharedPtr bnd_mat;
119 
120  int nel, i, j, k, n, cnt, gId;
121  int meshVertId, meshEdgeId, meshFaceId;
122 
123  const int nExp = expList->GetExpSize();
124  const int nDirBnd = m_locToGloMap->GetNumGlobalDirBndCoeffs();
125 
126  // Grab periodic geometry information.
128  PeriodicMap periodicVerts, periodicEdges, periodicFaces;
129  expList->GetPeriodicEntities(
130  periodicVerts, periodicEdges, periodicFaces);
131 
132  // The vectors below are of size 3 to have separate storage for
133  // vertices, edges and faces.
134 
135  // Maps from geometry ID to the matrix representing the extracted
136  // portion of S_1. For example idMats[2] folds the S_1 face blocks.
137  vector<map<int, vector<NekDouble> > > idMats(3);
138 
139  // Maps from the global ID, as obtained from AssemblyMapCG's
140  // localToGlobalMap, to the geometry ID.
141  vector<map<int, int> > gidMeshIds(3);
142 
143  // Maps from the global ID to the number of degrees of freedom for
144  // this geometry object.
145  vector<map<int, int> > gidDofs(3);
146 
147  // Array containing maximum information needed for the universal
148  // numbering later. For i = 0,1,2 for each geometry dimension:
149  //
150  // maxVertIds[2*i] = maximum geometry ID at dimension i
151  // maxVertIds[2*i+1] = maximum number of degrees of freedom for all
152  // elements of dimension i.
153  Array<OneD, int> maxVertIds(6, -1);
154 
155  // Iterator for idMats members.
156  map<int, vector<NekDouble> >::iterator gIt;
157 
158  // Figure out mapping from each elemental contribution to offset in
159  // (vert,edge,face) triples.
160  for (cnt = n = 0; n < nExp; ++n)
161  {
162  nel = expList->GetOffset_Elmt_Id(n);
163  exp = expList->GetExp(nel);
164 
165  // Grab reference to local Schur complement matrix.
166  DNekScalMatSharedPtr schurMat =
167  m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0,0);
168 
169  // Process vertices to extract relevant portion of the Schur
170  // complement matrix.
171  for (i = 0; i < exp->GetNverts(); ++i)
172  {
173  meshVertId = exp->GetGeom()->GetVid(i);
174  int locId = exp->GetVertexMap(i);
175 
176  // Get the global ID of this vertex.
177  gId = m_locToGloMap->GetLocalToGlobalMap(
178  cnt + locId) - nDirBnd;
179 
180  // Ignore all Dirichlet vertices.
181  if (gId < 0)
182  {
183  continue;
184  }
185 
186  gidDofs[0][gId] = 1;
187 
188  // Extract vertex value from Schur complement matrix.
189  NekDouble vertVal = (*schurMat)(locId,locId);
190 
191  // See if we have processed this vertex from another
192  // element.
193  gIt = idMats[0].find(gId);
194 
195  if (gIt == idMats[0].end())
196  {
197  // If not then put our 'matrix' inside idMats.
198  idMats[0][gId] = vector<NekDouble>(1, vertVal);
199  }
200  else
201  {
202  // Otherwise combine with the value that is already
203  // there (i.e. do assembly on this degree of freedom).
204  gIt->second[0] += vertVal;
205  }
206 
207  // Now check to see if the vertex is periodic. If it is,
208  // then we change meshVertId to be the minimum of all the
209  // other periodic vertices, so that we don't end up
210  // duplicating the matrix in our final block matrix.
211  pIt = periodicVerts.find(meshVertId);
212  if (pIt != periodicVerts.end())
213  {
214  for (j = 0; j < pIt->second.size(); ++j)
215  {
216  meshVertId = min(meshVertId, pIt->second[j].id);
217  }
218  }
219 
220  // Finally record the other information we need into the
221  // other matrices.
222  gidMeshIds[0][gId] = meshVertId;
223  maxVertIds[0] = max(maxVertIds[0], meshVertId);
224  maxVertIds[1] = 1;
225  }
226 
227  // Process edges. This logic is mostly the same as the previous
228  // block.
229  for (i = 0; i < exp->GetNedges(); ++i)
230  {
231  meshEdgeId = exp->GetGeom()->GetEid(i);
232 
233  Array<OneD, unsigned int> bmap, bmap2;
235  StdRegions::Orientation edgeOrient =
236  exp->GetGeom()->GetEorient(i);
237 
238  // Check if this edge is periodic. We may need to flip
239  // orientation if it is.
240  pIt = periodicEdges.find(meshEdgeId);
241  if (pIt != periodicEdges.end())
242  {
243  pair<int, StdRegions::Orientation> idOrient =
245  meshEdgeId, edgeOrient, pIt->second);
246  meshEdgeId = idOrient.first;
247  edgeOrient = idOrient.second;
248  }
249 
250  // Grab edge interior map, and the edge inverse boundary
251  // map, so that we can extract this edge from the Schur
252  // complement matrix.
253  exp->GetEdgeInteriorMap(i, edgeOrient, bmap, sign);
254  bmap2 = exp->GetEdgeInverseBoundaryMap(i);
255 
256  // Allocate temporary storage for the extracted edge matrix.
257  const int nEdgeCoeffs = bmap.num_elements();
258  vector<NekDouble> tmpStore(nEdgeCoeffs*nEdgeCoeffs);
259 
260  gId = m_locToGloMap->GetLocalToGlobalMap(cnt + bmap[0]);
261 
262  for (j = 0; j < nEdgeCoeffs; ++j)
263  {
264  // We record the minimum ID from the edge for our
265  // maps. This follows the logic that the assembly map
266  // ordering will always give us a contiguous ordering of
267  // global degrees of freedom for edge interior
268  // coefficients.
269  gId = min(gId,
270  m_locToGloMap->GetLocalToGlobalMap(
271  cnt + bmap[j])
272  - nDirBnd);
273 
274  // Ignore Dirichlet edges.
275  if (gId < 0)
276  {
277  continue;
278  }
279 
280  const NekDouble sign1 = sign[j];
281 
282  // Extract this edge, along with sign array for assembly
283  // later.
284  for (k = 0; k < nEdgeCoeffs; ++k)
285  {
286  tmpStore[k+j*nEdgeCoeffs] =
287  sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
288  }
289  }
290 
291  if (gId < 0)
292  {
293  continue;
294  }
295 
296  gidDofs[1][gId] = nEdgeCoeffs;
297 
298  // Assemble this edge matrix with another one, if it exists.
299  gIt = idMats[1].find(gId);
300  if (gIt == idMats[1].end())
301  {
302  idMats[1][gId] = tmpStore;
303  }
304  else
305  {
306  ASSERTL1(tmpStore.size() == gIt->second.size(),
307  "Number of modes mismatch");
308  Vmath::Vadd(nEdgeCoeffs*nEdgeCoeffs, &gIt->second[0], 1,
309  &tmpStore[0], 1, &gIt->second[0], 1);
310  }
311 
312  gidMeshIds[1][gId] = meshEdgeId;
313  maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
314  maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
315  }
316 
317  // Process faces. This logic is mostly the same as the previous
318  // block.
319  for (i = 0; i < exp->GetNfaces(); ++i)
320  {
321  meshFaceId = exp->GetGeom()->GetFid(i);
322 
323  Array<OneD, unsigned int> bmap, bmap2;
325  StdRegions::Orientation faceOrient =
326  exp->GetGeom()->GetForient(i);
327 
328  // Check if this face is periodic. We may need to flip
329  // orientation if it is.
330  pIt = periodicFaces.find(meshFaceId);
331  if (pIt != periodicFaces.end())
332  {
333  meshFaceId = min(meshFaceId, pIt->second[0].id);
334  faceOrient = DeterminePeriodicFaceOrient(
335  faceOrient, pIt->second[0].orient);
336  }
337 
338  exp->GetFaceInteriorMap(i, faceOrient, bmap, sign);
339  bmap2 = exp->GetFaceInverseBoundaryMap(i);
340 
341  // Allocate temporary storage for the extracted face matrix.
342  const int nFaceCoeffs = bmap.num_elements();
343  vector<NekDouble> tmpStore(nFaceCoeffs*nFaceCoeffs);
344 
345  gId = m_locToGloMap->GetLocalToGlobalMap(cnt + bmap[0]);
346 
347  for (j = 0; j < nFaceCoeffs; ++j)
348  {
349  gId = min(gId,
350  m_locToGloMap->GetLocalToGlobalMap(
351  cnt + bmap[j])
352  - nDirBnd);
353 
354  // Ignore Dirichlet faces.
355  if (gId < 0)
356  {
357  continue;
358  }
359 
360  const NekDouble sign1 = sign[j];
361 
362  // Extract this face, along with sign array for assembly
363  // later.
364  for (k = 0; k < nFaceCoeffs; ++k)
365  {
366  tmpStore[k+j*nFaceCoeffs] =
367  sign1*sign[k]*(*schurMat)(bmap2[j], bmap2[k]);
368  }
369  }
370 
371  if (gId < 0)
372  {
373  continue;
374  }
375 
376  gidDofs[2][gId] = nFaceCoeffs;
377 
378  // Assemble this face matrix with another one, if it exists.
379  gIt = idMats[2].find(gId);
380  if (gIt == idMats[2].end())
381  {
382  idMats[2][gId] = tmpStore;
383  }
384  else
385  {
386  ASSERTL1(tmpStore.size() == gIt->second.size(),
387  "Number of modes mismatch");
388  Vmath::Vadd(nFaceCoeffs*nFaceCoeffs, &gIt->second[0], 1,
389  &tmpStore[0], 1, &gIt->second[0], 1);
390  }
391 
392  gidMeshIds[2][gId] = meshFaceId;
393  maxVertIds[4] = max(maxVertIds[4], meshFaceId);
394  maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
395  }
396 
397  cnt += exp->GetNcoeffs();
398  }
399 
400  // Perform a reduction to find maximum vertex, edge and face
401  // geometry IDs.
402  m_comm = expList->GetSession()->GetComm()->GetRowComm();
403  m_comm->AllReduce(maxVertIds, LibUtilities::ReduceMax);
404 
405  // Concatenate all matrices into contiguous storage and figure out
406  // universal ID numbering.
407  vector<NekDouble> storageBuf;
408  vector<long> globalToUniversal;
409 
410  for (i = 0, cnt = 1; i < 3; ++i)
411  {
412  const int maxDofs = maxVertIds[2*i+1];
413 
414  // Note that iterating over the map uses the property that keys
415  // are accessed in order of ascending order, putting everything
416  // in the right order for the global system.
417  for (gIt = idMats[i].begin(); gIt != idMats[i].end(); ++gIt)
418  {
419  // Copy matrix into storage.
420  storageBuf.insert(storageBuf.end(),
421  gIt->second.begin(), gIt->second.end());
422 
423  // Get mesh ID from global ID number.
424  ASSERTL1(gidMeshIds[i].count(gIt->first) > 0,
425  "Unable to find global ID " +
426  boost::lexical_cast<string>(gIt->first) +
427  " inside map");
428  meshVertId = gidMeshIds[i][gIt->first];
429 
430  for (j = 0; j < gIt->second.size(); ++j)
431  {
432  globalToUniversal.push_back(
433  cnt + meshVertId*maxDofs*maxDofs + j);
434  }
435 
436  // Free up the temporary storage.
437  gIt->second.clear();
438  }
439 
440  cnt += (maxVertIds[2*i]+1)*maxDofs*maxDofs;
441  }
442 
443  ASSERTL1(storageBuf.size() == globalToUniversal.size(),
444  "Storage buffer and global to universal map size does "
445  "not match");
446 
447  Array<OneD, NekDouble> storageData(
448  storageBuf.size(), &storageBuf[0]);
449  Array<OneD, long> globalToUniversalMap(
450  globalToUniversal.size(), &globalToUniversal[0]);
451 
452  // Use GS to assemble data between processors.
453  Gs::gs_data *tmpGs = Gs::Init(globalToUniversalMap, m_comm);
454  Gs::Gather(storageData, Gs::gs_add, tmpGs);
455 
456  // Figure out what storage we need in the block matrix.
458  1 + idMats[1].size() + idMats[2].size());
459 
460  // Vertex block is a diagonal matrix.
461  n_blks[0] = idMats[0].size();
462 
463  // Now extract number of rows in each edge and face block from the
464  // gidDofs map.
465  cnt = 1;
466  for (i = 1; i < 3; ++i)
467  {
468  for (gIt = idMats[i].begin(); gIt != idMats[i].end(); ++gIt)
469  {
470  ASSERTL1(gidDofs[i].count(gIt->first) > 0,
471  "Unable to find number of degrees of freedom for "
472  "global ID " + boost::lexical_cast<string>(
473  gIt->first));
474 
475  n_blks[cnt++] = gidDofs[i][gIt->first];
476  }
477  }
478 
479  // Allocate storage for the block matrix.
481  ::AllocateSharedPtr(n_blks, n_blks, eDIAGONAL);
482 
483  // We deal with the vertex matrix separately since all vertices can
484  // be condensed into a single, block-diagonal matrix.
486  ::AllocateSharedPtr(n_blks[0], n_blks[0], 0.0, eDIAGONAL);
487 
488  // Fill the vertex matrix with the inverse of each vertex value.
489  cnt = 0;
490  for (gIt = idMats[0].begin(); gIt != idMats[0].end(); ++gIt, ++cnt)
491  {
492  (*vertMat)(cnt, cnt) = 1.0/storageData[cnt];
493  }
494 
495  // Put the vertex matrix in the block matrix.
496  m_blkMat->SetBlock(0,0,vertMat);
497 
498  // Finally, grab the matrices from the block storage, invert them
499  // and place them in the correct position inside the block matrix.
500  int cnt2 = 1;
501  for (i = 1; i < 3; ++i)
502  {
503  for (gIt = idMats[i].begin(); gIt != idMats[i].end();
504  ++gIt, ++cnt2)
505  {
506  int nDofs = gidDofs[i][gIt->first];
507 
509  ::AllocateSharedPtr(nDofs, nDofs);
510 
511  for (j = 0; j < nDofs; ++j)
512  {
513  for (k = 0; k < nDofs; ++k)
514  {
515  (*tmp)(j,k) = storageData[k+j*nDofs + cnt];
516  }
517  }
518 
519  cnt += nDofs*nDofs;
520 
521  tmp->Invert();
522  m_blkMat->SetBlock(cnt2, cnt2, tmp);
523  }
524  }
525  }
LibUtilities::CommSharedPtr m_comm
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:223
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::map< int, vector< PeriodicEntity > > PeriodicMap
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:150
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
StdRegions::Orientation DeterminePeriodicFaceOrient(StdRegions::Orientation faceOrient, StdRegions::Orientation perFaceOrient)
Determine relative orientation between two faces.
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
double NekDouble
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...
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
boost::shared_ptr< AssemblyMap > m_locToGloMap
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
const boost::weak_ptr< GlobalLinSys > m_linsys
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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:285
void Nektar::MultiRegions::PreconditionerBlock::BlockPreconditionerHDG ( void  )
private

Definition at line 540 of file PreconditionerBlock.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eDIAGONAL, Gs::Gather(), Gs::gs_add, Gs::Init(), m_blkMat, Nektar::MultiRegions::Preconditioner::m_comm, m_linsys, m_locToGloMap, and Nektar::LibUtilities::ReduceMax.

Referenced by v_BuildPreconditioner().

541  {
542  boost::shared_ptr<MultiRegions::ExpList>
543  expList = ((m_linsys.lock())->GetLocMat()).lock();
544  boost::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
546  DNekScalBlkMatSharedPtr loc_mat;
547  DNekScalMatSharedPtr bnd_mat;
548 
549  AssemblyMapDGSharedPtr asmMap = boost::dynamic_pointer_cast<
550  AssemblyMapDG>(m_locToGloMap);
551 
552  int i, j, k, n, cnt, cnt2;
553 
554  // Figure out number of Dirichlet trace elements
555  int nTrace = expList->GetTrace()->GetExpSize();
556  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
557 
558  for (cnt = n = 0; n < nTrace; ++n)
559  {
560  if (cnt >= nDir)
561  {
562  break;
563  }
564 
565  cnt += trace->GetExp(n)->GetNcoeffs();
566  }
567 
568  nDir = n;
569 
570  // Allocate storage for block matrix. Need number of unique faces in
571  // trace space.
572  int maxTraceSize = -1;
573  map<int, int> arrayOffsets;
574 
575  for (cnt = 0, n = nDir; n < nTrace; ++n)
576  {
577  int nCoeffs = trace->GetExp(n)->GetNcoeffs();
578  int nCoeffs2 = nCoeffs * nCoeffs;
579 
580  arrayOffsets[n] = cnt;
581  cnt += nCoeffs2;
582 
583  if (maxTraceSize < nCoeffs2)
584  {
585  maxTraceSize = nCoeffs2;
586  }
587  }
588 
589  // Find maximum trace size.
590  m_comm = expList->GetSession()->GetComm()->GetRowComm();
591  m_comm->AllReduce(maxTraceSize, LibUtilities::ReduceMax);
592 
593  // Zero matrix storage.
594  Array<OneD, NekDouble> tmpStore(cnt, 0.0);
595 
596  // Assemble block matrices for each trace element.
597  for (cnt = n = 0; n < expList->GetExpSize(); ++n)
598  {
599  int elmt = expList->GetOffset_Elmt_Id(n);
600  locExpansion = expList->GetExp(elmt);
601 
603  asmMap->GetElmtToTrace()[elmt];
604 
605  // Block matrix (lambda)
606  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
607  bnd_mat = loc_mat->GetBlock(0,0);
608 
609  int nFacets = locExpansion->GetNumBases() == 2 ?
610  locExpansion->GetNedges() : locExpansion->GetNfaces();
611 
612  for (cnt2 = i = 0; i < nFacets; ++i)
613  {
614  int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
615  int elmtId = elmtToTraceMap[i]->GetElmtId ();
616 
617  // Ignore Dirichlet edges.
618  if (elmtId < nDir)
619  {
620  cnt += nCoeffs;
621  cnt2 += nCoeffs;
622  continue;
623  }
624 
625  NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
626 
627  for (j = 0; j < nCoeffs; ++j)
628  {
629  NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(
630  cnt + j);
631  for (k = 0; k < nCoeffs; ++k)
632  {
633  NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(
634  cnt + k);
635  off[j*nCoeffs + k] +=
636  (*bnd_mat)(cnt2+j, cnt2+k) * sign1 * sign2;
637  }
638  }
639 
640  cnt += nCoeffs;
641  cnt2 += nCoeffs;
642  }
643  }
644 
645  // Set up IDs for universal numbering.
646  Array<OneD, long> uniIds(tmpStore.num_elements());
647  for (cnt = 0, n = nDir; n < nTrace; ++n)
648  {
649  LocalRegions::ExpansionSharedPtr traceExp = trace->GetExp(n);
650  int nCoeffs = traceExp->GetNcoeffs();
651  int geomId = traceExp->GetGeom()->GetGlobalID();
652 
653  for (i = 0; i < nCoeffs*nCoeffs; ++i)
654  {
655  uniIds[cnt++] = geomId * maxTraceSize + i + 1;
656  }
657  }
658 
659  // Assemble matrices across partitions.
660  Gs::gs_data *gsh = Gs::Init(uniIds, m_comm);
661  Gs::Gather(tmpStore, Gs::gs_add, gsh);
662 
663  // Set up diagonal block matrix
664  Array<OneD, unsigned int> n_blks(nTrace - nDir);
665  for (n = 0; n < nTrace - nDir; ++n)
666  {
667  n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
668  }
669 
671  ::AllocateSharedPtr(n_blks, n_blks, eDIAGONAL);
672 
673  for (n = 0; n < nTrace - nDir; ++n)
674  {
675  int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
677  ::AllocateSharedPtr(nCoeffs, nCoeffs);
678  NekDouble *off = &tmpStore[arrayOffsets[n+nDir]];
679 
680  for (i = 0; i < nCoeffs; ++i)
681  {
682  for (j = 0; j < nCoeffs; ++j)
683  {
684  (*tmp)(i,j) = off[i*nCoeffs + j];
685  }
686  }
687 
688  tmp->Invert();
689  m_blkMat->SetBlock(n, n, tmp);
690  }
691  }
boost::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:49
LibUtilities::CommSharedPtr m_comm
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:223
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:150
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
boost::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:74
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
boost::shared_ptr< AssemblyMap > m_locToGloMap
const boost::weak_ptr< GlobalLinSys > m_linsys
static PreconditionerSharedPtr Nektar::MultiRegions::PreconditionerBlock::create ( const boost::shared_ptr< GlobalLinSys > &  plinsys,
const boost::shared_ptr< AssemblyMap > &  pLocToGloMap 
)
inlinestatic

Creates an instance of this class.

Definition at line 56 of file PreconditionerBlock.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

60  {
62  ::AllocateSharedPtr(plinsys,pLocToGloMap);
63  p->InitObject();
64  return p;
65  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< Preconditioner > PreconditionerSharedPtr
Definition: GlobalLinSys.h:62
void Nektar::MultiRegions::PreconditionerBlock::v_BuildPreconditioner ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 83 of file PreconditionerBlock.cpp.

References BlockPreconditionerCG(), BlockPreconditionerHDG(), Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), and m_linsys.

84  {
85  GlobalLinSysKey key = m_linsys.lock()->GetKey();
86 
87  // Different setup for HDG and CG.
88  if (key.GetMatrixType() == StdRegions::eHybridDGHelmBndLam)
89  {
91  }
92  else
93  {
95  }
96  }
void BlockPreconditionerCG(void)
Construct a block preconditioner from for the continuous Galerkin system.
const boost::weak_ptr< GlobalLinSys > m_linsys
void Nektar::MultiRegions::PreconditionerBlock::v_DoPreconditioner ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput 
)
privatevirtual

Apply preconditioner to pInput and store the result in pOutput.

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 697 of file PreconditionerBlock.cpp.

References Nektar::eWrapper, and m_locToGloMap.

700  {
701  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
702  int nGlobal = m_locToGloMap->GetNumGlobalBndCoeffs();
703  int nNonDir = nGlobal-nDir;
704  DNekBlkMat &M = (*m_blkMat);
705  NekVector<NekDouble> r(nNonDir,pInput,eWrapper);
706  NekVector<NekDouble> z(nNonDir,pOutput,eWrapper);
707  z = M * r;
708  }
boost::shared_ptr< AssemblyMap > m_locToGloMap
void Nektar::MultiRegions::PreconditionerBlock::v_InitObject ( )
privatevirtual

Reimplemented from Nektar::MultiRegions::Preconditioner.

Definition at line 75 of file PreconditionerBlock.cpp.

References ASSERTL0, Nektar::MultiRegions::eIterativeStaticCond, Nektar::MultiRegions::ePETScStaticCond, and m_locToGloMap.

76  {
77  GlobalSysSolnType solvertype=m_locToGloMap->GetGlobalSysSolnType();
79  solvertype == MultiRegions::ePETScStaticCond,
80  "Solver type not valid");
81  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
boost::shared_ptr< AssemblyMap > m_locToGloMap

Member Data Documentation

string Nektar::MultiRegions::PreconditionerBlock::className
static
Initial value:
"Block",
"Block Diagonal Preconditioning")

Name of class.

Registers the class with the Factory.

Definition at line 68 of file PreconditionerBlock.h.

DNekBlkMatSharedPtr Nektar::MultiRegions::PreconditionerBlock::m_blkMat
protected

Definition at line 80 of file PreconditionerBlock.h.

Referenced by BlockPreconditionerCG(), and BlockPreconditionerHDG().

const boost::weak_ptr<GlobalLinSys> Nektar::MultiRegions::PreconditionerBlock::m_linsys
protected
boost::shared_ptr<AssemblyMap> Nektar::MultiRegions::PreconditionerBlock::m_locToGloMap
protected
PreconditionerType Nektar::MultiRegions::PreconditionerBlock::m_preconType
protected

Definition at line 79 of file PreconditionerBlock.h.