Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PreconditionerDiagonal.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PreconditionerDiagonal.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: PreconditionerDiagonal definition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
39 #include <cmath>
40 
41 using namespace std;
42 
43 namespace Nektar
44 {
45 namespace MultiRegions
46 {
47 /**
48  * Registers the class with the Factory.
49  */
50 string PreconditionerDiagonal::className =
52  "Diagonal", PreconditionerDiagonal::create, "Diagonal Preconditioning");
53 /**
54  * @class Preconditioner
55  *
56  * This class implements diagonal preconditioning for the conjugate
57  * gradient matrix solver.
58  */
59 
60 PreconditionerDiagonal::PreconditionerDiagonal(
61  const std::shared_ptr<GlobalLinSys> &plinsys,
62  const AssemblyMapSharedPtr &pLocToGloMap)
63  : Preconditioner(plinsys, pLocToGloMap)
64 {
65 }
66 
68 {
69 }
70 
72 {
73  GlobalSysSolnType solvertype = m_locToGloMap.lock()->GetGlobalSysSolnType();
74  if (solvertype == eIterativeFull)
75  {
77  }
78  else if (solvertype == eIterativeStaticCond ||
79  solvertype == eIterativeMultiLevelStaticCond ||
80  solvertype == ePETScStaticCond ||
81  solvertype == ePETScMultiLevelStaticCond)
82  {
84  }
85  else
86  {
87  ASSERTL0(0, "Unsupported solver type");
88  }
89 }
90 
91 /**
92  * Diagonal preconditioner computed by summing the relevant elements of
93  * the local matrix system.
94  */
96 {
97  std::shared_ptr<MultiRegions::ExpList> expList =
98  ((m_linsys.lock())->GetLocMat()).lock();
99 
101 
102  auto asmMap = m_locToGloMap.lock();
103 
104  int i, j, n, cnt, gid1, gid2;
105  NekDouble sign1, sign2, value;
106  int nGlobal = asmMap->GetNumGlobalCoeffs();
107  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
108  int nInt = nGlobal - nDir;
109 
110  // fill global matrix
111  DNekScalMatSharedPtr loc_mat;
112  Array<OneD, NekDouble> vOutput(nGlobal, 0.0);
113 
114  int loc_row;
115  int nElmt = expList->GetNumElmts();
116  for (n = cnt = 0; n < nElmt; ++n)
117  {
118  loc_mat = (m_linsys.lock())->GetBlock(n);
119  loc_row = loc_mat->GetRows();
120 
121  for (i = 0; i < loc_row; ++i)
122  {
123  gid1 = asmMap->GetLocalToGlobalMap(cnt + i) - nDir;
124  sign1 = asmMap->GetLocalToGlobalSign(cnt + i);
125 
126  if (gid1 >= 0)
127  {
128  for (j = 0; j < loc_row; ++j)
129  {
130  gid2 = asmMap->GetLocalToGlobalMap(cnt + j) - nDir;
131  sign2 = asmMap->GetLocalToGlobalSign(cnt + j);
132  if (gid2 == gid1)
133  {
134  // When global matrix is symmetric,
135  // only add the value for the upper
136  // triangular part in order to avoid
137  // entries to be entered twice
138  value = vOutput[gid1 + nDir] +
139  sign1 * sign2 * (*loc_mat)(i, j);
140  vOutput[gid1 + nDir] = value;
141  }
142  }
143  }
144  }
145  cnt += loc_row;
146  }
147 
148  // Assemble diagonal contributions across processes
149  asmMap->UniversalAssemble(vOutput);
150 
152  Vmath::Sdiv(nInt, 1.0, &vOutput[nDir], 1, &m_diagonals[0], 1);
153 }
154 
155 /**
156  * Diagonal preconditioner defined as the inverse of the main
157  * diagonal of the Schur complement
158  *
159  */
161 {
162  auto asmMap = m_locToGloMap.lock();
163 
164  int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
165  int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
166  int rows = nGlobalBnd - nDirBnd;
167 
168  Array<OneD, NekDouble> vOutput(nGlobalBnd, 0.0);
169 
170  // Extract diagonal contributions
172  for (unsigned int i = 0; i < rows; ++i)
173  {
174  vOutput[nDirBnd + i] = diagonals[i];
175  }
176 
177  // Assemble diagonal contributions across processes
178  asmMap->UniversalAssembleBnd(vOutput);
179 
181  Vmath::Sdiv(rows, 1.0, &vOutput[nDirBnd], 1, &m_diagonals[0], 1);
182 }
183 
184 /**
185  *
186  */
188  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
189 {
190  auto asmMap = m_locToGloMap.lock();
191 
192  GlobalSysSolnType solvertype = asmMap->GetGlobalSysSolnType();
193 
194  int nGlobal = solvertype == eIterativeFull
195  ? asmMap->GetNumGlobalCoeffs()
196  : asmMap->GetNumGlobalBndCoeffs();
197  int nDir = asmMap->GetNumGlobalDirBndCoeffs();
198  int nNonDir = nGlobal - nDir;
199  Vmath::Vmul(nNonDir, &pInput[0], 1, &m_diagonals[0], 1, &pOutput[0], 1);
200 }
201 
204  "Null", PreconditionerNull::create, "No Preconditioning");
205 
206 /**
207  * @class Null Preconditioner
208  *
209  * This class implements no preconditioning for the conjugate
210  * gradient matrix solver.
211  */
213  const std::shared_ptr<GlobalLinSys> &plinsys,
214  const AssemblyMapSharedPtr &pLocToGloMap)
215  : Preconditioner(plinsys, pLocToGloMap)
216 {
217 }
218 
219 /**
220  *
221  */
223 {
224 }
225 
226 /**
227  *
228  */
230 {
231 }
232 
233 /**
234  *
235  */
237  const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput)
238 {
239  Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
240 }
241 
242 } // namespace MultiRegions
243 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply a preconditioner to the conjugate gradient method.
const std::weak_ptr< GlobalLinSys > m_linsys
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
std::weak_ptr< AssemblyMap > m_locToGloMap
static std::string className
Name of class.
PreconditionerNull(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply a preconditioner to the conjugate gradient method.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:51
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:324
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255