Nektar++
PreconditionerLinear.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PreconditionerLinear.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: PreconditionerLinear definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
41 
42 #ifdef NEKTAR_USING_PETSC
44 #endif
45 
46 #include <LocalRegions/MatrixKey.h>
47 #include <cmath>
48 
49 using namespace std;
50 
51 namespace Nektar
52 {
53 namespace MultiRegions
54 {
55 /**
56  * Registers the class with the Factory.
57  */
58 
59 string PreconditionerLinear::className1 =
61  "FullLinearSpace", PreconditionerLinear::create,
62  "Full Linear space inverse Preconditioning");
63 
64 std::string PreconditionerLinear::solveType =
65  LibUtilities::SessionReader::RegisterDefaultSolverInfo("LinearPreconSolver",
66  "Xxt");
67 std::string PreconditionerLinear::solveTypeIds[] = {
68  LibUtilities::SessionReader::RegisterEnumValue(
69  "LinearPreconSolver", "PETSc", MultiRegions::eLinearPreconPETSc),
70  LibUtilities::SessionReader::RegisterEnumValue(
71  "LinearPreconSolver", "Xxt", MultiRegions::eLinearPreconXxt)};
72 
73 /**
74  * @class PreconditionerLinear
75  *
76  * This class implements preconditioning for the conjugate
77  * gradient matrix solver.
78  */
79 
80 PreconditionerLinear::PreconditionerLinear(
81  const std::shared_ptr<GlobalLinSys> &plinsys,
82  const AssemblyMapSharedPtr &pLocToGloMap)
83  : Preconditioner(plinsys, pLocToGloMap)
84 {
85 }
86 
88 {
89 }
90 
92 {
93  GlobalSysSolnType sType = m_locToGloMap.lock()->GetGlobalSysSolnType();
94  ASSERTL0(sType == eIterativeStaticCond || sType == ePETScStaticCond,
95  "This type of preconditioning is not implemented "
96  "for this solver");
97 
98  std::shared_ptr<MultiRegions::ExpList> expList =
99  ((m_linsys.lock())->GetLocMat()).lock();
100 
102  expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
103  "LinearPreconSolver");
104 
105  GlobalSysSolnType linSolveType;
106 
107  switch (solveType)
108  {
109  case eLinearPreconPETSc:
110  {
111  linSolveType = ePETScFullMatrix;
112 #ifndef NEKTAR_USING_PETSC
113  NEKERROR(ErrorUtil::efatal, "Nektar++ has not been compiled with "
114  "PETSc support.");
115 #endif
116  break;
117  }
118  case eLinearPreconXxt:
119  default:
120  {
121  linSolveType = eXxtFullMatrix;
122  break;
123  }
124  }
125 
127  m_locToGloMap.lock()->LinearSpaceMap(*expList, linSolveType);
128 
129  // Generate linear solve system.
130  StdRegions::MatrixType mType =
131  m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass
134 
135  GlobalLinSysKey preconKey(mType, m_vertLocToGloMap,
136  (m_linsys.lock())->GetKey().GetConstFactors());
137 
138  switch (solveType)
139  {
140  case eLinearPreconXxt:
141  {
142  m_vertLinsys =
144  preconKey, expList, m_vertLocToGloMap);
145  break;
146  }
147  case eLinearPreconPETSc:
148  {
149 #ifdef NEKTAR_USING_PETSC
150  m_vertLinsys =
152  preconKey, expList, m_vertLocToGloMap);
153 #else
154  ASSERTL0(false, "Nektar++ has not been compiled with "
155  "PETSc support.");
156 #endif
157  }
158  }
159 }
160 
161 /**
162  *
163  */
165  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
166 {
169 }
170 
171 /**
172  *
173  */
175  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
176  const Array<OneD, NekDouble> &pNonVertOutput,
177  Array<OneD, NekDouble> &pVertForce)
178 {
179  GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
180  switch (solvertype)
181  {
184  {
185  int i, val;
186  int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
187  int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
188  // mapping from full space to vertices
189  Array<OneD, int> LocToGloBnd =
190  m_vertLocToGloMap->GetLocalToGlobalBndMap();
191 
192  // Global to local for linear solver (different from above)
193  Array<OneD, int> LocToGlo =
194  m_vertLocToGloMap->GetLocalToGlobalMap();
195 
196  // number of Dir coeffs in from full problem
197  int nDirFull = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
198 
199  Array<OneD, NekDouble> In(nglo, 0.0);
200  Array<OneD, NekDouble> Out(nglo, 0.0);
201 
202  // Gather rhs
203  for (i = 0; i < nloc; ++i)
204  {
205  val = LocToGloBnd[i];
206  if (val >= nDirFull)
207  {
208  In[LocToGlo[i]] = pInput[val - nDirFull];
209  }
210  }
211 
212  // Do solve without enforcing any boundary conditions.
213  m_vertLinsys->SolveLinearSystem(
214  m_vertLocToGloMap->GetNumGlobalCoeffs(), In, Out,
216  m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
217 
218  if (pNonVertOutput != NullNekDouble1DArray)
219  {
220  ASSERTL1(pNonVertOutput.size() >= pOutput.size(),
221  "Non Vert output is not of sufficient length");
222  Vmath::Vcopy(pOutput.size(), pNonVertOutput, 1, pOutput, 1);
223  }
224  else
225  {
226  // Copy input to output as a unit preconditioner on
227  // any other value
228  Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
229  }
230 
231  if (pVertForce != NullNekDouble1DArray)
232  {
233  Vmath::Zero(pVertForce.size(), pVertForce, 1);
234  // Scatter back soln from linear solve
235  for (i = 0; i < nloc; ++i)
236  {
237  val = LocToGloBnd[i];
238  if (val >= nDirFull)
239  {
240  pOutput[val - nDirFull] = Out[LocToGlo[i]];
241  // copy vertex forcing into this vector
242  pVertForce[val - nDirFull] = In[LocToGlo[i]];
243  }
244  }
245  }
246  else
247  {
248  // Scatter back soln from linear solve
249  for (i = 0; i < nloc; ++i)
250  {
251  val = LocToGloBnd[i];
252  if (val >= nDirFull)
253  {
254  pOutput[val - nDirFull] = Out[LocToGlo[i]];
255  }
256  }
257  }
258  }
259  break;
260  default:
261  ASSERTL0(0, "Unsupported solver type");
262  break;
263  }
264 }
265 } // namespace MultiRegions
266 } // 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.
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply a preconditioner to the conjugate gradient method.
std::shared_ptr< AssemblyMap > m_vertLocToGloMap
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce) override
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:51
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255