Nektar++
PreconditionerDiagonal.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PreconditionerDiagonal.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: PreconditionerDiagonal definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
39#include <cmath>
40
41using namespace std;
42
44{
45/**
46 * Registers the class with the Factory.
47 */
50 "Diagonal", PreconditionerDiagonal::create, "Diagonal Preconditioning");
51/**
52 * @class Preconditioner
53 *
54 * This class implements diagonal preconditioning for the conjugate
55 * gradient matrix solver.
56 */
57
59 const std::shared_ptr<GlobalLinSys> &plinsys,
60 const AssemblyMapSharedPtr &pLocToGloMap)
61 : Preconditioner(plinsys, pLocToGloMap)
62{
63}
64
66{
67}
68
70{
71 GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
72 if (solvertype == eIterativeFull)
73 {
75 }
76 else if (solvertype == eIterativeStaticCond ||
77 solvertype == eIterativeMultiLevelStaticCond ||
78 solvertype == ePETScStaticCond ||
79 solvertype == ePETScMultiLevelStaticCond)
80 {
82 }
83 else
84 {
85 ASSERTL0(0, "Unsupported solver type");
86 }
87}
88
89/**
90 * Diagonal preconditioner computed by summing the relevant elements of
91 * the local matrix system.
92 */
94{
95 std::shared_ptr<MultiRegions::ExpList> expList =
96 ((m_linsys.lock())->GetLocMat()).lock();
97
98 auto asmMap = m_locToGloMap.lock();
99 int nGlobal = asmMap->GetNumGlobalCoeffs();
100 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
101 int nInt = nGlobal - nDir;
102 int nElmt = expList->GetNumElmts();
103
104 Array<OneD, NekDouble> vOutput(nGlobal, 0.0);
105 for (int n = 0, cnt = 0; n < nElmt; ++n)
106 {
107 auto loc_mat = (m_linsys.lock())->GetBlock(n);
108 int loc_row = loc_mat->GetRows();
109 for (int i = 0; i < loc_row; ++i)
110 {
111 int gid1 = asmMap->GetLocalToGlobalMap(cnt + i);
112 vOutput[gid1] += (*loc_mat)(i, i);
113 }
114 cnt += loc_row;
115 }
116
117 // Assemble diagonal contributions across processes
118 asmMap->UniversalAssemble(vOutput);
119
121 Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
122}
123
124/**
125 * Diagonal preconditioner defined as the inverse of the main
126 * diagonal of the Schur complement
127 *
128 */
130{
131 auto asmMap = m_locToGloMap.lock();
132
133 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
134 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
135 int rows = nGlobalBnd - nDirBnd;
136
137 Array<OneD, NekDouble> vOutput(nGlobalBnd, 0.0);
138
139 // Extract diagonal contributions
141 for (unsigned int i = 0; i < rows; ++i)
142 {
143 vOutput[nDirBnd + i] = diagonals[i];
144 }
145
146 // Assemble diagonal contributions across processes
147 asmMap->UniversalAssembleBnd(vOutput);
148
150 Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
151}
152
153/**
154 *
155 */
157 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
158 const bool &IsLocal)
159{
160 auto asmMap = m_locToGloMap.lock();
161
162 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
163
164 bool isFull = solvertype == eIterativeFull ? true : false;
165
166 int nGlobal = (isFull) ? asmMap->GetNumGlobalCoeffs()
167 : asmMap->GetNumGlobalBndCoeffs();
168 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
169 int nNonDir = nGlobal - nDir;
170
171 if (IsLocal)
172 {
173 Array<OneD, NekDouble> wk(nGlobal);
174 (isFull) ? asmMap->Assemble(pInput, wk)
175 : asmMap->AssembleBnd(pInput, wk);
176 Vmath::Vmul(nNonDir, wk.data() + nDir, 1, m_diagonals.data(), 1,
177 wk.data() + nDir, 1);
178 Vmath::Zero(nDir, wk, 1);
179 (isFull) ? asmMap->GlobalToLocal(wk, pOutput)
180 : asmMap->GlobalToLocalBnd(wk, pOutput);
181 }
182 else
183 {
184 Vmath::Vmul(nNonDir, &pInput[0], 1, &m_diagonals[0], 1, &pOutput[0], 1);
185 }
186}
187
190 "Null", PreconditionerNull::create, "No Preconditioning");
191
192/**
193 * @class Null Preconditioner
194 *
195 * This class implements no preconditioning for the conjugate
196 * gradient matrix solver.
197 */
199 const std::shared_ptr<GlobalLinSys> &plinsys,
200 const AssemblyMapSharedPtr &pLocToGloMap)
201 : Preconditioner(plinsys, pLocToGloMap)
202{
203}
204
205/**
206 *
207 */
209{
210}
211
212/**
213 *
214 */
216{
217}
218
219/**
220 *
221 */
223 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
224 const bool &isLocal)
225
226{
227 auto asmMap = m_locToGloMap.lock();
228
229 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
230
231 bool isFull = solvertype == eIterativeFull ? true : false;
232
233 int nGlobal = (isFull) ? asmMap->GetNumGlobalCoeffs()
234 : asmMap->GetNumGlobalBndCoeffs();
235 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
236
237 if (isLocal)
238 {
239 Array<OneD, NekDouble> wk(nGlobal);
240 (isFull) ? asmMap->Assemble(pInput, wk)
241 : asmMap->AssembleBnd(pInput, wk);
242 Vmath::Zero(nDir, wk, 1);
243 (isFull) ? asmMap->GlobalToLocal(wk, pOutput)
244 : asmMap->GlobalToLocalBnd(wk, pOutput);
245 }
246 else
247 {
248 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
249 }
250}
251
252/**
253 * Registers the class with the Factory.
254 */
257 "Jacobi", PreconditionerJacobi::create, "Jacobi Preconditioning");
258/**
259 * @class PreconditionerJacobi
260 *
261 * This class implements jacobi preconditioning building on diagonal version
262 * above
263 */
264
266 const std::shared_ptr<GlobalLinSys> &plinsys,
267 const AssemblyMapSharedPtr &pLocToGloMap)
268 : PreconditionerDiagonal(plinsys, pLocToGloMap)
269{
270}
271
273{
274}
275
277{
279
280 auto expList = ((m_linsys.lock())->GetLocMat()).lock();
281 std::shared_ptr<LibUtilities::SessionReader> session =
282 expList->GetSession();
283
284 std::string var = m_locToGloMap.lock()->GetVariable();
285
286 if (session->DefinesGlobalSysSolnInfo(var, "JacobiIterations"))
287 {
288 m_niter = boost::lexical_cast<int>(
289 session->GetGlobalSysSolnInfo(var, "JacobiIterations").c_str());
290 }
291 else
292 {
293 m_niter = 3;
294 }
295}
296
297/**
298 *
299 */
301 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
302 const bool &IsLocal)
303{
304 auto asmMap = m_locToGloMap.lock();
305
306 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
307
308 int nGlobal = solvertype == eIterativeFull
309 ? asmMap->GetNumGlobalCoeffs()
310 : asmMap->GetNumGlobalBndCoeffs();
311 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
312 int nNonDir = nGlobal - nDir;
313 Array<OneD, NekDouble> wk(nGlobal);
314
315 if (IsLocal)
316 {
317 int nLocal = solvertype == eIterativeFull
318 ? asmMap->GetNumLocalCoeffs()
319 : asmMap->GetNumLocalBndCoeffs();
320
321 Array<OneD, NekDouble> wk1(nLocal);
322 asmMap->Assemble(pInput, wk);
323 Vmath::Vmul(nNonDir, wk.data() + nDir, 1, m_diagonals.data(), 1,
324 wk.data() + nDir, 1);
325 Vmath::Zero(nDir, wk, 1);
326
327 for (int n = 1; n < m_niter; ++n)
328 {
329 asmMap->GlobalToLocal(wk, pOutput);
330
331 // do Ax operator
332 std::dynamic_pointer_cast<GlobalLinSysIterative>(m_linsys.lock())
333 ->DoMatrixMultiply(pOutput, wk1);
334
335 Vmath::Vsub(nLocal, pInput, 1, wk1, 1, wk1, 1);
336
337 asmMap->Assemble(wk1, pOutput);
338 Vmath::Vvtvp(nNonDir, pOutput.data() + nDir, 1, m_diagonals.data(),
339 1, wk.data() + nDir, 1, wk.data() + nDir, 1);
340 }
341
342 asmMap->GlobalToLocal(wk, pOutput);
343 }
344 else
345 {
346 Array<OneD, NekDouble> wk1(nGlobal);
347 Vmath::Vmul(nNonDir, pInput.data(), 1, m_diagonals.data(), 1,
348 wk.data() + nDir, 1);
349 Vmath::Zero(nDir, wk, 1);
350
351 for (int n = 1; n < m_niter; ++n)
352 {
353 // do Ax operator
354 std::dynamic_pointer_cast<GlobalLinSysIterative>(m_linsys.lock())
355 ->DoMatrixMultiply(wk, wk1);
356
357 // b - Ax
358 Vmath::Vsub(nNonDir, pInput.data(), 1, wk1.data() + nDir, 1,
359 wk1.data() + nDir, 1);
360
361 // new sol = 1/diag (b-Ax) + old sol
362 Vmath::Vvtvp(nNonDir, wk1.data() + nDir, 1, m_diagonals.data(), 1,
363 wk.data() + nDir, 1, wk.data() + nDir, 1);
364 }
365
366 Vmath::Vcopy(nNonDir, wk.data() + nDir, 1, pOutput.data(), 1);
367 }
368}
369
370} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconditionerDiagonal(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
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.
const std::weak_ptr< GlobalLinSys > m_linsys
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
std::weak_ptr< AssemblyMap > m_locToGloMap
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconditionerJacobi(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
static std::string className
Name of class.
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.
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.
static std::string className
Name of class.
PreconditionerNull(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
Definition: Vmath.hpp:154
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
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