Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // 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: GlobalLinSys definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <map>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44  namespace MultiRegions
45  {
46  /**
47  * @class GlobalLinSysIterativeCG
48  *
49  * This class implements a conjugate gradient matrix solver.
50  * Preconditioning is implemented using a Jacobi (diagonal)
51  * preconditioner.
52  */
53 
54  /**
55  * Registers the class with the Factory.
56  */
57  string GlobalLinSysIterativeFull::className
59  "IterativeFull",
60  GlobalLinSysIterativeFull::create,
61  "Iterative solver for full matrix system.");
62 
63 
64  /**
65  * Constructor for full direct matrix solve.
66  * @param pKey Key specifying matrix to solve.
67  * @param pExp Shared pointer to expansion list for applying
68  * matrix evaluations.
69  * @param pLocToGloMap Local to global mapping.
70  */
71  GlobalLinSysIterativeFull::GlobalLinSysIterativeFull(
72  const GlobalLinSysKey &pKey,
73  const boost::weak_ptr<ExpList> &pExp,
74  const boost::shared_ptr<AssemblyMap> &pLocToGloMap)
75  : GlobalLinSys (pKey, pExp, pLocToGloMap),
76  GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
77  {
79  "This routine should only be used when using an Iterative "
80  "conjugate gradient matrix solve.");
81  }
82 
83 
84  /**
85  *
86  */
88  {
89 
90  }
91 
92 
93  /**
94  * Solve a global linear system with Dirichlet forcing using a
95  * conjugate gradient method. This routine performs handling of the
96  * Dirichlet forcing terms and wraps the underlying iterative solver
97  * used for the remaining degrees of freedom.
98  *
99  * Consider solving for \f$x\f$, the matrix system \f$Ax=b\f$, where
100  * \f$b\f$ is known. To enforce the Dirichlet terms we instead solve
101  * \f[A(x-x_0) = b - Ax_0 \f]
102  * where \f$x_0\f$ is the Dirichlet forcing.
103  *
104  * @param pInput RHS of linear system, \f$b\f$.
105  * @param pOutput On input, values of dirichlet degrees
106  * of freedom with initial guess on other values.
107  * On output, the solution \f$x\f$.
108  * @param pLocToGloMap Local to global mapping.
109  * @param pDirForcing Precalculated Dirichlet forcing.
110  */
112  const Array<OneD, const NekDouble> &pInput,
113  Array<OneD, NekDouble> &pOutput,
114  const AssemblyMapSharedPtr &pLocToGloMap,
115  const Array<OneD, const NekDouble> &pDirForcing)
116  {
117  boost::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
118  bool vCG;
119  if ((m_locToGloMap = boost::dynamic_pointer_cast<AssemblyMapCG>(
120  pLocToGloMap)))
121  {
122  vCG = true;
123  }
124  else if ((m_locToGloMap = boost::dynamic_pointer_cast<
125  AssemblyMapDG>(pLocToGloMap)))
126  {
127  vCG = false;
128  }
129  else
130  {
131  ASSERTL0(false, "Unknown map type");
132  }
133 
134  bool dirForcCalculated = (bool) pDirForcing.num_elements();
135  int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
136  int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
137  int nDirTotal = nDirDofs;
138 
139  expList->GetComm()->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  boost::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  // retrieve robin boundary condition information and apply robin
199  // boundary conditions to the solution.
200  const std::map<int, RobinBCInfoSharedPtr> vRobinBCInfo
201  = expList->GetRobinBCInfo();
202  if(vRobinBCInfo.size() > 0)
203  {
204  ASSERTL0(false,
205  "Robin boundaries not set up in IterativeFull solver.");
206  int nGlobal = m_locToGloMap->GetNumGlobalCoeffs();
207  int nLocal = m_locToGloMap->GetNumLocalCoeffs();
208  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
209  int nNonDir = nGlobal - nDir;
210  Array<OneD, NekDouble> robin_A(nGlobal, 0.0);
211  Array<OneD, NekDouble> robin_l(nLocal, 0.0);
213  NekVector<NekDouble> robin(nNonDir,
214  tmp = robin_A + nDir, eWrapper);
215 
216  // Operation: p_A = A * d_A
217  // First map d_A to local solution
218  m_locToGloMap->GlobalToLocal(pInput, robin_l);
219 
220  // Iterate over all the elements computing Robin BCs where
221  // necessary
222  for (int n = 0; n < expList->GetNumElmts(); ++n)
223  {
224  int nel = expList->GetOffset_Elmt_Id(n);
225  int offset = expList->GetCoeff_Offset(n);
226  int ncoeffs = expList->GetExp(nel)->GetNcoeffs();
227 
228  if(vRobinBCInfo.count(nel) != 0) // add robin mass matrix
229  {
232  StdRegions::StdExpansionSharedPtr vExp = expList->GetExp(nel);
233 
234  // add local matrix contribution
235  for(rBC = vRobinBCInfo.find(nel)->second;rBC; rBC = rBC->next)
236  {
237  vExp->AddRobinEdgeContribution(rBC->m_robinID,rBC->m_robinPrimitiveCoeffs, tmp = robin_l + offset);
238  }
239  }
240  else
241  {
242  Vmath::Zero(ncoeffs, &robin_l[offset], 1);
243  }
244  }
245 
246  // Map local Robin contribution back to global coefficients
247  m_locToGloMap->LocalToGlobal(robin_l, robin_A);
248  // Add them to the output of the GeneralMatrixOp
249  Vmath::Vadd(nGlobal, pOutput, 1, robin_A, 1, pOutput, 1);
250  }
251 
252  }
253 
254  /**
255  *
256  */
258  {
259  m_map = m_locToGloMap->GetGlobalToUniversalMapUnique();
260  }
261 
262  }
263 }
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
boost::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
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:199
STL namespace.
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:127
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:329
A global linear system.
Definition: GlobalLinSys.h:74
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:359
GlobalLinSysFactory & GetGlobalLinSysFactory()
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
const boost::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:129