Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PreconditionerLinear.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Preconditioner.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: Preconditioner definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
42 
43 #ifdef NEKTAR_USING_PETSC
45 #endif
46 
47 #include <LocalRegions/MatrixKey.h>
48 #include <math.h>
49 
50 using namespace std;
51 
52 namespace Nektar
53 {
54  namespace MultiRegions
55  {
56  /**
57  * Registers the class with the Factory.
58  */
59 
60  string PreconditionerLinear::className1
62  "FullLinearSpace",
63  PreconditionerLinear::create,
64  "Full Linear space inverse Preconditioning");
65 
66  std::string PreconditionerLinear::solveType =
67  LibUtilities::SessionReader::RegisterDefaultSolverInfo(
68  "LinearPreconSolver",
69  "Xxt");
70  std::string PreconditionerLinear::solveTypeIds[] = {
71  LibUtilities::SessionReader::RegisterEnumValue(
72  "LinearPreconSolver",
73  "PETSc",
75  LibUtilities::SessionReader::RegisterEnumValue(
76  "LinearPreconSolver",
77  "Xxt",
79  };
80 
81  /**
82  * @class PreconditionerLinear
83  *
84  * This class implements preconditioning for the conjugate
85  * gradient matrix solver.
86  */
87 
88  PreconditionerLinear::PreconditionerLinear(
89  const boost::shared_ptr<GlobalLinSys> &plinsys,
90  const AssemblyMapSharedPtr &pLocToGloMap)
91  : Preconditioner(plinsys, pLocToGloMap)
92  {
93  }
94 
96  {
97  }
98 
100  {
101  GlobalSysSolnType sType = m_locToGloMap->GetGlobalSysSolnType();
102  ASSERTL0(sType == eIterativeStaticCond || sType == ePETScStaticCond,
103  "This type of preconditioning is not implemented "
104  "for this solver");
105 
106  boost::shared_ptr<MultiRegions::ExpList>
107  expList=((m_linsys.lock())->GetLocMat()).lock();
108 
110  expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
111  "LinearPreconSolver");
112 
113  GlobalSysSolnType linSolveType;
114 
115  switch(solveType)
116  {
117  case eLinearPreconPETSc:
118  {
119  linSolveType = ePETScFullMatrix;
120 #ifndef NEKTAR_USING_PETSC
121  ASSERTL0(false, "Nektar++ has not been compiled with "
122  "PETSc support.");
123 #endif
124  }
125  case eLinearPreconXxt:
126  default:
127  {
128  linSolveType = eXxtFullMatrix;
129  break;
130  }
131  }
132 
133  m_vertLocToGloMap = m_locToGloMap->LinearSpaceMap(*expList, linSolveType);
134 
135  // Generate linear solve system.
136  StdRegions::MatrixType mType =
137  m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass ?
140 
141  GlobalLinSysKey preconKey(mType, m_vertLocToGloMap, (m_linsys.lock())->GetKey().GetConstFactors());
142 
143  switch(solveType)
144  {
145  case eLinearPreconXxt:
146  {
148  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
149  break;
150  }
151  case eLinearPreconPETSc:
152  {
153 #ifdef NEKTAR_USING_PETSC
155  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
156 #else
157  ASSERTL0(false, "Nektar++ has not been compiled with "
158  "PETSc support.");
159 #endif
160  }
161  }
162  }
163 
164  /**
165  *
166  */
168  const Array<OneD, NekDouble>& pInput,
169  Array<OneD, NekDouble>& pOutput)
170  {
173  }
174 
175  /**
176  *
177  */
179  const Array<OneD, NekDouble>& pInput,
180  Array<OneD, NekDouble>& pOutput,
181  const Array<OneD, NekDouble>& pNonVertOutput,
182  Array<OneD, NekDouble>& pVertForce)
183  {
184  GlobalSysSolnType solvertype=m_locToGloMap->GetGlobalSysSolnType();
185  switch(solvertype)
186  {
189  {
190  int i,val;
191  int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
192  int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
193  // mapping from full space to vertices
194  Array<OneD, int> LocToGloBnd = m_vertLocToGloMap->GetLocalToGlobalBndMap();
195 
196  // Global to local for linear solver (different from above)
197  Array<OneD, int> LocToGlo = m_vertLocToGloMap->GetLocalToGlobalMap();
198 
199  // number of Dir coeffs in from full problem
200  int nDirFull = m_locToGloMap->GetNumGlobalDirBndCoeffs();
201 
202  Array<OneD,NekDouble> In(nglo,0.0);
203  Array<OneD,NekDouble> Out(nglo,0.0);
204 
205  // Gather rhs
206  for(i = 0; i < nloc; ++i)
207  {
208  val = LocToGloBnd[i];
209  if(val >= nDirFull)
210  {
211  In[LocToGlo[i]] = pInput[val-nDirFull];
212  }
213  }
214 
215  // Do solve without enforcing any boundary conditions.
216  m_vertLinsys->SolveLinearSystem(
217  m_vertLocToGloMap->GetNumGlobalCoeffs(),
218  In,Out,m_vertLocToGloMap,
219  m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
220 
221 
222  if(pNonVertOutput != NullNekDouble1DArray)
223  {
224  ASSERTL1(pNonVertOutput.num_elements() >= pOutput.num_elements(),"Non Vert output is not of sufficient length");
225  Vmath::Vcopy(pOutput.num_elements(),pNonVertOutput,1,pOutput,1);
226  }
227  else
228  {
229  //Copy input to output as a unit preconditioner on
230  //any other value
231  Vmath::Vcopy(pInput.num_elements(),pInput,1,pOutput,1);
232  }
233 
234  if(pVertForce != NullNekDouble1DArray)
235  {
236  Vmath::Zero(pVertForce.num_elements(),pVertForce,1);
237  // Scatter back soln from linear solve
238  for(i = 0; i < nloc; ++i)
239  {
240  val = LocToGloBnd[i];
241  if(val >= nDirFull)
242  {
243  pOutput[val-nDirFull] = Out[LocToGlo[i]];
244  // copy vertex forcing into this vector
245  pVertForce[val-nDirFull] = In[LocToGlo[i]];
246  }
247  }
248  }
249  else
250  {
251  // Scatter back soln from linear solve
252  for(i = 0; i < nloc; ++i)
253  {
254  val = LocToGloBnd[i];
255  if(val >= nDirFull)
256  {
257  pOutput[val-nDirFull] = Out[LocToGlo[i]];
258  }
259  }
260  }
261  }
262  break;
263  default:
264  ASSERTL0(0,"Unsupported solver type");
265  break;
266  }
267  }
268  }
269 }
270 
271 
272 
273 
274 
275 
boost::shared_ptr< AssemblyMap > m_locToGloMap
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.
static Array< OneD, NekDouble > NullNekDouble1DArray
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
PreconFactory & GetPreconFactory()
STL namespace.
boost::shared_ptr< AssemblyMap > m_vertLocToGloMap
const boost::weak_ptr< GlobalLinSys > m_linsys
Describe a linear system.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215