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 
39 #include <LocalRegions/MatrixKey.h>
40 #include <cmath>
41 
42 namespace Nektar
43 {
44  namespace MultiRegions
45  {
46  std::string Preconditioner::lookupIds[8] = {
48  "Preconditioner", "Null", eNull),
50  "Preconditioner", "Diagonal", eDiagonal),
52  "Preconditioner", "FullLinearSpaceWithDiagonal",
55  "Preconditioner", "FullLinearSpace",eLinear),
57  "Preconditioner", "LowEnergyBlock",eLowEnergy),
59  "Preconditioner", "FullLinearSpaceWithLowEnergyBlock",
62  "Preconditioner", "Block",eBlock),
64  "Preconditioner", "FullLinearSpaceWithBlock",eLinearWithBlock),
65  };
66  std::string Preconditioner::def =
68  "Preconditioner", "Diagonal");
69 
70  /**
71  * @class Preconditioner
72  *
73  * This class implements preconditioning for the conjugate
74  * gradient matrix solver.
75  */
76 
78  const std::shared_ptr<GlobalLinSys> &plinsys,
79  const AssemblyMapSharedPtr &pLocToGloMap)
80  : m_linsys(plinsys),
81  m_preconType(pLocToGloMap->GetPreconType()),
82  m_locToGloMap(pLocToGloMap)
83  {
84  }
85 
86  /**
87  *
88  */
90  {
91  static PreconFactory instance;
92  return instance;
93  }
94 
96  {
97  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
98  }
99 
100  /**
101  * \brief Apply a preconditioner to the conjugate gradient method
102  */
104  const Array<OneD, NekDouble>& pInput,
105  Array<OneD, NekDouble>& pOutput)
106  {
107  boost::ignore_unused(pInput, pOutput);
108  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
109  }
110 
111  /**
112  * \brief Apply a preconditioner to the conjugate gradient method with
113  * an output for non-vertex degrees of freedom.
114  */
116  const Array<OneD, NekDouble> &pInput,
117  Array<OneD, NekDouble> &pOutput,
118  const Array<OneD, NekDouble> &pNonVertOutput,
119  Array<OneD, NekDouble>& pVertForce)
120  {
121  boost::ignore_unused(pInput, pOutput, pNonVertOutput, pVertForce);
122  NEKERROR(ErrorUtil::efatal, "Method does not exist");
123  }
124 
125  /**
126  * \brief Transform from original basis to low energy basis
127  */
129  Array<OneD, NekDouble>& pInOut)
130  {
131  boost::ignore_unused(pInOut);
132  }
133 
134  /**
135  * \brief Transform from low energy coeffs to orignal basis
136  */
138  Array<OneD, NekDouble>& pInOut)
139  {
140  boost::ignore_unused(pInOut);
141  }
142 
143  /**
144  * \brief Multiply by the block inverse transformation matrix
145  */
147  const Array<OneD, NekDouble> &pInput,
148  Array<OneD, NekDouble> &pOutput)
149  {
150  boost::ignore_unused(pInput, pOutput);
151  }
152 
153  /**
154  * \brief Multiply by the block transposed inverse transformation matrix
155  */
157  const Array<OneD, NekDouble> &pInput,
158  Array<OneD, NekDouble> &pOutput)
159  {
160  boost::ignore_unused(pInput, pOutput);
161  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
162  }
163 
165  {
166  }
167 
168  /**
169  * \brief Get block elemental transposed transformation matrix
170  * \f$\mathbf{R}^{T}\f$
171  */
173  int offset, int bnd_offset,
174  const std::shared_ptr<DNekScalMat> &loc_mat)
175  {
176  boost::ignore_unused(offset, bnd_offset);
177  return loc_mat;
178  }
179 
180  /**
181  * @brief Performs global assembly of diagonal entries to global Schur
182  * complement matrix.
183  */
186  {
187  auto asmMap = m_locToGloMap.lock();
188 
189  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
190  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
191  int rows = nGlobalBnd - nDirBnd;
192 
193  DNekScalBlkMatSharedPtr loc_mat;
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)
214  - nDirBnd;
215  sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
216 
217  if (gid1 < 0)
218  {
219  continue;
220  }
221 
222  for (j = 0; j < bnd_row; ++j)
223  {
224  gid2 = asmMap->GetLocalToGlobalBndMap (cnt + j)
225  - nDirBnd;
226  sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
227 
228  if (gid2 == gid1)
229  {
230  diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
231  }
232  }
233  }
234  cnt += bnd_row;
235  }
236 
237  return diagonals;
238  }
239  }
240 }
241 
242 
243 
244 
245 
246 
#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:52
@ eNull
No Solution type specified.
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73