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 <math.h>
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  int offset)
131  {
132  boost::ignore_unused(pInOut, offset);
133  }
134 
135  /**
136  * \brief Transform from original basis to low energy basis
137  */
139  const Array<OneD, NekDouble> &pInOut,
140  Array<OneD, NekDouble> &pOutput)
141  {
142  boost::ignore_unused(pInOut, pOutput);
143  }
144 
145  /**
146  * \brief Transform from low energy basis to orignal basis
147  */
149  Array<OneD, NekDouble>& pInput)
150  {
151  Vmath::Smul(pInput.num_elements(), 1.0, pInput, 1, pInput, 1);
152  }
153 
154  /**
155  * \brief Multiply by the block inverse transformation matrix
156  */
158  const Array<OneD, NekDouble> &pInput,
159  Array<OneD, NekDouble> &pOutput)
160  {
161  boost::ignore_unused(pInput, pOutput);
162  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
163  }
164 
165  /**
166  * \brief Multiply by the block transposed inverse transformation matrix
167  */
169  const Array<OneD, NekDouble> &pInput,
170  Array<OneD, NekDouble> &pOutput)
171  {
172  boost::ignore_unused(pInput, pOutput);
173  NEKERROR(ErrorUtil::efatal,"Method does not exist" );
174  }
175 
177  {
178  }
179 
180  /**
181  * \brief Get block elemental transposed transformation matrix
182  * \f$\mathbf{R}^{T}\f$
183  */
185  int offset, int bnd_offset,
186  const std::shared_ptr<DNekScalMat> &loc_mat)
187  {
188  boost::ignore_unused(offset, bnd_offset);
189  return loc_mat;
190  }
191 
192  /**
193  * @brief Performs global assembly of diagonal entries to global Schur
194  * complement matrix.
195  */
198  {
199  auto asmMap = m_locToGloMap.lock();
200 
201  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
202  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
203  int rows = nGlobalBnd - nDirBnd;
204 
205  DNekScalBlkMatSharedPtr loc_mat;
206  DNekScalMatSharedPtr bnd_mat;
207  int sign1, sign2, gid1, gid2, i, j, n, cnt;
208  Array<OneD, NekDouble> diagonals(rows,0.0);
209 
210  // Extract diagonal contributions of globally assembled
211  // schur complement matrix
212  for (cnt = n = 0; n < m_linsys.lock()->GetNumBlocks(); ++n)
213  {
214  // Get statically condensed local matrix.
215  loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
216 
217  // Extract boundary block.
218  bnd_mat = loc_mat->GetBlock(0, 0);
219 
220  // Offset by number of rows.
221  int bnd_row = bnd_mat->GetRows();
222 
223  for (i = 0; i < bnd_row; ++i)
224  {
225  gid1 = asmMap->GetLocalToGlobalBndMap (cnt + i)
226  - nDirBnd;
227  sign1 = asmMap->GetLocalToGlobalBndSign(cnt + i);
228 
229  if (gid1 < 0)
230  {
231  continue;
232  }
233 
234  for (j = 0; j < bnd_row; ++j)
235  {
236  gid2 = asmMap->GetLocalToGlobalBndMap (cnt + j)
237  - nDirBnd;
238  sign2 = asmMap->GetLocalToGlobalBndSign(cnt + j);
239 
240  if (gid2 == gid1)
241  {
242  diagonals[gid1] += sign1 * sign2 * (*bnd_mat)(i, j);
243  }
244  }
245  }
246  cnt += bnd_row;
247  }
248 
249  return diagonals;
250  }
251  }
252 }
253 
254 
255 
256 
257 
258 
Preconditioner(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual DNekScalMatSharedPtr v_TransformedSchurCompl(int offset, int bndoffset, const std::shared_ptr< DNekScalMat > &loc_mat)
Get block elemental transposed transformation matrix .
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
virtual void v_DoMultiplybyInverseTransposedTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block transposed inverse transformation matrix.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:73
const std::weak_ptr< GlobalLinSys > m_linsys
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Apply a preconditioner to the conjugate gradient method.
PreconFactory & GetPreconFactory()
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...
virtual void v_DoMultiplybyInverseTransformationMatrix(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput)
Multiply by the block inverse transformation matrix.
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
virtual void v_DoTransformFromLowEnergy(Array< OneD, NekDouble > &pInput)
Transform from low energy basis to orignal basis.
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
No Solution type specified.
virtual void v_DoTransformToLowEnergy(Array< OneD, NekDouble > &pInOut, int offset)
Transform from original basis to low energy basis.
static std::string RegisterDefaultSolverInfo(const std::string &pName, const std::string &pValue)
Registers the default string value of a solver info property.
std::weak_ptr< AssemblyMap > m_locToGloMap
Provides a generic Factory class.
Definition: NekFactory.hpp:103