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
35#include <boost/core/ignore_unused.hpp>
36
40#include <cmath>
41
42namespace Nektar
43{
44namespace MultiRegions
45{
46
47// register default solver value as diagonal
48std::string Preconditioner::def =
50 "Diagonal");
51/**
52 * @class Preconditioner
53 *
54 * This class implements preconditioning for the conjugate
55 * gradient matrix solver.
56 */
57
58Preconditioner::Preconditioner(const std::shared_ptr<GlobalLinSys> &plinsys,
59 const AssemblyMapSharedPtr &pLocToGloMap)
60 : m_linsys(plinsys), m_preconType(pLocToGloMap->GetPreconType()),
61 m_locToGloMap(pLocToGloMap)
62{
63}
64
65/**
66 *
67 */
69{
70 static PreconFactory instance;
71 return instance;
72}
73
75{
76 NEKERROR(ErrorUtil::efatal, "Method does not exist");
77}
78
79/**
80 * \brief Apply a preconditioner to the conjugate gradient method
81 */
84 const bool &IsLocal)
85{
86 boost::ignore_unused(pInput, pOutput, IsLocal);
87 NEKERROR(ErrorUtil::efatal, "Method does not exist");
88}
89
90/**
91 * \brief Apply an assembly and scatter back to lcoal array
92 */
95 const bool &ZeroDir)
96{
97 auto assMap = m_locToGloMap.lock();
98 GlobalSysSolnType solvertype = assMap->GetGlobalSysSolnType();
99
100 if (solvertype == eIterativeFull)
101 {
102 assMap->Assemble(pInput, pOutput);
103 if (ZeroDir)
104 {
105 int nDir = assMap->GetNumGlobalDirBndCoeffs();
106 Vmath::Zero(nDir, pOutput, 1);
107 }
108 assMap->GlobalToLocal(pOutput, pOutput);
109 }
110 else // bnd version.
111 {
112 assMap->AssembleBnd(pInput, pOutput);
113 if (ZeroDir)
114 {
115 int nDir = assMap->GetNumGlobalDirBndCoeffs();
116 Vmath::Zero(nDir, pOutput, 1);
117 }
118 assMap->GlobalToLocalBnd(pOutput, pOutput);
119 }
120}
121
122/**
123 * \brief Apply a preconditioner to the conjugate gradient method with
124 * an output for non-vertex degrees of freedom.
125 */
127 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
128 const Array<OneD, NekDouble> &pNonVertOutput,
129 Array<OneD, NekDouble> &pVertForce)
130{
131 boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
132 NEKERROR(ErrorUtil::efatal, "Method does not exist");
133}
134
135/**
136 * \brief Transform from original basis to low energy basis
137 */
140{
141 boost::ignore_unused(pInOut);
142}
143
144/**
145 * \brief Transform from low energy coeffs to orignal basis
146 */
149{
150 boost::ignore_unused(pInOut);
151}
152
153/**
154 * \brief Multiply by the block inverse transformation matrix
155 */
157 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
158{
159 boost::ignore_unused(pInput, pOutput);
160}
161
162/**
163 * \brief Multiply by the block transposed inverse transformation matrix
164 */
166 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
167{
168 boost::ignore_unused(pInput, pOutput);
169 NEKERROR(ErrorUtil::efatal, "Method does not exist");
170}
171
173{
174}
175
176/**
177 * \brief Get block elemental transposed transformation matrix
178 * \f$\mathbf{R}^{T}\f$
179 */
181 int offset, int bnd_offset, const std::shared_ptr<DNekScalMat> &loc_mat)
182{
183 boost::ignore_unused(offset, bnd_offset);
184 return loc_mat;
185}
186
187/**
188 * @brief Performs global assembly of diagonal entries to global Schur
189 * complement matrix.
190 */
192{
193 auto asmMap = m_locToGloMap.lock();
194
195 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
196 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
197 int rows = nGlobalBnd - nDirBnd;
198
200 DNekScalMatSharedPtr bnd_mat;
201 int sign1, sign2, gid1, gid2, i, j, n, cnt;
202 Array<OneD, NekDouble> diagonals(rows, 0.0);
203
204 // Extract diagonal contributions of globally assembled
205 // schur complement matrix
206 for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
207 {
208 // Get statically condensed local matrix.
209 loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
210
211 // Extract boundary block.
212 bnd_mat = loc_mat->GetBlock(0, 0);
213
214 // Offset by number of rows.
215 int bnd_row = bnd_mat->GetRows();
216
217 for (i = 0; i < bnd_row; ++i)
218 {
219 gid1 = asmMap->GetLocalToGlobalBndMap(cnt + i) - nDirBnd;
220 sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
221
222 if (gid1 < 0)
223 {
224 continue;
225 }
226
227 for (j = 0; j < bnd_row; ++j)
228 {
229 gid2 = asmMap->GetLocalToGlobalBndMap(cnt + j) - nDirBnd;
230 sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
231
232 if (gid2 == gid1)
233 {
234 diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
235 }
236 }
237 }
238 cnt += bnd_row;
239 }
240
241 return diagonals;
242}
243} // namespace MultiRegions
244} // namespace Nektar
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
Provides a generic Factory class.
Definition: NekFactory.hpp:105
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:52
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:487