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 
37 #include <LocalRegions/MatrixKey.h>
40 #include <cmath>
41 
42 namespace Nektar
43 {
44 namespace MultiRegions
45 {
46 std::string Preconditioner::lookupIds[8] = {
47  LibUtilities::SessionReader::RegisterEnumValue("Preconditioner", "Null",
48  eNull),
49  LibUtilities::SessionReader::RegisterEnumValue("Preconditioner", "Diagonal",
50  eDiagonal),
52  "Preconditioner", "FullLinearSpaceWithDiagonal", eLinearWithDiagonal),
54  "FullLinearSpace", eLinear),
56  "Preconditioner", "LowEnergyBlock", eLowEnergy),
58  "Preconditioner", "FullLinearSpaceWithLowEnergyBlock",
60  LibUtilities::SessionReader::RegisterEnumValue("Preconditioner", "Block",
61  eBlock),
63  "Preconditioner", "FullLinearSpaceWithBlock", eLinearWithBlock),
64 };
65 std::string Preconditioner::def =
67  "Diagonal");
68 
69 /**
70  * @class Preconditioner
71  *
72  * This class implements preconditioning for the conjugate
73  * gradient matrix solver.
74  */
75 
76 Preconditioner::Preconditioner(const std::shared_ptr<GlobalLinSys> &plinsys,
77  const AssemblyMapSharedPtr &pLocToGloMap)
78  : m_linsys(plinsys), m_preconType(pLocToGloMap->GetPreconType()),
79  m_locToGloMap(pLocToGloMap)
80 {
81 }
82 
83 /**
84  *
85  */
87 {
88  static PreconFactory instance;
89  return instance;
90 }
91 
93 {
94  NEKERROR(ErrorUtil::efatal, "Method does not exist");
95 }
96 
97 /**
98  * \brief Apply a preconditioner to the conjugate gradient method
99  */
101  Array<OneD, NekDouble> &pOutput)
102 {
103  boost::ignore_unused(pInput, pOutput);
104  NEKERROR(ErrorUtil::efatal, "Method does not exist");
105 }
106 
107 /**
108  * \brief Apply a preconditioner to the conjugate gradient method with
109  * an output for non-vertex degrees of freedom.
110  */
112  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
113  const Array<OneD, NekDouble> &pNonVertOutput,
114  Array<OneD, NekDouble> &pVertForce)
115 {
116  boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
117  NEKERROR(ErrorUtil::efatal, "Method does not exist");
118 }
119 
120 /**
121  * \brief Transform from original basis to low energy basis
122  */
124  Array<OneD, NekDouble> &pInOut)
125 {
126  boost::ignore_unused(pInOut);
127 }
128 
129 /**
130  * \brief Transform from low energy coeffs to orignal basis
131  */
133  Array<OneD, NekDouble> &pInOut)
134 {
135  boost::ignore_unused(pInOut);
136 }
137 
138 /**
139  * \brief Multiply by the block inverse transformation matrix
140  */
142  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
143 {
144  boost::ignore_unused(pInput, pOutput);
145 }
146 
147 /**
148  * \brief Multiply by the block transposed inverse transformation matrix
149  */
151  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
152 {
153  boost::ignore_unused(pInput, pOutput);
154  NEKERROR(ErrorUtil::efatal, "Method does not exist");
155 }
156 
158 {
159 }
160 
161 /**
162  * \brief Get block elemental transposed transformation matrix
163  * \f$\mathbf{R}^{T}\f$
164  */
166  int offset, int bnd_offset, const std::shared_ptr<DNekScalMat> &loc_mat)
167 {
168  boost::ignore_unused(offset, bnd_offset);
169  return loc_mat;
170 }
171 
172 /**
173  * @brief Performs global assembly of diagonal entries to global Schur
174  * complement matrix.
175  */
177 {
178  auto asmMap = m_locToGloMap.lock();
179 
180  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
181  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
182  int rows = nGlobalBnd - nDirBnd;
183 
184  DNekScalBlkMatSharedPtr loc_mat;
185  DNekScalMatSharedPtr bnd_mat;
186  int sign1, sign2, gid1, gid2, i, j, n, cnt;
187  Array<OneD, NekDouble> diagonals(rows, 0.0);
188 
189  // Extract diagonal contributions of globally assembled
190  // schur complement matrix
191  for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
192  {
193  // Get statically condensed local matrix.
194  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
195 
196  // Extract boundary block.
197  bnd_mat = loc_mat->GetBlock(0, 0);
198 
199  // Offset by number of rows.
200  int bnd_row = bnd_mat->GetRows();
201 
202  for (i = 0; i < bnd_row; ++i)
203  {
204  gid1 = asmMap->GetLocalToGlobalBndMap(cnt + i) - nDirBnd;
205  sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
206 
207  if (gid1 < 0)
208  {
209  continue;
210  }
211 
212  for (j = 0; j < bnd_row; ++j)
213  {
214  gid2 = asmMap->GetLocalToGlobalBndMap(cnt + j) - nDirBnd;
215  sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
216 
217  if (gid2 == gid1)
218  {
219  diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
220  }
221  }
222  }
223  cnt += bnd_row;
224  }
225 
226  return diagonals;
227 }
228 } // namespace MultiRegions
229 } // 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 RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
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_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.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.
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:51
@ eNull
No Solution type specified.
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