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
43namespace Nektar
44{
45namespace MultiRegions
46{
47/**
48 * Registers the class with the Factory.
49 */
52 "Diagonal", PreconditionerDiagonal::create, "Diagonal Preconditioning");
53/**
54 * @class Preconditioner
55 *
56 * This class implements diagonal preconditioning for the conjugate
57 * gradient matrix solver.
58 */
59
61 const std::shared_ptr<GlobalLinSys> &plinsys,
62 const AssemblyMapSharedPtr &pLocToGloMap)
63 : Preconditioner(plinsys, pLocToGloMap)
64{
65}
66
68{
69}
70
72{
73 GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
74 if (solvertype == eIterativeFull)
75 {
77 }
78 else if (solvertype == eIterativeStaticCond ||
79 solvertype == eIterativeMultiLevelStaticCond ||
80 solvertype == ePETScStaticCond ||
81 solvertype == ePETScMultiLevelStaticCond)
82 {
84 }
85 else
86 {
87 ASSERTL0(0, "Unsupported solver type");
88 }
89}
90
91/**
92 * Diagonal preconditioner computed by summing the relevant elements of
93 * the local matrix system.
94 */
96{
97 std::shared_ptr<MultiRegions::ExpList> expList =
98 ((m_linsys.lock())->GetLocMat()).lock();
99
101
102 auto asmMap = m_locToGloMap.lock();
103
104 int i, j, n, cnt, gid1, gid2;
105 NekDouble sign1, sign2, value;
106 int nGlobal = asmMap->GetNumGlobalCoeffs();
107 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
108 int nInt = nGlobal - nDir;
109
110 // fill global matrix
111 DNekScalMatSharedPtr loc_mat;
112 Array<OneD, NekDouble> vOutput(nGlobal, 0.0);
113
114 int loc_row;
115 int nElmt = expList->GetNumElmts();
116 for (n = cnt = 0; n < nElmt; ++n)
117 {
118 loc_mat = (m_linsys.lock())->GetBlock(n);
119 loc_row = loc_mat->GetRows();
120
121 for (i = 0; i < loc_row; ++i)
122 {
123 gid1 = asmMap->GetLocalToGlobalMap(cnt + i) - nDir;
124 sign1 = asmMap->GetLocalToGlobalSign(cnt + i);
125
126 if (gid1 >= 0)
127 {
128 for (j = 0; j < loc_row; ++j)
129 {
130 gid2 = asmMap->GetLocalToGlobalMap(cnt + j) - nDir;
131 sign2 = asmMap->GetLocalToGlobalSign(cnt + j);
132 if (gid2 == gid1)
133 {
134 // When global matrix is symmetric,
135 // only add the value for the upper
136 // triangular part in order to avoid
137 // entries to be entered twice
138 value = vOutput[gid1 + nDir] +
139 sign1 * sign2 * (*loc_mat)(i, j);
140 vOutput[gid1 + nDir] = value;
141 }
142 }
143 }
144 }
145 cnt += loc_row;
146 }
147
148 // Assemble diagonal contributions across processes
149 asmMap->UniversalAssemble(vOutput);
150
152 Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
153}
154
155/**
156 * Diagonal preconditioner defined as the inverse of the main
157 * diagonal of the Schur complement
158 *
159 */
161{
162 auto asmMap = m_locToGloMap.lock();
163
164 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
165 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
166 int rows = nGlobalBnd - nDirBnd;
167
168 Array<OneD, NekDouble> vOutput(nGlobalBnd, 0.0);
169
170 // Extract diagonal contributions
172 for (unsigned int i = 0; i < rows; ++i)
173 {
174 vOutput[nDirBnd + i] = diagonals[i];
175 }
176
177 // Assemble diagonal contributions across processes
178 asmMap->UniversalAssembleBnd(vOutput);
179
181 Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
182}
183
184/**
185 *
186 */
188 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
189 const bool &IsLocal)
190{
191 auto asmMap = m_locToGloMap.lock();
192
193 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
194
195 bool isFull = solvertype == eIterativeFull ? true : false;
196
197 int nGlobal = (isFull) ? asmMap->GetNumGlobalCoeffs()
198 : asmMap->GetNumGlobalBndCoeffs();
199 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
200 int nNonDir = nGlobal - nDir;
201
202 if (IsLocal)
203 {
204 Array<OneD, NekDouble> wk(nGlobal);
205 (isFull) ? asmMap->Assemble(pInput, wk)
206 : asmMap->AssembleBnd(pInput, wk);
207 Vmath::Vmul(nNonDir, wk.data() + nDir, 1, m_diagonals.data(), 1,
208 wk.data() + nDir, 1);
209 Vmath::Zero(nDir, wk, 1);
210 (isFull) ? asmMap->GlobalToLocal(wk, pOutput)
211 : asmMap->GlobalToLocalBnd(wk, pOutput);
212 }
213 else
214 {
215 Vmath::Vmul(nNonDir, &pInput[0], 1, &m_diagonals[0], 1, &pOutput[0], 1);
216 }
217}
218
221 "Null", PreconditionerNull::create, "No Preconditioning");
222
223/**
224 * @class Null Preconditioner
225 *
226 * This class implements no preconditioning for the conjugate
227 * gradient matrix solver.
228 */
230 const std::shared_ptr<GlobalLinSys> &plinsys,
231 const AssemblyMapSharedPtr &pLocToGloMap)
232 : Preconditioner(plinsys, pLocToGloMap)
233{
234}
235
236/**
237 *
238 */
240{
241}
242
243/**
244 *
245 */
247{
248}
249
250/**
251 *
252 */
254 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
255 const bool &isLocal)
256
257{
258 boost::ignore_unused(isLocal);
259 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
260}
261
262/**
263 * Registers the class with the Factory.
264 */
267 "Jacobi", PreconditionerJacobi::create, "Jacobi Preconditioning");
268/**
269 * @class PreconditionerJacobi
270 *
271 * This class implements jacobi preconditioning building on diagonal version
272 * above
273 */
274
276 const std::shared_ptr<GlobalLinSys> &plinsys,
277 const AssemblyMapSharedPtr &pLocToGloMap)
278 : PreconditionerDiagonal(plinsys, pLocToGloMap)
279{
280}
281
283{
284}
285
287{
289
290 auto expList = ((m_linsys.lock())->GetLocMat()).lock();
291 std::shared_ptr<LibUtilities::SessionReader> session =
292 expList->GetSession();
293
294 std::string var = m_locToGloMap.lock()->GetVariable();
295
296 if (session->DefinesGlobalSysSolnInfo(var, "JacobiIterations"))
297 {
298 m_niter = boost::lexical_cast<int>(
299 session->GetGlobalSysSolnInfo(var, "JacobiIterations").c_str());
300 }
301 else
302 {
303 m_niter = 3;
304 }
305}
306
307/**
308 *
309 */
311 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
312 const bool &IsLocal)
313{
314 auto asmMap = m_locToGloMap.lock();
315
316 GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
317
318 int nGlobal = solvertype == eIterativeFull
319 ? asmMap->GetNumGlobalCoeffs()
320 : asmMap->GetNumGlobalBndCoeffs();
321 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
322 int nNonDir = nGlobal - nDir;
323 Array<OneD, NekDouble> wk(nGlobal);
324
325 if (IsLocal)
326 {
327 int nLocal = solvertype == eIterativeFull
328 ? asmMap->GetNumLocalCoeffs()
329 : asmMap->GetNumLocalBndCoeffs();
330
331 Array<OneD, NekDouble> wk1(nLocal);
332 asmMap->Assemble(pInput, wk);
333 Vmath::Vmul(nNonDir, wk.data() + nDir, 1, m_diagonals.data(), 1,
334 wk.data() + nDir, 1);
335 Vmath::Zero(nDir, wk, 1);
336
337 for (int n = 1; n < m_niter; ++n)
338 {
339 asmMap->GlobalToLocal(wk, pOutput);
340
341 // do Ax operator
342 std::dynamic_pointer_cast<GlobalLinSysIterative>(m_linsys.lock())
343 ->DoMatrixMultiply(pOutput, wk1);
344
345 Vmath::Vsub(nLocal, pInput, 1, wk1, 1, wk1, 1);
346
347 asmMap->Assemble(wk1, pOutput);
348 Vmath::Vvtvp(nNonDir, pOutput.data() + nDir, 1, m_diagonals.data(),
349 1, wk.data() + nDir, 1, wk.data() + nDir, 1);
350 }
351
352 asmMap->GlobalToLocal(wk, pOutput);
353 }
354 else
355 {
356 Array<OneD, NekDouble> wk1(nGlobal);
357 Vmath::Vmul(nNonDir, pInput.data(), 1, m_diagonals.data(), 1,
358 wk.data() + nDir, 1);
359 Vmath::Zero(nDir, wk, 1);
360
361 for (int n = 1; n < m_niter; ++n)
362 {
363 // do Ax operator
364 std::dynamic_pointer_cast<GlobalLinSysIterative>(m_linsys.lock())
365 ->DoMatrixMultiply(wk, wk1);
366
367 // b - Ax
368 Vmath::Vsub(nNonDir, pInput.data(), 1, wk1.data() + nDir, 1,
369 wk1.data() + nDir, 1);
370
371 // new sol = 1/diag (b-Ax) + old sol
372 Vmath::Vvtvp(nNonDir, wk1.data() + nDir, 1, m_diagonals.data(), 1,
373 wk.data() + nDir, 1, wk.data() + nDir, 1);
374 }
375
376 Vmath::Vcopy(nNonDir, wk.data() + nDir, 1, pOutput.data(), 1);
377 }
378}
379
380} // namespace MultiRegions
381} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
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)
virtual 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.
virtual 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.
virtual 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.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
double NekDouble
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.cpp:207
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.cpp:569
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.cpp:319
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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