Nektar++
Preconditioner.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Preconditioner.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: Preconditioner definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
38#include <cmath>
39
41{
42
43// register default solver value as diagonal
44std::string Preconditioner::def =
46 "Diagonal");
47/**
48 * @class Preconditioner
49 *
50 * This class implements preconditioning for the conjugate
51 * gradient matrix solver.
52 */
53
54Preconditioner::Preconditioner(const std::shared_ptr<GlobalLinSys> &plinsys,
55 const AssemblyMapSharedPtr &pLocToGloMap)
56 : m_linsys(plinsys), m_preconType(pLocToGloMap->GetPreconType()),
57 m_locToGloMap(pLocToGloMap)
58{
59}
60
61/**
62 *
63 */
65{
66 static PreconFactory instance;
67 return instance;
68}
69
71{
72 NEKERROR(ErrorUtil::efatal, "Method does not exist");
73}
74
75/**
76 * \brief Apply a preconditioner to the conjugate gradient method
77 */
79 [[maybe_unused]] const Array<OneD, NekDouble> &pInput,
80 [[maybe_unused]] Array<OneD, NekDouble> &pOutput,
81 [[maybe_unused]] const bool &IsLocal)
82{
83 NEKERROR(ErrorUtil::efatal, "Method does not exist");
84}
85
86/**
87 * \brief Apply an assembly and scatter back to lcoal array
88 */
91 const bool &ZeroDir)
92{
93 auto assMap = m_locToGloMap.lock();
94 GlobalSysSolnType solvertype = assMap->GetGlobalSysSolnType();
95
96 if (solvertype == eIterativeFull)
97 {
98 assMap->Assemble(pInput, pOutput);
99 if (ZeroDir)
100 {
101 int nDir = assMap->GetNumGlobalDirBndCoeffs();
102 Vmath::Zero(nDir, pOutput, 1);
103 }
104 assMap->GlobalToLocal(pOutput, pOutput);
105 }
106 else // bnd version.
107 {
108 assMap->AssembleBnd(pInput, pOutput);
109 if (ZeroDir)
110 {
111 int nDir = assMap->GetNumGlobalDirBndCoeffs();
112 Vmath::Zero(nDir, pOutput, 1);
113 }
114 assMap->GlobalToLocalBnd(pOutput, pOutput);
115 }
116}
117
118/**
119 * \brief Apply a preconditioner to the conjugate gradient method with
120 * an output for non-vertex degrees of freedom.
121 */
123 [[maybe_unused]] const Array<OneD, NekDouble> &pInput,
124 [[maybe_unused]] Array<OneD, NekDouble> &pOutput,
125 [[maybe_unused]] const Array<OneD, NekDouble> &pNonVertOutput,
126 [[maybe_unused]] Array<OneD, NekDouble> &pVertForce)
127{
128 NEKERROR(ErrorUtil::efatal, "Method does not exist");
129}
130
131/**
132 * \brief Transform from original basis to low energy basis
133 */
135 [[maybe_unused]] Array<OneD, NekDouble> &pInOut)
136{
137}
138
139/**
140 * \brief Transform from low energy coeffs to orignal basis
141 */
143 [[maybe_unused]] Array<OneD, NekDouble> &pInOut)
144{
145}
146
147/**
148 * \brief Multiply by the block inverse transformation matrix
149 */
151 [[maybe_unused]] const Array<OneD, NekDouble> &pInput,
152 [[maybe_unused]] Array<OneD, NekDouble> &pOutput)
153{
154}
155
156/**
157 * \brief Multiply by the block transposed inverse transformation matrix
158 */
160 [[maybe_unused]] const Array<OneD, NekDouble> &pInput,
161 [[maybe_unused]] Array<OneD, NekDouble> &pOutput)
162{
163 NEKERROR(ErrorUtil::efatal, "Method does not exist");
164}
165
167{
168}
169
170/**
171 * \brief Get block elemental transposed transformation matrix
172 * \f$\mathbf{R}^{T}\f$
173 */
175 [[maybe_unused]] int offset, [[maybe_unused]] int bnd_offset,
176 const std::shared_ptr<DNekScalMat> &loc_mat)
177{
178 return loc_mat;
179}
180
181/**
182 * @brief Performs global assembly of diagonal entries to global Schur
183 * complement matrix.
184 */
186{
187 auto asmMap = m_locToGloMap.lock();
188
189 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
190 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
191 int rows = nGlobalBnd - nDirBnd;
192
194 DNekScalMatSharedPtr bnd_mat;
195 int sign1, sign2, gid1, gid2, i, j, n, cnt;
196 Array<OneD, NekDouble> diagonals(rows, 0.0);
197
198 // Extract diagonal contributions of globally assembled
199 // schur complement matrix
200 for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
201 {
202 // Get statically condensed local matrix.
203 loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
204
205 // Extract boundary block.
206 bnd_mat = loc_mat->GetBlock(0, 0);
207
208 // Offset by number of rows.
209 int bnd_row = bnd_mat->GetRows();
210
211 for (i = 0; i < bnd_row; ++i)
212 {
213 gid1 = asmMap->GetLocalToGlobalBndMap(cnt + i) - nDirBnd;
214 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
215
216 if (gid1 < 0)
217 {
218 continue;
219 }
220
221 for (j = 0; j < bnd_row; ++j)
222 {
223 gid2 = asmMap->GetLocalToGlobalBndMap(cnt + j) - nDirBnd;
224 sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
225
226 if (gid2 == gid1)
227 {
228 diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
229 }
230 }
231 }
232 cnt += bnd_row;
233 }
234
235 return diagonals;
236}
237} // namespace Nektar::MultiRegions
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
Provides a generic Factory class.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false)
Apply a preconditioner to the conjugate gradient method.
virtual void v_DoTransformBasisToLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from original basis to low energy basis.
const std::weak_ptr< GlobalLinSys > m_linsys
virtual void v_DoPreconditionerWithNonVertOutput(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const Array< OneD, NekDouble > &pNonVertOutput, Array< OneD, NekDouble > &pVertForce)
Apply a preconditioner to the conjugate gradient method with an output for non-vertex degrees of free...
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
void DoAssembleLoc(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &ZeroDir)
Apply an assembly and scatter back to lcoal array.
std::weak_ptr< AssemblyMap > m_locToGloMap
virtual void v_DoTransformBasisFromLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.
Preconditioner(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual void v_DoTransformCoeffsFromLowEnergy(Array< OneD, NekDouble > &pInOut)
Transform from low energy coeffs to orignal basis.
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
Get block elemental transposed transformation matrix .
virtual void v_DoTransformCoeffsToLowEnergy(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273