Nektar++
GlobalLinSysIterativeFull.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GlobalLinSysIterativeFull.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: GlobalLinSysIterativeFull definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
37#include <map>
38
39using namespace std;
40
42{
43/**
44 * @class GlobalLinSysIterativeCG
45 *
46 * This class implements a conjugate gradient matrix solver.
47 * Preconditioning is implemented using a Jacobi (diagonal)
48 * preconditioner.
49 */
50
51/**
52 * Registers the class with the Factory.
53 */
56 "IterativeFull", GlobalLinSysIterativeFull::create,
57 "Iterative solver for full matrix system.");
58
59/**
60 * Constructor for full direct matrix solve.
61 * @param pKey Key specifying matrix to solve.
62 * @param pExp Shared pointer to expansion list for applying
63 * matrix evaluations.
64 * @param pLocToGloMap Local to global mapping.
65 */
67 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExp,
68 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
69 : GlobalLinSys(pKey, pExp, pLocToGloMap),
70 GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
71{
73 "This routine should only be used when using an Iterative "
74 "conjugate gradient matrix solve.");
75}
76
77/**
78 * Solve a global linear system with Dirichlet forcing using a
79 * conjugate gradient method. This routine performs handling of the
80 * Dirichlet forcing terms and wraps the underlying iterative solver
81 * used for the remaining degrees of freedom.
82 *
83 * Consider solving for \f$x\f$, the matrix system \f$Ax=b\f$, where
84 * \f$b\f$ is known. To enforce the Dirichlet terms we instead solve
85 * \f[A(x-x_0) = b - Ax_0 \f]
86 * where \f$x_0\f$ is the Dirichlet forcing.
87 *
88 * @param pInput RHS of linear system, \f$b\f$.
89 * @param pOutput On input, values of dirichlet degrees
90 * of freedom with initial guess on other values.
91 * On output, the solution \f$x\f$.
92 * @param pLocToGloMap Local to global mapping.
93 * @param pDirForcing Precalculated Dirichlet forcing.
94 */
96 const Array<OneD, const NekDouble> &pLocInput,
97 Array<OneD, NekDouble> &pLocOutput,
98 const AssemblyMapSharedPtr &pLocToGloMap,
99 const Array<OneD, const NekDouble> &pDirForcing)
100{
101 m_locToGloMap = pLocToGloMap;
102
103 bool dirForcCalculated = (bool)pDirForcing.size();
104 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
105 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
106 int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
107
108 int nDirTotal = nDirDofs;
109 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
110 expList->GetComm()->GetRowComm()->AllReduce(nDirTotal,
112
113 if (nDirTotal)
114 {
115 Array<OneD, NekDouble> rhs(nLocDofs);
116
117 // Calculate the Dirichlet forcing
118 if (dirForcCalculated)
119 {
120 // Assume pDirForcing is in local space
121 ASSERTL0(
122 pDirForcing.size() >= nLocDofs,
123 "DirForcing is not of sufficient size. Is it in local space?");
124 Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, rhs, 1);
125 }
126 else
127 {
128 // Calculate initial condition and Dirichlet forcing and subtract it
129 // from the rhs
130 expList->GeneralMatrixOp(m_linSysKey, pLocOutput, rhs);
131
132 // Iterate over all the elements computing Robin BCs where
133 // necessary
134 for (auto &r : m_robinBCInfo) // add robin mass matrix
135 {
138
139 int n = r.first;
140 int offset = expList->GetCoeff_Offset(n);
141
142 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
143 // Add local matrix contribution
144 for (rBC = r.second; rBC; rBC = rBC->next)
145 {
146 vExp->AddRobinTraceContribution(
147 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
148 pLocOutput + offset, rhsloc = rhs + offset);
149 }
150 }
151 Vmath::Vsub(nLocDofs, pLocInput, 1, rhs, 1, rhs, 1);
152 }
153
154 if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
155 {
156 Array<OneD, NekDouble> diff(nLocDofs);
157
158 // Solve for perturbation from initial guess in pOutput
159 SolveLinearSystem(nGlobDofs, rhs, diff, pLocToGloMap, nDirDofs);
160
161 // Add back initial and boundary condition
162 Vmath::Vadd(nLocDofs, diff, 1, pLocOutput, 1, pLocOutput, 1);
163 }
164 else
165 {
166 ASSERTL0(false, "Need DG solve if using Dir BCs");
167 }
168 }
169 else
170 {
171 SolveLinearSystem(nGlobDofs, pLocInput, pLocOutput, pLocToGloMap,
172 nDirDofs);
173 }
174}
175
176/**
177 *
178 */
180 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
181{
182 bool isLocal = m_linsol->IsLocal();
183
184 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
185
186 AssemblyMapSharedPtr asmMap = m_locToGloMap.lock();
187
188 int ncoeffs = expList->GetNcoeffs();
189 Array<OneD, NekDouble> InputLoc, OutputLoc;
190
191 if (isLocal)
192 {
193 InputLoc = pInput;
194 OutputLoc = pOutput;
195
196 // Perform matrix-vector operation A*d_i
197 expList->GeneralMatrixOp(m_linSysKey, InputLoc, OutputLoc);
198 }
199 else
200 {
201 InputLoc = Array<OneD, NekDouble>(ncoeffs);
202 OutputLoc = Array<OneD, NekDouble>(ncoeffs);
203
204 asmMap->GlobalToLocal(pInput, InputLoc);
205 // Perform matrix-vector operation A*d_i
206 expList->GeneralMatrixOp(m_linSysKey, InputLoc, OutputLoc);
207 }
208
209 // Apply robin boundary conditions to the solution.
210 for (auto &r : m_robinBCInfo) // add robin mass matrix
211 {
214
215 int n = r.first;
216
217 int offset = expList->GetCoeff_Offset(n);
218 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
219
220 // add local matrix contribution
221 for (rBC = r.second; rBC; rBC = rBC->next)
222 {
223 vExp->AddRobinTraceContribution(
224 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, InputLoc + offset,
225 tmp = OutputLoc + offset);
226 }
227 }
228
229 // put back in global coeffs
230 if (isLocal == false)
231 {
232 asmMap->Assemble(OutputLoc, pOutput);
233 }
234}
235
236/**
237 *
238 */
240{
241 m_map = m_locToGloMap.lock()->GetGlobalToUniversalMapUnique();
242}
243
244/**
245 *
246 */
248 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
249 Array<OneD, NekDouble> &pOutput, const AssemblyMapSharedPtr &plocToGloMap,
250 const int nDir)
251{
252
253 if (!m_linsol)
254 {
256 m_expList.lock()->GetComm()->GetRowComm();
258 m_expList.lock()->GetSession();
259
260 // Check such a module exists for this equation.
263 "NekLinSysIter '" + m_linSysIterSolver +
264 "' is not defined.\n");
265
266 // Create the key to hold solver settings
267 auto sysKey = LibUtilities::NekSysKey();
268 string variable = plocToGloMap->GetVariable();
269
270 // Either get the solnInfo from <GlobalSysSolInfo> or from
271 // <Parameters>
272 if (pSession->DefinesGlobalSysSolnInfo(variable,
273 "IterativeSolverTolerance"))
274 {
275 sysKey.m_NekLinSysTolerance = boost::lexical_cast<double>(
276 pSession
277 ->GetGlobalSysSolnInfo(variable, "IterativeSolverTolerance")
278 .c_str());
279 }
280 else
281 {
282 pSession->LoadParameter("IterativeSolverTolerance",
283 sysKey.m_NekLinSysTolerance,
285 }
286
287 if (pSession->DefinesGlobalSysSolnInfo(variable,
288 "NekLinSysMaxIterations"))
289 {
290 sysKey.m_NekLinSysMaxIterations = boost::lexical_cast<int>(
291 pSession
292 ->GetGlobalSysSolnInfo(variable, "NekLinSysMaxIterations")
293 .c_str());
294 }
295 else
296 {
297 pSession->LoadParameter("NekLinSysMaxIterations",
298 sysKey.m_NekLinSysMaxIterations, 5000);
299 }
300
301 if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysMaxStorage"))
302 {
303 sysKey.m_LinSysMaxStorage = boost::lexical_cast<int>(
304 pSession->GetGlobalSysSolnInfo(variable, "LinSysMaxStorage")
305 .c_str());
306 }
307 else
308 {
309 pSession->LoadParameter("LinSysMaxStorage",
310 sysKey.m_LinSysMaxStorage, 100);
311 }
312
313 if (pSession->DefinesGlobalSysSolnInfo(variable, "GMRESMaxHessMatBand"))
314 {
315 sysKey.m_KrylovMaxHessMatBand = boost::lexical_cast<int>(
316 pSession->GetGlobalSysSolnInfo(variable, "GMRESMaxHessMatBand")
317 .c_str());
318 }
319 else
320 {
321 pSession->LoadParameter("GMRESMaxHessMatBand",
322 sysKey.m_KrylovMaxHessMatBand,
323 sysKey.m_LinSysMaxStorage + 1);
324 }
325
326 // The following settings have no correponding tests and are rarely
327 // used.
328 pSession->MatchSolverInfo("GMRESLeftPrecon", "True",
329 sysKey.m_NekLinSysLeftPrecon, false);
330 pSession->MatchSolverInfo("GMRESRightPrecon", "True",
331 sysKey.m_NekLinSysRightPrecon, true);
332
334 m_linSysIterSolver, pSession, vRowComm, nGlobal - nDir, sysKey);
335
336 m_linsol->SetSysOperators(m_NekSysOp);
337 v_UniqueMap();
338 m_linsol->SetUniversalUniqueMap(m_map);
339 }
340
341 if (!m_precon)
342 {
343 m_precon = CreatePrecon(plocToGloMap);
344 m_precon->BuildPreconditioner();
345 }
346
347 m_linsol->SetRhsMagnitude(m_isAbsoluteTolerance ? 1.0 : m_rhs_magnitude);
348
349 if (m_useProjection)
350 {
351 Array<OneD, NekDouble> gloIn(nGlobal);
352 Array<OneD, NekDouble> gloOut(nGlobal, 0.0);
353 plocToGloMap->Assemble(pInput, gloIn);
354 DoProjection(nGlobal, gloIn, gloOut, nDir, m_isAconjugate);
355 plocToGloMap->GlobalToLocal(gloOut, pOutput);
356 }
357 else
358 {
359 if (m_linsol->IsLocal())
360 {
361 int nLocDofs = plocToGloMap->GetNumLocalCoeffs();
362 Vmath::Zero(nLocDofs, pOutput, 1);
363 m_linsol->SolveSystem(nLocDofs, pInput, pOutput, nDir);
364 }
365 else
366 {
367 Array<OneD, NekDouble> gloIn(nGlobal);
368 Array<OneD, NekDouble> gloOut(nGlobal, 0.0);
369 plocToGloMap->Assemble(pInput, gloIn);
370 m_linsol->SolveSystem(nGlobal, gloIn, gloOut, nDir);
371 plocToGloMap->GlobalToLocal(gloOut, pOutput);
372 }
373 }
374}
375
376} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
A global linear system.
Definition: GlobalLinSys.h:70
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:122
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:124
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:190
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:120
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map.
static GlobalLinSysSharedPtr create(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
GlobalLinSysIterativeFull(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray) override
Solve the linear system for given input and output vectors using a specified local to global map.
void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir) override
Solve the matrix system.
LibUtilities::NekLinSysIterSharedPtr m_linsol
void DoProjection(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir, const bool isAconjugate)
projection technique
bool m_useProjection
Whether to apply projection technique.
Array< OneD, int > m_map
Global to universal unique map.
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
static const NekDouble kNekIterativeTol
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.hpp:180
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.hpp:220
STL namespace.