Nektar++
GlobalLinSysIterativeFull.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 
35 #include <map>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  /**
46  * @class GlobalLinSysIterativeCG
47  *
48  * This class implements a conjugate gradient matrix solver.
49  * Preconditioning is implemented using a Jacobi (diagonal)
50  * preconditioner.
51  */
52 
53  /**
54  * Registers the class with the Factory.
55  */
56  string GlobalLinSysIterativeFull::className
58  "IterativeFull",
59  GlobalLinSysIterativeFull::create,
60  "Iterative solver for full matrix system.");
61 
62 
63  /**
64  * Constructor for full direct matrix solve.
65  * @param pKey Key specifying matrix to solve.
66  * @param pExp Shared pointer to expansion list for applying
67  * matrix evaluations.
68  * @param pLocToGloMap Local to global mapping.
69  */
70  GlobalLinSysIterativeFull::GlobalLinSysIterativeFull(
71  const GlobalLinSysKey &pKey,
72  const std::weak_ptr<ExpList> &pExp,
73  const std::shared_ptr<AssemblyMap> &pLocToGloMap)
74  : GlobalLinSys (pKey, pExp, pLocToGloMap),
75  GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
76  {
78  "This routine should only be used when using an Iterative "
79  "conjugate gradient matrix solve.");
80  }
81 
82 
83  /**
84  *
85  */
87  {
88 
89  }
90 
91 
92  /**
93  * Solve a global linear system with Dirichlet forcing using a
94  * conjugate gradient method. This routine performs handling of the
95  * Dirichlet forcing terms and wraps the underlying iterative solver
96  * used for the remaining degrees of freedom.
97  *
98  * Consider solving for \f$x\f$, the matrix system \f$Ax=b\f$, where
99  * \f$b\f$ is known. To enforce the Dirichlet terms we instead solve
100  * \f[A(x-x_0) = b - Ax_0 \f]
101  * where \f$x_0\f$ is the Dirichlet forcing.
102  *
103  * @param pInput RHS of linear system, \f$b\f$.
104  * @param pOutput On input, values of dirichlet degrees
105  * of freedom with initial guess on other values.
106  * On output, the solution \f$x\f$.
107  * @param pLocToGloMap Local to global mapping.
108  * @param pDirForcing Precalculated Dirichlet forcing.
109  */
111  const Array<OneD, const NekDouble> &pInput,
112  Array<OneD, NekDouble> &pOutput,
113  const AssemblyMapSharedPtr &pLocToGloMap,
114  const Array<OneD, const NekDouble> &pDirForcing)
115  {
116  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
117  bool vCG = false;
118  m_locToGloMap = pLocToGloMap;
119 
120  if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
121  {
122  vCG = true;
123  }
124  else if (std::dynamic_pointer_cast<AssemblyMapDG>(pLocToGloMap))
125  {
126  vCG = false;
127  }
128  else
129  {
130  NEKERROR(ErrorUtil::efatal, "Unknown map type");
131  }
132 
133  bool dirForcCalculated = (bool) pDirForcing.num_elements();
134  int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
135  int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
136  int nDirTotal = nDirDofs;
137 
138  expList->GetComm()->GetRowComm()
139  ->AllReduce(nDirTotal, LibUtilities::ReduceSum);
140 
141  Array<OneD, NekDouble> tmp(nGlobDofs), tmp2;
142 
143  if(nDirTotal)
144  {
145  // calculate the Dirichlet forcing
146  if(dirForcCalculated)
147  {
148  Vmath::Vsub(nGlobDofs, pInput.get(), 1,
149  pDirForcing.get(), 1,
150  tmp.get(), 1);
151  }
152  else
153  {
154  // Calculate the dirichlet forcing B_b (== X_b) and
155  // substract it from the rhs
156  expList->GeneralMatrixOp(
157  m_linSysKey, pOutput, tmp, eGlobal);
158 
159  Vmath::Vsub(nGlobDofs, pInput.get(), 1,
160  tmp.get(), 1,
161  tmp.get(), 1);
162  }
163  if (vCG)
164  {
165  Array<OneD, NekDouble> out(nGlobDofs,0.0);
166 
167  // solve for perturbation from intiial guess in pOutput
169  nGlobDofs, tmp, out, pLocToGloMap, nDirDofs);
170  Vmath::Vadd(nGlobDofs-nDirDofs, &out [nDirDofs], 1,
171  &pOutput[nDirDofs], 1, &pOutput[nDirDofs], 1);
172  }
173  else
174  {
175  ASSERTL0(false, "Need DG solve if using Dir BCs");
176  }
177  }
178  else
179  {
180  Vmath::Vcopy(nGlobDofs, pInput, 1, tmp, 1);
181  SolveLinearSystem(nGlobDofs, tmp, pOutput, pLocToGloMap);
182  }
183  }
184 
185 
186  /**
187  *
188  */
190  const Array<OneD, NekDouble>& pInput,
191  Array<OneD, NekDouble>& pOutput)
192  {
193  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
194  // Perform matrix-vector operation A*d_i
195  expList->GeneralMatrixOp(m_linSysKey,
196  pInput, pOutput, eGlobal);
197 
198  AssemblyMapSharedPtr asmMap = m_locToGloMap.lock();
199 
200  // Apply robin boundary conditions to the solution.
201  if(m_robinBCInfo.size() > 0)
202  {
203  ASSERTL0(false,
204  "Robin boundaries not set up in IterativeFull solver.");
205  int nGlobal = asmMap->GetNumGlobalCoeffs();
206  int nLocal = asmMap->GetNumLocalCoeffs();
207  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
208  int nNonDir = nGlobal - nDir;
209  Array<OneD, NekDouble> robin_A(nGlobal, 0.0);
210  Array<OneD, NekDouble> robin_l(nLocal, 0.0);
212  NekVector<NekDouble> robin(nNonDir,
213  tmp = robin_A + nDir, eWrapper);
214 
215  // Operation: p_A = A * d_A
216  // First map d_A to local solution
217  asmMap->GlobalToLocal(pInput, robin_l);
218 
219  // Iterate over all the elements computing Robin BCs where
220  // necessary
221  for (int n = 0; n < expList->GetNumElmts(); ++n)
222  {
223  int nel = n;
224  int offset = expList->GetCoeff_Offset(n);
225  int ncoeffs = expList->GetExp(nel)->GetNcoeffs();
226 
227  if(m_robinBCInfo.count(nel) != 0) // add robin mass matrix
228  {
231  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(nel);
232 
233  // add local matrix contribution
234  for(rBC = m_robinBCInfo.find(nel)->second;rBC; rBC = rBC->next)
235  {
236  vExp->AddRobinEdgeContribution(rBC->m_robinID,rBC->m_robinPrimitiveCoeffs, tmp = robin_l + offset);
237  }
238  }
239  else
240  {
241  Vmath::Zero(ncoeffs, &robin_l[offset], 1);
242  }
243  }
244 
245  // Map local Robin contribution back to global coefficients
246  asmMap->LocalToGlobal(robin_l, robin_A);
247  // Add them to the output of the GeneralMatrixOp
248  Vmath::Vadd(nGlobal, pOutput, 1, robin_A, 1, pOutput, 1);
249  }
250 
251  }
252 
253  /**
254  *
255  */
257  {
258  m_map = m_locToGloMap.lock()->GetGlobalToUniversalMapUnique();
259  }
260 
261  }
262 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
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
STL namespace.
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
Global coefficients.
Array< OneD, int > m_map
Global to universal unique map.
Describe a linear system.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:125
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:346
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
A global linear system.
Definition: GlobalLinSys.h:72
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:127
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, 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...
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
GlobalLinSysFactory & GetGlobalLinSysFactory()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:129
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
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)