Nektar++
GlobalLinSysDirectFull.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: GlobalLinSysDirectFull.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: GlobalLinSysDirectFull definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <MultiRegions/ExpList.h>
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42 namespace MultiRegions
43 {
44 /**
45  * @class GlobalLinSysDirect
46  *
47  * Consider a linear system
48  * \f$\boldsymbol{M\hat{u}}_g=\boldsymbol{\hat{f}}\f$
49  * to be solved, where \f$\boldsymbol{M}\f$ is a matrix of type
50  * specified by \a mkey. This function assembles the global system
51  * matrix \f$\boldsymbol{M}\f$ out of the elemental submatrices
52  * \f$\boldsymbol{M}^e\f$. This is equivalent to:
53  * \f[ \boldsymbol{M}=\boldsymbol{\mathcal{A}}^T
54  * \underline{\boldsymbol{M}}^e\boldsymbol{\mathcal{A}}.\f]
55  * where the matrix \f$\boldsymbol{\mathcal{A}}\f$ is a sparse
56  * permutation matrix of size \f$N_{\mathrm{eof}}\times
57  * N_{\mathrm{dof}}\f$. However, due to the size and sparsity of the
58  * matrix \f$\boldsymbol{\mathcal{A}}\f$, it is more efficient to
59  * assemble the global matrix using the mapping array \a
60  * map\f$[e][i]\f$ contained in the input argument \a locToGloMap.
61  * The global assembly is then evaluated as:
62  * \f[ \boldsymbol{M}\left[\mathrm{\texttt{map}}[e][i]\right]
63  * \left[\mathrm{\texttt{map}}[e][j]\right]
64  * =\mathrm{\texttt{sign}}[e][i]\cdot
65  * \mathrm{\texttt{sign}}[e][j] \cdot\boldsymbol{M}^e[i][j]\f]
66  * where the values \a sign\f$[e][i]\f$ ensure the correct connectivity.
67  */
68 
69 /**
70  * Registers the class with the Factory.
71  */
72 string GlobalLinSysDirectFull::className =
74  "DirectFull", GlobalLinSysDirectFull::create, "Direct Full.");
75 
76 /// Constructor for full direct matrix solve.
77 GlobalLinSysDirectFull::GlobalLinSysDirectFull(
78  const GlobalLinSysKey &pLinSysKey, const std::weak_ptr<ExpList> &pExp,
79  const std::shared_ptr<AssemblyMap> &pLocToGloMap)
80  : GlobalLinSys(pLinSysKey, pExp, pLocToGloMap),
81  GlobalLinSysDirect(pLinSysKey, pExp, pLocToGloMap)
82 {
83 
85  "This routine should only be used when using a Full Direct"
86  " matrix solve");
87  ASSERTL1(pExp.lock()->GetComm()->GetSize() == 1,
88  "Direct full matrix solve can only be used in serial.");
89 
90  AssembleFullMatrix(pLocToGloMap);
91 }
92 
94 {
95 }
96 
97 /**
98  * Solve the linear system using a full global matrix system.
99  */
101  const Array<OneD, const NekDouble> &pLocInput,
102  Array<OneD, NekDouble> &pLocOutput,
103  const AssemblyMapSharedPtr &pLocToGloMap,
104  const Array<OneD, const NekDouble> &pDirForcing)
105 {
106  bool dirForcCalculated = (bool)pDirForcing.size();
107 
108  int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
109  int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
110  int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
111 
112  Array<OneD, NekDouble> tmp(nLocDofs);
113  Array<OneD, NekDouble> tmp1(nLocDofs);
114  Array<OneD, NekDouble> global(nGlobDofs, 0.0);
115 
116  if (nDirDofs)
117  {
118  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
119 
120  // calculate the dirichlet forcing
121  if (dirForcCalculated)
122  {
123  // assume pDirForcing is in local space
124  ASSERTL0(
125  pDirForcing.size() >= nLocDofs,
126  "DirForcing is not of sufficient size. Is it in local space?");
127  Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, tmp1, 1);
128  }
129  else
130  {
131 
132  // Calculate Dirichlet forcing and subtract it from the rhs
133  expList->GeneralMatrixOp(m_linSysKey, pLocOutput, tmp);
134 
135  // Iterate over all the elements computing Robin BCs where
136  // necessary
137  for (auto &r : m_robinBCInfo) // add robin mass matrix
138  {
140  Array<OneD, NekDouble> tmploc;
141 
142  int n = r.first;
143  int offset = expList->GetCoeff_Offset(n);
144 
145  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
146  // add local matrix contribution
147  for (rBC = r.second; rBC; rBC = rBC->next)
148  {
149  vExp->AddRobinTraceContribution(
150  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
151  pLocOutput + offset, tmploc = tmp + offset);
152  }
153  }
154  Vmath::Vsub(nLocDofs, pLocInput, 1, tmp, 1, tmp1, 1);
155  }
156 
157  pLocToGloMap->Assemble(tmp1, tmp);
158 
159  SolveLinearSystem(nGlobDofs, tmp, global, pLocToGloMap, nDirDofs);
160  pLocToGloMap->GlobalToLocal(global, tmp);
161 
162  // Add back initial condition
163  Vmath::Vadd(nLocDofs, tmp, 1, pLocOutput, 1, pLocOutput, 1);
164  }
165  else
166  {
167  pLocToGloMap->Assemble(pLocInput, tmp);
168  SolveLinearSystem(nGlobDofs, tmp, global, pLocToGloMap, nDirDofs);
169  pLocToGloMap->GlobalToLocal(global, pLocOutput);
170  }
171 }
172 
173 /**
174  * Assemble a full matrix from the block matrix stored in
175  * #m_blkMatrices and the given local to global mapping information.
176  * @param locToGloMap Local to global mapping information.
177  */
179  const AssemblyMapSharedPtr &pLocToGloMap)
180 {
181  int i, j, n, cnt, gid1, gid2;
182  NekDouble sign1, sign2, value;
183  int totDofs = pLocToGloMap->GetNumGlobalCoeffs();
184  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
185 
186  unsigned int rows = totDofs - NumDirBCs;
187  unsigned int cols = totDofs - NumDirBCs;
188  NekDouble zero = 0.0;
189 
190  DNekMatSharedPtr Gmat;
191  int bwidth = pLocToGloMap->GetFullSystemBandWidth();
192  MatrixStorage matStorage = eFULL;
193 
194  switch (m_linSysKey.GetMatrixType())
195  {
196  // case for all symmetric matices
197  case StdRegions::eMass:
201  {
202  if ((2 * (bwidth + 1)) < rows)
203  {
206  rows, cols, zero, matStorage, bwidth, bwidth);
207  }
208  else
209  {
210  matStorage = ePOSITIVE_DEFINITE_SYMMETRIC;
212  rows, cols, zero, matStorage);
213  }
214  break;
215  }
218  {
219  matStorage = eFULL;
220  Gmat = MemoryManager<DNekMat>::AllocateSharedPtr(rows, cols, zero,
221  matStorage);
222  break;
223  }
224  default:
225  {
226  NEKERROR(ErrorUtil::efatal, "Add MatrixType to switch "
227  "statement");
228  }
229  }
230 
231  // fill global matrix
232  DNekScalMatSharedPtr loc_mat;
233 
234  int loc_lda;
235  for (n = cnt = 0; n < m_expList.lock()->GetNumElmts(); ++n)
236  {
237  loc_mat = GetBlock(n);
238  loc_lda = loc_mat->GetRows();
239 
240  for (i = 0; i < loc_lda; ++i)
241  {
242  gid1 = pLocToGloMap->GetLocalToGlobalMap(cnt + i) - NumDirBCs;
243  sign1 = pLocToGloMap->GetLocalToGlobalSign(cnt + i);
244  if (gid1 >= 0)
245  {
246  for (j = 0; j < loc_lda; ++j)
247  {
248  gid2 =
249  pLocToGloMap->GetLocalToGlobalMap(cnt + j) - NumDirBCs;
250  sign2 = pLocToGloMap->GetLocalToGlobalSign(cnt + j);
251  if (gid2 >= 0)
252  {
253  // When global matrix is symmetric,
254  // only add the value for the upper
255  // triangular part in order to avoid
256  // entries to be entered twice
257  if ((matStorage == eFULL) || (gid2 >= gid1))
258  {
259  value = Gmat->GetValue(gid1, gid2) +
260  sign1 * sign2 * (*loc_mat)(i, j);
261  Gmat->SetValue(gid1, gid2, value);
262  }
263  }
264  }
265  }
266  }
267  cnt += loc_lda;
268  }
269 
270  if (rows)
271  {
274  }
275 }
276 } // namespace MultiRegions
277 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual void v_Solve(const Array< OneD, const NekDouble > &pLocInput, Array< OneD, NekDouble > &pLocalOutput, 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.
void AssembleFullMatrix(const std::shared_ptr< AssemblyMap > &locToGloMap)
DNekLinSysSharedPtr m_linSys
Basic linear system object.
A global linear system.
Definition: GlobalLinSys.h:72
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:124
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:126
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
DNekScalMatSharedPtr GetBlock(unsigned int n)
Definition: GlobalLinSys.h:211
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:51
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
@ ePOSITIVE_DEFINITE_SYMMETRIC_BANDED
@ ePOSITIVE_DEFINITE_SYMMETRIC
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
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:359
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:419