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
41namespace Nektar
42{
43namespace MultiRegions
44{
45/**
46 * @class GlobalLinSysIterativeCG
47 *
48 * This class implements a conjugate gradient matrix solver.
49 * Preconditioning is implemented using a Jacobi (diagonal)
50 * preconditioner.
51 */
52
53/**
54 * Registers the class with the Factory.
55 */
58 "IterativeFull", GlobalLinSysIterativeFull::create,
59 "Iterative solver for full matrix system.");
60
61/**
62 * Constructor for full direct matrix solve.
63 * @param pKey Key specifying matrix to solve.
64 * @param pExp Shared pointer to expansion list for applying
65 * matrix evaluations.
66 * @param pLocToGloMap Local to global mapping.
67 */
69 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExp,
70 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
71 : GlobalLinSys(pKey, pExp, pLocToGloMap),
72 GlobalLinSysIterative(pKey, pExp, pLocToGloMap)
73{
75 "This routine should only be used when using an Iterative "
76 "conjugate gradient matrix solve.");
77}
78
79/**
80 *
81 */
83{
84}
85
86/**
87 * Solve a global linear system with Dirichlet forcing using a
88 * conjugate gradient method. This routine performs handling of the
89 * Dirichlet forcing terms and wraps the underlying iterative solver
90 * used for the remaining degrees of freedom.
91 *
92 * Consider solving for \f$x\f$, the matrix system \f$Ax=b\f$, where
93 * \f$b\f$ is known. To enforce the Dirichlet terms we instead solve
94 * \f[A(x-x_0) = b - Ax_0 \f]
95 * where \f$x_0\f$ is the Dirichlet forcing.
96 *
97 * @param pInput RHS of linear system, \f$b\f$.
98 * @param pOutput On input, values of dirichlet degrees
99 * of freedom with initial guess on other values.
100 * On output, the solution \f$x\f$.
101 * @param pLocToGloMap Local to global mapping.
102 * @param pDirForcing Precalculated Dirichlet forcing.
103 */
105 const Array<OneD, const NekDouble> &pLocInput,
106 Array<OneD, NekDouble> &pLocOutput,
107 const AssemblyMapSharedPtr &pLocToGloMap,
108 const Array<OneD, const NekDouble> &pDirForcing)
109{
110 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
111 bool vCG = false;
112 m_locToGloMap = pLocToGloMap;
113
114 if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
115 {
116 vCG = true;
117 }
118 else if (std::dynamic_pointer_cast<AssemblyMapDG>(pLocToGloMap))
119 {
120 vCG = false;
121 }
122 else
123 {
124 NEKERROR(ErrorUtil::efatal, "Unknown map type");
125 }
126
127 bool dirForcCalculated = (bool)pDirForcing.size();
128 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
129 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
130 int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
131
132 int nDirTotal = nDirDofs;
133 expList->GetComm()->GetRowComm()->AllReduce(nDirTotal,
135
136 Array<OneD, NekDouble> tmp(nLocDofs);
137 Array<OneD, NekDouble> tmp1(nLocDofs);
138
139 if (nDirTotal)
140 {
141 // calculate the Dirichlet forcing
142 if (dirForcCalculated)
143 {
144 Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, tmp1, 1);
145 }
146 else
147 {
148 // Calculate the dirichlet forcing B_b (== X_b) and
149 // substract it from the rhs
150 expList->GeneralMatrixOp(m_linSysKey, pLocOutput, tmp);
151
152 // Iterate over all the elements computing Robin BCs where
153 // necessary
154 for (auto &r : m_robinBCInfo) // add robin mass matrix
155 {
158
159 int n = r.first;
160 int offset = expList->GetCoeff_Offset(n);
161
162 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
163 // add local matrix contribution
164 for (rBC = r.second; rBC; rBC = rBC->next)
165 {
166 vExp->AddRobinTraceContribution(
167 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
168 pLocOutput + offset, tmploc = tmp + offset);
169 }
170 }
171
172 Vmath::Vsub(nLocDofs, pLocInput, 1, tmp, 1, tmp1, 1);
173 }
174 if (vCG)
175 {
176 // solve for perturbation from initial guess in pOutput
177 SolveLinearSystem(nGlobDofs, tmp1, tmp, pLocToGloMap, nDirDofs);
178
179 // Add back initial condition
180 Vmath::Vadd(nLocDofs, tmp, 1, pLocOutput, 1, pLocOutput, 1);
181 }
182 else
183 {
184 ASSERTL0(false, "Need DG solve if using Dir BCs");
185 }
186 }
187 else
188 {
189 SolveLinearSystem(nGlobDofs, pLocInput, pLocOutput, pLocToGloMap,
190 nDirDofs);
191 }
192}
193
194/**
195 *
196 */
198 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
199{
200 bool isLocal = m_linsol->IsLocal();
201
202 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
203
204 AssemblyMapSharedPtr asmMap = m_locToGloMap.lock();
205
206 int ncoeffs = expList->GetNcoeffs();
207 Array<OneD, NekDouble> InputLoc, OutputLoc;
208
209 if (isLocal)
210 {
211 InputLoc = pInput;
212 OutputLoc = pOutput;
213
214 // Perform matrix-vector operation A*d_i
215 expList->GeneralMatrixOp(m_linSysKey, InputLoc, OutputLoc);
216 }
217 else
218 {
219 InputLoc = Array<OneD, NekDouble>(ncoeffs);
220 OutputLoc = Array<OneD, NekDouble>(ncoeffs);
221
222 asmMap->GlobalToLocal(pInput, InputLoc);
223 // Perform matrix-vector operation A*d_i
224 expList->GeneralMatrixOp(m_linSysKey, InputLoc, OutputLoc);
225 }
226
227 // Apply robin boundary conditions to the solution.
228 for (auto &r : m_robinBCInfo) // add robin mass matrix
229 {
232
233 int n = r.first;
234
235 int offset = expList->GetCoeff_Offset(n);
236 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
237
238 // add local matrix contribution
239 for (rBC = r.second; rBC; rBC = rBC->next)
240 {
241 vExp->AddRobinTraceContribution(
242 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, InputLoc + offset,
243 tmp = OutputLoc + offset);
244 }
245 }
246
247 // put back in global coeffs
248 if (isLocal == false)
249 {
250 asmMap->Assemble(OutputLoc, pOutput);
251 }
252}
253
254/**
255 *
256 */
258{
259 m_map = m_locToGloMap.lock()->GetGlobalToUniversalMapUnique();
260}
261
262/**
263 *
264 */
266 const int nGlobal, const Array<OneD, const NekDouble> &pInput,
267 Array<OneD, NekDouble> &pOutput, const AssemblyMapSharedPtr &plocToGloMap,
268 const int nDir)
269{
270
271 if (!m_linsol)
272 {
274 m_expList.lock()->GetComm()->GetRowComm();
276 m_expList.lock()->GetSession();
277
278 // Check such a module exists for this equation.
281 "NekLinSysIter '" + m_linSysIterSolver +
282 "' is not defined.\n");
284 m_linSysIterSolver, pSession, vRowComm, nGlobal - nDir,
286
287 m_linsol->SetSysOperators(m_NekSysOp);
288 v_UniqueMap();
289 m_linsol->setUniversalUniqueMap(m_map);
290 }
291
292 if (!m_precon)
293 {
294 m_precon = CreatePrecon(plocToGloMap);
295 m_precon->BuildPreconditioner();
296 }
297
298 m_linsol->setRhsMagnitude(m_rhs_magnitude);
299
300 if (m_useProjection)
301 {
302 Array<OneD, NekDouble> gloIn(nGlobal);
303 Array<OneD, NekDouble> gloOut(nGlobal, 0.0);
304 plocToGloMap->Assemble(pInput, gloIn);
305 DoProjection(nGlobal, gloIn, gloOut, nDir, m_tolerance, m_isAconjugate);
306 plocToGloMap->GlobalToLocal(gloOut, pOutput);
307 }
308 else
309 {
310 if (m_linsol->IsLocal())
311 {
312 int nLocDofs = plocToGloMap->GetNumLocalCoeffs();
313 Vmath::Zero(nLocDofs, pOutput, 1);
314 m_linsol->SolveSystem(nLocDofs, pInput, pOutput, nDir, m_tolerance);
315 }
316 else
317 {
318 Array<OneD, NekDouble> gloIn(nGlobal);
319 Array<OneD, NekDouble> gloOut(nGlobal, 0.0);
320 plocToGloMap->Assemble(pInput, gloIn);
321 m_linsol->SolveSystem(nGlobal, gloIn, gloOut, nDir, m_tolerance);
322 plocToGloMap->GlobalToLocal(gloOut, pOutput);
323 }
324 }
325}
326} // namespace MultiRegions
327} // 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
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
A global linear system.
Definition: GlobalLinSys.h:72
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:124
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:126
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:192
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:122
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.
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
virtual 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.
virtual 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 NekDouble tol, const bool isAconjugate)
projection technique
bool m_useProjection
Whether to apply projection technique.
NekDouble m_tolerance
Tolerance of iterative solver.
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:57
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:354
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
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.cpp:414