Nektar++
GlobalLinSysDirectFull.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File GlobalLinSys.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: GlobalLinSys definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
36 #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",
75  GlobalLinSysDirectFull::create,
76  "Direct Full.");
77 
78 
79  /// Constructor for full direct matrix solve.
80  GlobalLinSysDirectFull::GlobalLinSysDirectFull(
81  const GlobalLinSysKey &pLinSysKey,
82  const std::weak_ptr<ExpList> &pExp,
83  const std::shared_ptr<AssemblyMap> &pLocToGloMap)
84  : GlobalLinSys(pLinSysKey, pExp, pLocToGloMap),
85  GlobalLinSysDirect(pLinSysKey, pExp, pLocToGloMap)
86  {
87 
89  "This routine should only be used when using a Full Direct"
90  " matrix solve");
91  ASSERTL1(pExp.lock()->GetComm()->GetSize() == 1,
92  "Direct full matrix solve can only be used in serial.");
93 
94  AssembleFullMatrix(pLocToGloMap);
95  }
96 
97 
99  {
100 
101  }
102 
103 
104  /**
105  * Solve the linear system using a full global matrix system.
106  */
108  const Array<OneD, const NekDouble> &pLocInput,
109  Array<OneD, NekDouble> &pLocOutput,
110  const AssemblyMapSharedPtr &pLocToGloMap,
111  const Array<OneD, const NekDouble> &pDirForcing)
112  {
113  bool dirForcCalculated = (bool) pDirForcing.size();
114 
115  int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
116  int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
117  int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
118 
119  Array<OneD, NekDouble> tmp (nLocDofs);
120  Array<OneD, NekDouble> tmp1(nLocDofs);
121  Array<OneD, NekDouble> global(nGlobDofs,0.0);
122 
123  if(nDirDofs)
124  {
125  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
126 
127  // calculate the dirichlet forcing
128  if(dirForcCalculated)
129  {
130  // assume pDirForcing is in local space
131  ASSERTL0(pDirForcing.size() >= nLocDofs,
132  "DirForcing is not of sufficient size. Is it in local space?");
133  Vmath::Vsub(nLocDofs, pLocInput, 1,
134  pDirForcing, 1, tmp1, 1);
135  }
136  else
137  {
138 
139  // Calculate Dirichlet forcing and subtract it from the rhs
140  expList->GeneralMatrixOp(m_linSysKey, pLocOutput, tmp);
141 
142  // Iterate over all the elements computing Robin BCs where
143  // necessary
144  for(auto &r : m_robinBCInfo) // add robin mass matrix
145  {
147  Array<OneD, NekDouble> tmploc;
148 
149  int n = r.first;
150  int offset = expList->GetCoeff_Offset(n);
151 
152  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
153  // add local matrix contribution
154  for(rBC = r.second;rBC; rBC = rBC->next)
155  {
156  vExp->AddRobinTraceContribution(rBC->m_robinID,
157  rBC->m_robinPrimitiveCoeffs,
158  pLocOutput + offset,
159  tmploc = tmp + offset);
160  }
161  }
162  Vmath::Vsub(nLocDofs, pLocInput, 1, tmp, 1, tmp1, 1);
163  }
164 
165  pLocToGloMap->Assemble(tmp1,tmp);
166 
167  SolveLinearSystem(nGlobDofs, tmp, global, pLocToGloMap, nDirDofs);
168  pLocToGloMap->GlobalToLocal(global,tmp);
169 
170  // Add back initial condition
171  Vmath::Vadd(nLocDofs, tmp, 1, pLocOutput, 1, pLocOutput, 1);
172  }
173  else
174  {
175  pLocToGloMap->Assemble(pLocInput,tmp);
176  SolveLinearSystem(nGlobDofs,tmp,global,pLocToGloMap,nDirDofs);
177  pLocToGloMap->GlobalToLocal(global,pLocOutput);
178  }
179  }
180 
181 
182  /**
183  * Assemble a full matrix from the block matrix stored in
184  * #m_blkMatrices and the given local to global mapping information.
185  * @param locToGloMap Local to global mapping information.
186  */
188  const AssemblyMapSharedPtr& pLocToGloMap)
189  {
190  int i,j,n,cnt,gid1,gid2;
191  NekDouble sign1,sign2,value;
192  int totDofs = pLocToGloMap->GetNumGlobalCoeffs();
193  int NumDirBCs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
194 
195  unsigned int rows = totDofs - NumDirBCs;
196  unsigned int cols = totDofs - NumDirBCs;
197  NekDouble zero = 0.0;
198 
199  DNekMatSharedPtr Gmat;
200  int bwidth = pLocToGloMap->GetFullSystemBandWidth();
201  MatrixStorage matStorage = eFULL;
202 
203  switch(m_linSysKey.GetMatrixType())
204  {
205  // case for all symmetric matices
206  case StdRegions::eMass:
210  {
211  if( (2*(bwidth+1)) < rows)
212  {
215  ::AllocateSharedPtr(rows, cols, zero,
216  matStorage,
217  bwidth, bwidth);
218  }
219  else
220  {
221  matStorage = ePOSITIVE_DEFINITE_SYMMETRIC;
223  ::AllocateSharedPtr(rows, cols, zero,
224  matStorage);
225  }
226  break;
227  }
230  {
231  matStorage = eFULL;
233  ::AllocateSharedPtr(rows, cols, zero,
234  matStorage);
235  break;
236  }
237  default:
238  {
239  NEKERROR(ErrorUtil::efatal, "Add MatrixType to switch "
240  "statement");
241  }
242  }
243 
244  // fill global matrix
245  DNekScalMatSharedPtr loc_mat;
246 
247  int loc_lda;
248  for(n = cnt = 0; n < m_expList.lock()->GetNumElmts(); ++n)
249  {
250  loc_mat = GetBlock(n);
251  loc_lda = loc_mat->GetRows();
252 
253  for(i = 0; i < loc_lda; ++i)
254  {
255  gid1 = pLocToGloMap->GetLocalToGlobalMap(cnt + i)-NumDirBCs;
256  sign1 = pLocToGloMap->GetLocalToGlobalSign(cnt + i);
257  if(gid1 >= 0)
258  {
259  for(j = 0; j < loc_lda; ++j)
260  {
261  gid2 = pLocToGloMap->GetLocalToGlobalMap(cnt + j)
262  - NumDirBCs;
263  sign2 = pLocToGloMap->GetLocalToGlobalSign(cnt + j);
264  if(gid2 >= 0)
265  {
266  // When global matrix is symmetric,
267  // only add the value for the upper
268  // triangular part in order to avoid
269  // entries to be entered twice
270  if((matStorage == eFULL)||(gid2 >= gid1))
271  {
272  value = Gmat->GetValue(gid1,gid2)
273  + sign1*sign2*(*loc_mat)(i,j);
274  Gmat->SetValue(gid1,gid2,value);
275  }
276  }
277  }
278  }
279  }
280  cnt += loc_lda;
281  }
282 
283  if(rows)
284  {
287  }
288  }
289  }
290 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#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:250
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void AssembleFullMatrix(const std::shared_ptr< AssemblyMap > &locToGloMap)
virtual void v_Solve(const Array< OneD, const NekDouble > &pLocInput, Array< OneD, NekDouble > &pLocalOutput, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve the linear system for given input and output vectors using a specified local to global map.
DNekLinSysSharedPtr m_linSys
Basic linear system object.
A global linear system.
Definition: GlobalLinSys.h:73
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:127
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:129
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:201
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:125
DNekScalMatSharedPtr GetBlock(unsigned int n)
Definition: GlobalLinSys.h:222
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:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
@ ePOSITIVE_DEFINITE_SYMMETRIC_BANDED
@ ePOSITIVE_DEFINITE_SYMMETRIC
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
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:322
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:372