Nektar++
GlobalLinSysIterativeFull.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: GlobalLinSysIterativeFull.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: GlobalLinSysIterativeFull definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 #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", GlobalLinSysIterativeFull::create,
59  "Iterative solver for full matrix system.");
60 
61 /**
62  * Constructor for full direct matrix solve.
63  * @param pKey Key specifying matrix to solve.
64  * @param pExp Shared pointer to expansion list for applying
65  * matrix evaluations.
66  * @param pLocToGloMap Local to global mapping.
67  */
68 GlobalLinSysIterativeFull::GlobalLinSysIterativeFull(
69  const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExp,
70  const std::shared_ptr<AssemblyMap> &pLocToGloMap)
71  : GlobalLinSys(pKey, pExp, pLocToGloMap),
72  GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
73 {
75  "This routine should only be used when using an Iterative "
76  "conjugate gradient matrix solve.");
77 }
78 
79 /**
80  *
81  */
83 {
84 }
85 
86 /**
87  * Solve a global linear system with Dirichlet forcing using a
88  * conjugate gradient method. This routine performs handling of the
89  * Dirichlet forcing terms and wraps the underlying iterative solver
90  * used for the remaining degrees of freedom.
91  *
92  * Consider solving for \f$x\f$, the matrix system \f$Ax=b\f$, where
93  * \f$b\f$ is known. To enforce the Dirichlet terms we instead solve
94  * \f[A(x-x_0) = b - Ax_0 \f]
95  * where \f$x_0\f$ is the Dirichlet forcing.
96  *
97  * @param pInput RHS of linear system, \f$b\f$.
98  * @param pOutput On input, values of dirichlet degrees
99  * of freedom with initial guess on other values.
100  * On output, the solution \f$x\f$.
101  * @param pLocToGloMap Local to global mapping.
102  * @param pDirForcing Precalculated Dirichlet forcing.
103  */
105  const Array<OneD, const NekDouble> &pLocInput,
106  Array<OneD, NekDouble> &pLocOutput,
107  const AssemblyMapSharedPtr &pLocToGloMap,
108  const Array<OneD, const NekDouble> &pDirForcing)
109 {
110  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
111  bool vCG = false;
112  m_locToGloMap = pLocToGloMap;
113 
114  if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
115  {
116  vCG = true;
117  }
118  else if (std::dynamic_pointer_cast<AssemblyMapDG>(pLocToGloMap))
119  {
120  vCG = false;
121  }
122  else
123  {
124  NEKERROR(ErrorUtil::efatal, "Unknown map type");
125  }
126 
127  bool dirForcCalculated = (bool)pDirForcing.size();
128  int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
129  int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
130  int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
131 
132  int nDirTotal = nDirDofs;
133  expList->GetComm()->GetRowComm()->AllReduce(nDirTotal,
135 
136  Array<OneD, NekDouble> tmp(nLocDofs);
137  Array<OneD, NekDouble> tmp1(nLocDofs);
138 
139  Array<OneD, NekDouble> global(nGlobDofs, 0.0);
140 
141  if (nDirTotal)
142  {
143  // calculate the Dirichlet forcing
144  if (dirForcCalculated)
145  {
146  Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, tmp1, 1);
147  }
148  else
149  {
150  // Calculate the dirichlet forcing B_b (== X_b) and
151  // substract it from the rhs
152  expList->GeneralMatrixOp(m_linSysKey, pLocOutput, tmp);
153 
154  // Iterate over all the elements computing Robin BCs where
155  // necessary
156  for (auto &r : m_robinBCInfo) // add robin mass matrix
157  {
159  Array<OneD, NekDouble> tmploc;
160 
161  int n = r.first;
162  int offset = expList->GetCoeff_Offset(n);
163 
164  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
165  // add local matrix contribution
166  for (rBC = r.second; rBC; rBC = rBC->next)
167  {
168  vExp->AddRobinTraceContribution(
169  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
170  pLocOutput + offset, tmploc = tmp + offset);
171  }
172  }
173 
174  Vmath::Vsub(nLocDofs, pLocInput, 1, tmp, 1, tmp1, 1);
175  }
176  if (vCG)
177  {
178  pLocToGloMap->Assemble(tmp1, tmp);
179 
180  // solve for perturbation from initial guess in pOutput
181  SolveLinearSystem(nGlobDofs, tmp, global, pLocToGloMap, nDirDofs);
182 
183  pLocToGloMap->GlobalToLocal(global, tmp);
184 
185  // Add back initial condition
186  Vmath::Vadd(nLocDofs, tmp, 1, pLocOutput, 1, pLocOutput, 1);
187  }
188  else
189  {
190  ASSERTL0(false, "Need DG solve if using Dir BCs");
191  }
192  }
193  else
194  {
195  pLocToGloMap->Assemble(pLocInput, tmp);
196  SolveLinearSystem(nGlobDofs, tmp, global, pLocToGloMap, nDirDofs);
197  pLocToGloMap->GlobalToLocal(global, pLocOutput);
198  }
199 }
200 
201 /**
202  *
203  */
205  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
206 {
207  std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
208 
209  AssemblyMapSharedPtr asmMap = m_locToGloMap.lock();
210 
211  int ncoeffs = expList->GetNcoeffs();
212 
213  Array<OneD, NekDouble> InputLoc(ncoeffs);
214  Array<OneD, NekDouble> OutputLoc(ncoeffs);
215  asmMap->GlobalToLocal(pInput, InputLoc);
216 
217  // Perform matrix-vector operation A*d_i
218  expList->GeneralMatrixOp(m_linSysKey, InputLoc, OutputLoc);
219 
220  // Apply robin boundary conditions to the solution.
221  for (auto &r : m_robinBCInfo) // add robin mass matrix
222  {
225 
226  int n = r.first;
227 
228  int offset = expList->GetCoeff_Offset(n);
229  LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
230 
231  // add local matrix contribution
232  for (rBC = r.second; rBC; rBC = rBC->next)
233  {
234  vExp->AddRobinTraceContribution(
235  rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, InputLoc + offset,
236  tmp = OutputLoc + offset);
237  }
238  }
239 
240  // put back in global coeffs
241  asmMap->Assemble(OutputLoc, pOutput);
242 }
243 
244 /**
245  *
246  */
248 {
249  m_map = m_locToGloMap.lock()->GetGlobalToUniversalMapUnique();
250 }
251 
252 } // namespace MultiRegions
253 } // 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
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
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, 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.
Array< OneD, int > m_map
Global to universal unique map.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution 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
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