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
40
41#ifdef NEKTAR_USING_PETSC
43#endif
44
46#include <cmath>
47
48using namespace std;
49
51{
52/**
53 * Registers the class with the Factory.
54 */
55
58 "FullLinearSpace", PreconditionerLinear::create,
59 "Full Linear space inverse Preconditioning");
60
63 "Xxt");
66 "LinearPreconSolver", "PETSc", MultiRegions::eLinearPreconPETSc),
68 "LinearPreconSolver", "Xxt", MultiRegions::eLinearPreconXxt)};
69
70/**
71 * @class PreconditionerLinear
72 *
73 * This class implements preconditioning for the conjugate
74 * gradient matrix solver.
75 */
76
78 const std::shared_ptr<GlobalLinSys> &plinsys,
79 const AssemblyMapSharedPtr &pLocToGloMap)
80 : Preconditioner(plinsys, pLocToGloMap)
81{
82}
83
85{
86}
87
89{
90 GlobalSysSolnType sType = m_locToGloMap.lock()->GetGlobalSysSolnType();
92 "This type of preconditioning is not implemented "
93 "for this solver");
94
95 std::shared_ptr<MultiRegions::ExpList> expList =
96 ((m_linsys.lock())->GetLocMat()).lock();
97
99 expList->GetSession()->GetSolverInfoAsEnum<LinearPreconSolver>(
100 "LinearPreconSolver");
101
102 GlobalSysSolnType linSolveType;
103
104 switch (solveType)
105 {
107 {
108 linSolveType = ePETScFullMatrix;
109#ifndef NEKTAR_USING_PETSC
110 NEKERROR(ErrorUtil::efatal, "Nektar++ has not been compiled with "
111 "PETSc support.");
112#endif
113 break;
114 }
115 case eLinearPreconXxt:
116 default:
117 {
118 linSolveType = eXxtFullMatrix;
119 break;
120 }
121 }
122
124 m_locToGloMap.lock()->LinearSpaceMap(*expList, linSolveType);
125
127
128 // Generate linear solve system.
130 m_linsys.lock()->GetKey().GetMatrixType() == StdRegions::eMass
133
134 GlobalLinSysKey preconKey(mType, m_vertLocToGloMap,
135 (m_linsys.lock())->GetKey().GetConstFactors());
136
137 switch (solveType)
138 {
139 case eLinearPreconXxt:
140 {
143 preconKey, expList, m_vertLocToGloMap);
144 break;
145 }
147 {
148#ifdef NEKTAR_USING_PETSC
151 preconKey, expList, m_vertLocToGloMap);
152#else
153 ASSERTL0(false, "Nektar++ has not been compiled with "
154 "PETSc support.");
155#endif
156 }
157 }
158}
159
160/**
161 *
162 */
164 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
165 const bool &isLocal)
166{
167 ASSERTL0(isLocal == false, "PreconditionerLinear is only set up for Global"
168 " iteratives sovles");
171}
172
173/**
174 *
175 */
177 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
178 const Array<OneD, NekDouble> &pNonVertOutput,
179 Array<OneD, NekDouble> &pVertForce)
180{
181 GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
182 switch (solvertype)
183 {
186 {
187 int i, val;
188 int nloc = m_vertLocToGloMap->GetNumLocalCoeffs();
189 // mapping from full space to vertices
190 Array<OneD, int> LocToGloBnd =
191 m_vertLocToGloMap->GetLocalToGlobalBndMap();
192
193 // Global to local for linear solver (different from above)
194 Array<OneD, int> LocToGlo =
195 m_vertLocToGloMap->GetLocalToGlobalMap();
196
197 // number of Dir coeffs in from full problem
198 int nDirFull = m_locToGloMap.lock()->GetNumGlobalDirBndCoeffs();
199
200#if 0
201 int nglo = m_vertLocToGloMap->GetNumGlobalCoeffs();
202 Array<OneD, NekDouble> In(nglo, 0.0);
203 Array<OneD, NekDouble> Out(nglo, 0.0);
204
205 // Gather rhs to local linear space.
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#else
215 Array<OneD, NekDouble> In(nloc, 0.0);
216 Array<OneD, NekDouble> Out(nloc, 0.0);
217
218 // Gather rhs to local linear space.
219 for (i = 0; i < nloc; ++i)
220 {
221 val = LocToGloBnd[i];
222 if (val >= nDirFull)
223 {
224 In[i] = pInput[val - nDirFull] * m_invMult[i];
225 }
226 }
227#endif
228
229 // Do solve without enforcing any boundary conditions.
230 m_vertLinsys->SolveLinearSystem(
231 m_vertLocToGloMap->GetNumGlobalCoeffs(), In, Out,
233 m_vertLocToGloMap->GetNumGlobalDirBndCoeffs());
234
235 if (pNonVertOutput != NullNekDouble1DArray)
236 {
237 ASSERTL1(pNonVertOutput.size() >= pOutput.size(),
238 "Non Vert output is not of sufficient length");
239 Vmath::Vcopy(pOutput.size(), pNonVertOutput, 1, pOutput, 1);
240 }
241 else
242 {
243 // Copy input to output as a unit preconditioner on
244 // any other value
245 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
246 }
247
248#if 0
249 if (pVertForce != NullNekDouble1DArray)
250 {
251 Vmath::Zero(pVertForce.size(), pVertForce, 1);
252 // Scatter back soln from linear solve
253 for (i = 0; i < nloc; ++i)
254 {
255 val = LocToGloBnd[i];
256 if (val >= nDirFull)
257 {
258 pOutput[val - nDirFull] = Out[LocToGlo[i]];
259 // copy vertex forcing into this vector
260 pVertForce[val - nDirFull] = In[LocToGlo[i]];
261 }
262 }
263 }
264 else
265 {
266 // Scatter back soln from linear solve
267 for (i = 0; i < nloc; ++i)
268 {
269 val = LocToGloBnd[i];
270 if (val >= nDirFull)
271 {
272 pOutput[val - nDirFull] = Out[LocToGlo[i]];
273 }
274 }
275 }
276#else
277 if (pVertForce != NullNekDouble1DArray)
278 {
279 Vmath::Zero(pVertForce.size(), pVertForce, 1);
280 // Scatter back soln from linear solve
281 for (i = 0; i < nloc; ++i)
282 {
283 val = LocToGloBnd[i];
284 if (val >= nDirFull)
285 {
286 pOutput[val - nDirFull] = Out[i];
287 // copy vertex forcing into this vector
288 pVertForce[val - nDirFull] = In[i];
289 }
290 }
291 }
292 else
293 {
294 // Scatter back soln from linear solve
295 for (i = 0; i < nloc; ++i)
296 {
297 val = LocToGloBnd[i];
298 if (val >= nDirFull)
299 {
300 pOutput[val - nDirFull] = Out[i];
301 }
302 }
303 }
304#endif
305 }
306 break;
307 default:
308 ASSERTL0(0, "Unsupported solver type");
309 break;
310 }
311}
312
313/**
314 * Create the inverse multiplicity map.
315 * @param locToGloMap Local to global mapping information.
316 */
318 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
319{
320 const Array<OneD, const int> &vMap = pLocToGloMap->GetLocalToGlobalMap();
321 unsigned int nGlo = pLocToGloMap->GetNumGlobalCoeffs();
322 unsigned int nEntries = pLocToGloMap->GetNumLocalCoeffs();
323 unsigned int i;
324
325 // Count the multiplicity of each global DOF on this process
326 Array<OneD, NekDouble> vCounts(nGlo, 0.0);
327 for (i = 0; i < nEntries; ++i)
328 {
329 vCounts[vMap[i]] += 1.0;
330 }
331
332 // Get universal multiplicity by globally assembling counts
333 pLocToGloMap->UniversalAssemble(vCounts);
334
335 // Construct a map of 1/multiplicity for use in XXT solve
337 for (i = 0; i < nEntries; ++i)
338 {
339 m_invMult[i] = 1.0 / vCounts[vMap[i]];
340 }
341}
342
343} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
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
void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false) override
Apply a preconditioner to the conjugate gradient method.
PreconditionerLinear(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
void SetupInvMult(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
std::shared_ptr< AssemblyMap > m_vertLocToGloMap
static std::string className1
Name of class.
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:50
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.