Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 namespace Nektar
51 {
52  namespace MultiRegions
53  {
54  /**
55  * Registers the class with the Factory.
56  */
57 
60  "FullLinearSpace",
62  "Full Linear space inverse Preconditioning");
63 
66  "LinearPreconSolver",
67  "Xxt");
68  std::string PreconditionerLinear::solveTypeIds[] = {
70  "LinearPreconSolver",
71  "PETSc",
74  "LinearPreconSolver",
75  "Xxt",
77  };
78 
79  /**
80  * @class PreconditionerLinear
81  *
82  * This class implements preconditioning for the conjugate
83  * gradient matrix solver.
84  */
85 
87  const boost::shared_ptr<GlobalLinSys> &plinsys,
88  const AssemblyMapSharedPtr &pLocToGloMap)
89  : Preconditioner(plinsys, pLocToGloMap)
90  {
91  }
92 
94  {
95  }
96 
98  {
99  GlobalSysSolnType sType = m_locToGloMap->GetGlobalSysSolnType();
100  ASSERTL0(sType == eIterativeStaticCond || sType == ePETScStaticCond,
101  "This type of preconditioning is not implemented "
102  "for this solver");
103 
104  boost::shared_ptr<MultiRegions::ExpList>
105  expList=((m_linsys.lock())->GetLocMat()).lock();
106 
108  expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
109  "LinearPreconSolver");
110 
111  GlobalSysSolnType linSolveType;
112 
113  switch(solveType)
114  {
115  case eLinearPreconPETSc:
116  {
117  linSolveType = ePETScFullMatrix;
118 #ifndef NEKTAR_USING_PETSC
119  ASSERTL0(false, "Nektar++ has not been compiled with "
120  "PETSc support.");
121 #endif
122  }
123  case eLinearPreconXxt:
124  default:
125  {
126  linSolveType = eXxtFullMatrix;
127  break;
128  }
129  }
130 
131  m_vertLocToGloMap = m_locToGloMap->LinearSpaceMap(*expList, linSolveType);
132 
133  // Generate linear solve system.
134  StdRegions::MatrixType mType =
135  m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass ?
138 
139  GlobalLinSysKey preconKey(mType, m_vertLocToGloMap, (m_linsys.lock())->GetKey().GetConstFactors());
140 
141  switch(solveType)
142  {
143  case eLinearPreconXxt:
144  {
146  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
147  break;
148  }
149  case eLinearPreconPETSc:
150  {
151 #ifdef NEKTAR_USING_PETSC
153  AllocateSharedPtr(preconKey,expList,m_vertLocToGloMap);
154 #else
155  ASSERTL0(false, "Nektar++ has not been compiled with "
156  "PETSc support.");
157 #endif
158  }
159  }
160  }
161 
162  /**
163  *
164  */
166  const Array<OneD, NekDouble>& pInput,
167  Array<OneD, NekDouble>& pOutput)
168  {
171  }
172 
173  /**
174  *
175  */
177  const Array<OneD, NekDouble>& pInput,
178  Array<OneD, NekDouble>& pOutput,
179  const Array<OneD, NekDouble>& pNonVertOutput,
180  Array<OneD, NekDouble>& pVertForce)
181  {
182  GlobalSysSolnType solvertype=m_locToGloMap->GetGlobalSysSolnType();
183  switch(solvertype)
184  {
187  {
188  int i,val;
189  int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
190  int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
191  // mapping from full space to vertices
192  Array<OneD, int> LocToGloBnd = m_vertLocToGloMap->GetLocalToGlobalBndMap();
193 
194  // Global to local for linear solver (different from above)
195  Array<OneD, int> LocToGlo = m_vertLocToGloMap->GetLocalToGlobalMap();
196 
197  // number of Dir coeffs in from full problem
198  int nDirFull = m_locToGloMap->GetNumGlobalDirBndCoeffs();
199 
200  Array<OneD,NekDouble> In(nglo,0.0);
201  Array<OneD,NekDouble> Out(nglo,0.0);
202 
203  // Gather rhs
204  for(i = 0; i < nloc; ++i)
205  {
206  val = LocToGloBnd[i];
207  if(val >= nDirFull)
208  {
209  In[LocToGlo[i]] = pInput[val-nDirFull];
210  }
211  }
212 
213  // Do solve without enforcing any boundary conditions.
214  m_vertLinsys->SolveLinearSystem(
215  m_vertLocToGloMap->GetNumGlobalCoeffs(),
216  In,Out,m_vertLocToGloMap,
217  m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
218 
219 
220  if(pNonVertOutput != NullNekDouble1DArray)
221  {
222  ASSERTL1(pNonVertOutput.num_elements() >= pOutput.num_elements(),"Non Vert output is not of sufficient length");
223  Vmath::Vcopy(pOutput.num_elements(),pNonVertOutput,1,pOutput,1);
224  }
225  else
226  {
227  //Copy input to output as a unit preconditioner on
228  //any other value
229  Vmath::Vcopy(pInput.num_elements(),pInput,1,pOutput,1);
230  }
231 
232  if(pVertForce != NullNekDouble1DArray)
233  {
234  Vmath::Zero(pVertForce.num_elements(),pVertForce,1);
235  // Scatter back soln from linear solve
236  for(i = 0; i < nloc; ++i)
237  {
238  val = LocToGloBnd[i];
239  if(val >= nDirFull)
240  {
241  pOutput[val-nDirFull] = Out[LocToGlo[i]];
242  // copy vertex forcing into this vector
243  pVertForce[val-nDirFull] = In[LocToGlo[i]];
244  }
245  }
246  }
247  else
248  {
249  // Scatter back soln from linear solve
250  for(i = 0; i < nloc; ++i)
251  {
252  val = LocToGloBnd[i];
253  if(val >= nDirFull)
254  {
255  pOutput[val-nDirFull] = Out[LocToGlo[i]];
256  }
257  }
258  }
259  }
260  break;
261  default:
262  ASSERTL0(0,"Unsupported solver type");
263  break;
264  }
265  }
266  }
267 }
268 
269 
270 
271 
272 
273 
boost::shared_ptr< AssemblyMap > m_locToGloMap
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
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...
static std::string className1
Name of class.
PreconFactory & GetPreconFactory()
boost::shared_ptr< AssemblyMap > m_vertLocToGloMap
PreconditionerLinear(const boost::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
const boost::weak_ptr< GlobalLinSys > m_linsys
static PreconditionerSharedPtr create(const boost::shared_ptr< GlobalLinSys > &plinsys, const boost::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
Describe a linear system.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
#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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215