Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GlobalLinSysIterativeStaticCond.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: GlobalLinSysIterativeStaticCond.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Implementation to linear solver using single-
33 // or multi-level static condensation
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
43 
44 namespace Nektar
45 {
46  namespace MultiRegions
47  {
48  /**
49  * @class GlobalLinSysIterativeStaticCond
50  *
51  * Solves a linear system iteratively using single- or multi-level
52  * static condensation.
53  */
54 
55  /**
56  * Registers the class with the Factory.
57  */
60  "IterativeStaticCond",
62  "Iterative static condensation.");
63 
66  "IterativeMultiLevelStaticCond",
68  "Iterative multi-level static condensation.");
69 
70 
73  "LocalMatrixStorageStrategy",
74  "Sparse");
77  "LocalMatrixStorageStrategy",
78  "Contiguous",
81  "LocalMatrixStorageStrategy",
82  "Non-contiguous",
85  "LocalMatrixStorageStrategy",
86  "Sparse",
88  };
89 
90  /**
91  * For a matrix system of the form @f[
92  * \left[ \begin{array}{cc}
93  * \boldsymbol{A} & \boldsymbol{B}\\
94  * \boldsymbol{C} & \boldsymbol{D}
95  * \end{array} \right]
96  * \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2}
97  * \end{array}\right]
98  * = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2}
99  * \end{array}\right],
100  * @f]
101  * where @f$\boldsymbol{D}@f$ and
102  * @f$(\boldsymbol{A-BD^{-1}C})@f$ are invertible, store and assemble
103  * a static condensation system, according to a given local to global
104  * mapping. #m_linSys is constructed by AssembleSchurComplement().
105  * @param mKey Associated matrix key.
106  * @param pLocMatSys LocalMatrixSystem
107  * @param locToGloMap Local to global mapping.
108  */
110  const GlobalLinSysKey &pKey,
111  const boost::weak_ptr<ExpList> &pExpList,
112  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
113  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
114  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
115  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
116  {
119  "This constructor is only valid when using static "
120  "condensation");
122  == pLocToGloMap->GetGlobalSysSolnType(),
123  "The local to global map is not set up for the requested "
124  "solution type");
125  }
126 
127 
128  /**
129  *
130  */
132  const GlobalLinSysKey &pKey,
133  const boost::weak_ptr<ExpList> &pExpList,
134  const DNekScalBlkMatSharedPtr pSchurCompl,
135  const DNekScalBlkMatSharedPtr pBinvD,
136  const DNekScalBlkMatSharedPtr pC,
137  const DNekScalBlkMatSharedPtr pInvD,
138  const boost::shared_ptr<AssemblyMap> &pLocToGloMap,
139  const PreconditionerSharedPtr pPrecon)
140  : GlobalLinSys (pKey, pExpList, pLocToGloMap),
141  GlobalLinSysIterative (pKey, pExpList, pLocToGloMap),
142  GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
143  {
144  m_schurCompl = pSchurCompl;
145  m_BinvD = pBinvD;
146  m_C = pC;
147  m_invD = pInvD;
148  m_precon = pPrecon;
149  }
150 
151 
153  {
155  = m_locToGloMap->GetPreconType();
156  std::string PreconType
159  PreconType,GetSharedThisPtr(),m_locToGloMap);
160 
161  // Allocate memory for top-level structure
163 
164  // Setup Block Matrix systems
165  int n, n_exp = m_expList.lock()->GetNumElmts();
166 
167  MatrixStorage blkmatStorage = eDIAGONAL;
168  const Array<OneD,const unsigned int>& nbdry_size
169  = m_locToGloMap->GetNumLocalBndCoeffsPerPatch();
170 
172  ::AllocateSharedPtr(nbdry_size, nbdry_size , blkmatStorage);
173 
174  // Preserve original matrix in m_S1Blk
175  for (n = 0; n < n_exp; ++n)
176  {
177  DNekScalMatSharedPtr mat = m_schurCompl->GetBlock(n, n);
178  m_S1Blk->SetBlock(n, n, mat);
179  }
180 
181  // Build preconditioner
182  m_precon->BuildPreconditioner();
183 
184  // Do transform of Schur complement matrix
185  for (n = 0; n < n_exp; ++n)
186  {
187  if (m_linSysKey.GetMatrixType() !=
189  {
190  DNekScalMatSharedPtr mat = m_S1Blk->GetBlock(n, n);
191  DNekScalMatSharedPtr t = m_precon->TransformedSchurCompl(
192  m_expList.lock()->GetOffset_Elmt_Id(n), mat);
193  m_schurCompl->SetBlock(n, n, t);
194  }
195  }
196 
197  // Construct this level
199  }
200 
201  /**
202  *
203  */
205  {
206 
207  }
208 
210  v_GetStaticCondBlock(unsigned int n)
211  {
212  DNekScalBlkMatSharedPtr schurComplBlock;
213  int scLevel = m_locToGloMap->GetStaticCondLevel();
214  DNekScalBlkMatSharedPtr sc = scLevel == 0 ? m_S1Blk : m_schurCompl;
215  DNekScalMatSharedPtr localMat = sc->GetBlock(n,n);
216  unsigned int nbdry = localMat->GetRows();
217  unsigned int nblks = 1;
218  unsigned int esize[1] = {nbdry};
219 
220  schurComplBlock = MemoryManager<DNekScalBlkMat>
221  ::AllocateSharedPtr(nblks, nblks, esize, esize);
222  schurComplBlock->SetBlock(0, 0, localMat);
223 
224  return schurComplBlock;
225  }
226 
227  /**
228  * Assemble the schur complement matrix from the block matrices stored
229  * in #m_blkMatrices and the given local to global mapping information.
230  * @param locToGloMap Local to global mapping information.
231  */
233  const AssemblyMapSharedPtr pLocToGloMap)
234  {
235  int i,j,n,cnt,gid1,gid2;
236  NekDouble sign1,sign2;
237 
238  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
239  DoGlobalMatOp(m_linSysKey.GetMatrixType());
240 
241  // Set up unique map
242  v_UniqueMap();
243 
244  // Build precon again if we in multi-level static condensation (a
245  // bit of a hack)
247  {
249  = m_locToGloMap->GetPreconType();
250  std::string PreconType
253  PreconType,GetSharedThisPtr(),m_locToGloMap);
254  m_precon->BuildPreconditioner();
255  }
256 
257  if (!doGlobalOp)
258  {
260  return;
261  }
262 
263  int nBndDofs = pLocToGloMap->GetNumGlobalBndCoeffs();
264  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
265  unsigned int rows = nBndDofs - NumDirBCs;
266  unsigned int cols = nBndDofs - NumDirBCs;
267 
268  // COO sparse storage to assist in assembly
269  COOMatType gmat_coo;
270 
271  // Get the matrix storage structure
272  // (whether to store only one triangular part, if symmetric)
273  MatrixStorage matStorage = eFULL;
274 
275  // assemble globally
276  DNekScalMatSharedPtr loc_mat;
277  int loc_lda;
278  for(n = cnt = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
279  {
280  loc_mat = m_schurCompl->GetBlock(n,n);
281  loc_lda = loc_mat->GetRows();
282 
283  // Set up Matrix;
284  for(i = 0; i < loc_lda; ++i)
285  {
286  gid1 = pLocToGloMap->GetLocalToGlobalBndMap (cnt + i)
287  - NumDirBCs;
288  sign1 = pLocToGloMap->GetLocalToGlobalBndSign(cnt + i);
289 
290  if(gid1 >= 0)
291  {
292  for(j = 0; j < loc_lda; ++j)
293  {
294  gid2 = pLocToGloMap->GetLocalToGlobalBndMap(cnt+j)
295  - NumDirBCs;
296  sign2 = pLocToGloMap->GetLocalToGlobalBndSign(cnt+j);
297 
298  if (gid2 >= 0)
299  {
300  gmat_coo[std::make_pair(gid1,gid2)] +=
301  sign1*sign2*(*loc_mat)(i,j);
302  }
303  }
304  }
305  }
306  cnt += loc_lda;
307  }
308 
310  sparseStorage (1);
311 
312  BCOMatType partMat;
313  convertCooToBco(rows, cols, 1, gmat_coo, partMat);
314 
315  sparseStorage[0] =
317  AllocateSharedPtr(rows, cols, 1, partMat, matStorage );
318 
319  // Create block diagonal matrix
321  AllocateSharedPtr(sparseStorage);
322  }
323 
324 
325  /**
326  * Populates sparse block-diagonal schur complement matrix from
327  * the block matrices stored in #m_blkMatrices.
328  */
330  {
331  LocalMatrixStorageStrategy storageStrategy =
332  m_expList.lock()->GetSession()->
333  GetSolverInfoAsEnum<LocalMatrixStorageStrategy>(
334  "LocalMatrixStorageStrategy");
335 
336  switch(storageStrategy)
337  {
340  {
341  size_t storageSize = 0;
342  int nBlk = m_schurCompl->GetNumberOfBlockRows();
343 
344  m_scale = Array<OneD, NekDouble> (nBlk, 1.0);
345  m_rows = Array<OneD, unsigned int> (nBlk, 0U);
346 
347  // Determine storage requirements for dense blocks.
348  for (int i = 0; i < nBlk; ++i)
349  {
350  m_rows[i] = m_schurCompl->GetBlock(i,i)->GetRows();
351  m_scale[i] = m_schurCompl->GetBlock(i,i)->Scale();
352  storageSize += m_rows[i] * m_rows[i];
353  }
354 
355  // Assemble dense storage blocks.
356  DNekScalMatSharedPtr loc_mat;
357  m_denseBlocks.resize(nBlk);
358  double *ptr = 0;
359 
360  if (MultiRegions::eContiguous == storageStrategy)
361  {
362  m_storage.resize (storageSize);
363  ptr = &m_storage[0];
364  }
365 
366  for (unsigned int n = 0; n < nBlk; ++n)
367  {
368  loc_mat = m_schurCompl->GetBlock(n,n);
369 
370  if (MultiRegions::eContiguous == storageStrategy)
371  {
372  int loc_lda = loc_mat->GetRows();
373  int blockSize = loc_lda * loc_lda;
374  m_denseBlocks[n] = ptr;
375  for(int i = 0; i < loc_lda; ++i)
376  {
377  for(int j = 0; j < loc_lda; ++j)
378  {
379  ptr[j*loc_lda+i] = (*loc_mat)(i,j);
380  }
381  }
382  ptr += blockSize;
384  m_expList.lock()->GetOffset_Elmt_Id(n));
385  }
386  else
387  {
388  m_denseBlocks[n] = loc_mat->GetRawPtr();
389  }
390  }
391  break;
392  }
394  {
395  DNekScalMatSharedPtr loc_mat;
396  int loc_lda;
397  int blockSize = 0;
398 
399  // First run through to split the set of local matrices into
400  // partitions of fixed block size, and count number of local
401  // matrices that belong to each partition.
402  std::vector<std::pair<int,int> > partitions;
403  for(int n = 0; n < m_schurCompl->GetNumberOfBlockRows(); ++n)
404  {
405  loc_mat = m_schurCompl->GetBlock(n,n);
406  loc_lda = loc_mat->GetRows();
407 
408  ASSERTL1(loc_lda >= 0,
409  boost::lexical_cast<std::string>(n) + "-th "
410  "matrix block in Schur complement has "
411  "rank 0!");
412 
413  if (blockSize == loc_lda)
414  {
415  partitions[partitions.size()-1].first++;
416  }
417  else
418  {
419  blockSize = loc_lda;
420  partitions.push_back(make_pair(1,loc_lda));
421  }
422  }
423 
424  MatrixStorage matStorage = eFULL;
425 
426  // Create a vector of sparse storage holders
428  sparseStorage (partitions.size());
429 
430  for (int part = 0, n = 0; part < partitions.size(); ++part)
431  {
432  BCOMatType partMat;
433 
434  for(int k = 0; k < partitions[part].first; ++k, ++n)
435  {
436  loc_mat = m_schurCompl->GetBlock(n,n);
437  loc_lda = loc_mat->GetRows();
438 
439  ASSERTL1(loc_lda == partitions[part].second,
440  boost::lexical_cast<std::string>(n) + "-th"
441  " matrix block in Schur complement has "
442  "unexpected rank");
443 
444  partMat[make_pair(k,k)] = BCOEntryType(
445  loc_lda*loc_lda, loc_mat->GetRawPtr());
446 
448  m_expList.lock()->GetOffset_Elmt_Id(n));
449  }
450 
451  sparseStorage[part] =
454  partitions[part].first, partitions[part].first,
455  partitions[part].second, partMat, matStorage );
456  }
457 
458  // Create block diagonal matrix
460  AllocateSharedPtr(sparseStorage);
461 
462  break;
463  }
464  default:
465  ErrorUtil::NekError("Solver info property \
466  LocalMatrixStorageStrategy takes values \
467  Contiguous, Non-contiguous and Sparse");
468  }
469  }
470 
471  /**
472  *
473  */
475  const Array<OneD, NekDouble>& pInput,
476  Array<OneD, NekDouble>& pOutput)
477  {
478  int nLocal = m_locToGloMap->GetNumLocalBndCoeffs();
479  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
480  bool doGlobalOp = m_expList.lock()->GetGlobalOptParam()->
481  DoGlobalMatOp(m_linSysKey.GetMatrixType());
482 
483  if(doGlobalOp)
484  {
485  // Do matrix multiply globally
486  Array<OneD, NekDouble> in = pInput + nDir;
487  Array<OneD, NekDouble> out = pOutput + nDir;
488 
489  m_sparseSchurCompl->Multiply(in,out);
490  m_locToGloMap->UniversalAssembleBnd(pOutput, nDir);
491  }
492  else if (m_sparseSchurCompl)
493  {
494  // Do matrix multiply locally using block-diagonal sparse matrix
495  Array<OneD, NekDouble> tmp = m_wsp + nLocal;
496 
497  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
498  m_sparseSchurCompl->Multiply(m_wsp,tmp);
499  m_locToGloMap->AssembleBnd(tmp, pOutput);
500  }
501  else
502  {
503  // Do matrix multiply locally, using direct BLAS calls
504  m_locToGloMap->GlobalToLocalBnd(pInput, m_wsp);
505  int i, cnt;
506  Array<OneD, NekDouble> tmpout = m_wsp + nLocal;
507  for (i = cnt = 0; i < m_denseBlocks.size(); cnt += m_rows[i], ++i)
508  {
509  const int rows = m_rows[i];
510  Blas::Dgemv('N', rows, rows,
511  m_scale[i], m_denseBlocks[i], rows,
512  m_wsp.get()+cnt, 1,
513  0.0, tmpout.get()+cnt, 1);
514  }
515  m_locToGloMap->AssembleBnd(tmpout, pOutput);
516  }
517  }
518 
520  {
521  m_map = m_locToGloMap->GetGlobalToUniversalBndMapUnique();
522  }
523 
525  int scLevel,
526  NekVector<NekDouble> &F_GlobBnd)
527  {
528  if (scLevel == 0)
529  {
530  Set_Rhs_Magnitude(F_GlobBnd);
531  return m_S1Blk;
532  }
533  else
534  {
535  return m_schurCompl;
536  }
537  }
538 
540  Array<OneD, NekDouble>& pInOut,
541  int offset)
542  {
543  m_precon->DoTransformToLowEnergy(pInOut, offset);
544  }
545 
547  Array<OneD, NekDouble>& pInOut)
548  {
549  m_precon->DoTransformFromLowEnergy(pInOut);
550  }
551 
553  const GlobalLinSysKey &mkey,
554  const boost::weak_ptr<ExpList> &pExpList,
555  const DNekScalBlkMatSharedPtr pSchurCompl,
556  const DNekScalBlkMatSharedPtr pBinvD,
557  const DNekScalBlkMatSharedPtr pC,
558  const DNekScalBlkMatSharedPtr pInvD,
559  const boost::shared_ptr<AssemblyMap> &l2gMap)
560  {
562  GlobalLinSysIterativeStaticCond>::AllocateSharedPtr(
563  mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap,
564  m_precon);
565  sys->Initialise(l2gMap);
566  return sys;
567  }
568  }
569 }