Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PreconCfsBRJ.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PreconCfsBRJ.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: PreconCfsBRJ definition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42 /**
43  * @class PreconCfsBRJ
44  *
45  * Solves a linear system using iterative methods.
46  */
47 std::string PreconCfsBRJ::className =
49  "PreconCfsBRJ", PreconCfsBRJ::create,
50  "Block Relaxed Jacobi Preconditioner for CFS.");
51 
52 PreconCfsBRJ::PreconCfsBRJ(
55  const LibUtilities::CommSharedPtr &vComm)
56  : PreconCfsOp(pFields, pSession, vComm)
57 {
58  pSession->LoadParameter("PreconItsStep", m_PreconItsStep, 7);
59  pSession->LoadParameter("BRJRelaxParam", m_BRJRelaxParam, 1.0);
60 
61  int nvariables = pFields.size();
63 
65  for (int i = 0; i < nvariables; i++)
66  {
68  }
70 
71  int nelmts = pFields[0]->GetNumElmts();
72  int nelmtcoef;
73  Array<OneD, unsigned int> nelmtmatdim(nelmts);
74  for (int i = 0; i < nelmts; i++)
75  {
76  nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
77  nelmtmatdim[i] = nelmtcoef * nvariables;
78  }
79  AllocateNekBlkMatDig(m_PreconMatSingle, nelmtmatdim, nelmtmatdim);
80 }
81 
83 {
85 }
86 
89  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
90  const bool &flag)
91 {
92  boost::ignore_unused(flag);
93 
94  int nBRJIterTot = m_PreconItsStep;
95  if (0 == nBRJIterTot)
96  {
97  DoNullPrecon(inarray, outarray, flag);
98  }
99  else
100  {
101  const NekDouble BRJParam = m_BRJRelaxParam;
102  const NekDouble OmBRJParam = 1.0 - BRJParam;
103 
104  unsigned int nvariables = pFields.size();
105  unsigned int npoints = pFields[0]->GetNcoeffs();
106  unsigned int ntotpnt = inarray.size();
107 
108  ASSERTL0(nvariables * npoints == ntotpnt,
109  "nvariables*npoints!=ntotpnt in PreconCoeff");
110 
111  Array<OneD, NekDouble> rhs(ntotpnt);
112 
113  Array<OneD, NekDouble> outN(ntotpnt);
114  Array<OneD, NekDouble> outTmp(ntotpnt);
115  Array<OneD, Array<OneD, NekDouble>> rhs2d(nvariables);
116  Array<OneD, Array<OneD, NekDouble>> out_2d(nvariables);
117  Array<OneD, Array<OneD, NekDouble>> outTmp_2d(nvariables);
118  for (int m = 0; m < nvariables; m++)
119  {
120  int moffset = m * npoints;
121  rhs2d[m] = rhs + moffset;
122  out_2d[m] = outarray + moffset;
123  outTmp_2d[m] = outTmp + moffset;
124  pFields[m]->MultiplyByMassMatrix(inarray + moffset, rhs2d[m]);
125  }
126 
127  int nphysic = pFields[0]->GetNpoints();
128  int nTracePts = pFields[0]->GetTrace()->GetNpoints();
130  for (int i = 0; i < m_spacedim; i++)
131  {
132  qfield[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
133  for (int j = 0; j < nvariables; j++)
134  {
135  qfield[i][j] = Array<OneD, NekDouble>(nphysic);
136  }
137  }
138  int ntmpTrace = 4 + 2 * m_spacedim;
139  TensorOfArray3D<NekDouble> tmpTrace(ntmpTrace);
140  for (int i = 0; i < ntmpTrace; i++)
141  {
142  tmpTrace[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
143  for (int j = 0; j < nvariables; j++)
144  {
145  tmpTrace[i][j] = Array<OneD, NekDouble>(nTracePts);
146  }
147  }
148  Array<OneD, Array<OneD, NekDouble>> FwdFluxDeriv(nvariables);
149  Array<OneD, Array<OneD, NekDouble>> BwdFluxDeriv(nvariables);
150  for (int j = 0; j < nvariables; j++)
151  {
152  FwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
153  BwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
154  }
155 
156  bool flagUpdateDervFlux = false;
157 
158  const int nwspTraceDataType = nvariables + 1;
159  NekSingle tmpSingle;
160  Array<OneD, Array<OneD, NekSingle>> wspTraceDataType(nwspTraceDataType);
161  for (int m = 0; m < nwspTraceDataType; m++)
162  {
163  wspTraceDataType[m] = Array<OneD, NekSingle>(nTracePts);
164  }
165 
166  PreconBlkDiag(pFields, rhs, outarray, m_PreconMatSingle, tmpSingle);
167 
168  for (int nrelax = 0; nrelax < nBRJIterTot - 1; nrelax++)
169  {
170  Vmath::Smul(ntotpnt, OmBRJParam, outarray, 1, outN, 1);
171 
173  pFields, nvariables, npoints, rhs2d, out_2d, flagUpdateDervFlux,
174  FwdFluxDeriv, BwdFluxDeriv, qfield, tmpTrace, wspTraceDataType,
177 
178  PreconBlkDiag(pFields, outarray, outTmp, m_PreconMatSingle,
179  tmpSingle);
180  Vmath::Svtvp(ntotpnt, BRJParam, outTmp, 1, outN, 1, outarray, 1);
181  }
182  }
183 }
184 
187  const Array<OneD, const Array<OneD, NekDouble>> &intmp,
188  const NekDouble time, const NekDouble lambda)
189 {
190  if (0 < m_PreconItsStep)
191  {
197 
198  if (m_verbose && m_root)
199  {
200  cout << " ## CalcuPreconMat " << endl;
201  }
202  }
203 
204  m_BndEvaluateTime = time;
205  m_DtLambdaPreconMat = lambda;
206 
207  m_CalcPreconMatFlag = false;
209 }
210 
212  const NekDouble dtLambda)
213 {
214  boost::ignore_unused(res);
215 
216  bool flag = false;
217 
218  if (m_CalcPreconMatFlag || (m_DtLambdaPreconMat != dtLambda))
219  {
220  flag = true;
221  }
222 
224  {
225  flag = true;
226  }
227 
228  m_CalcPreconMatFlag = flag;
229  return flag;
230 }
231 
232 template <typename DataType, typename TypeNekBlkMatSharedPtr>
235  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
236  const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
237 {
238  boost::ignore_unused(tmpDataType);
239 
240  unsigned int nvariables = pFields.size();
241  unsigned int npoints = pFields[0]->GetNcoeffs();
242  unsigned int npointsVar = nvariables * npoints;
243  Array<OneD, DataType> Sinarray(npointsVar);
244  Array<OneD, DataType> Soutarray(npointsVar);
245  NekVector<DataType> tmpVect(npointsVar, Sinarray, eWrapper);
246  NekVector<DataType> outVect(npointsVar, Soutarray, eWrapper);
247 
248  std::shared_ptr<LocalRegions::ExpansionVector> expvect =
249  pFields[0]->GetExp();
250  int nTotElmt = (*expvect).size();
251 
252  for (int m = 0; m < nvariables; m++)
253  {
254  int nVarOffset = m * npoints;
255  for (int ne = 0; ne < nTotElmt; ne++)
256  {
257  int nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
258  int nElmtCoef = pFields[0]->GetNcoeffs(ne);
259  int inOffset = nVarOffset + nCoefOffset;
260  int outOffset = nCoefOffset * nvariables + m * nElmtCoef;
261  for (int i = 0; i < nElmtCoef; i++)
262  {
263  Sinarray[outOffset + i] = DataType(inarray[inOffset + i]);
264  }
265  }
266  }
267 
268  outVect = (*PreconMatVars) * tmpVect;
269 
270  for (int m = 0; m < nvariables; m++)
271  {
272  int nVarOffset = m * npoints;
273  for (int ne = 0; ne < nTotElmt; ne++)
274  {
275  int nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
276  int nElmtCoef = pFields[0]->GetNcoeffs(ne);
277  int inOffset = nVarOffset + nCoefOffset;
278  int outOffset = nCoefOffset * nvariables + m * nElmtCoef;
279  for (int i = 0; i < nElmtCoef; i++)
280  {
281  outarray[inOffset + i] = NekDouble(Soutarray[outOffset + i]);
282  }
283  }
284  }
285 }
286 
287 template <typename DataType>
290  const int nvariables, const int nCoeffs,
291  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
292  Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
293  Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
294  Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
296  Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
297  const TensorOfArray4D<DataType> &TraceJacArray,
298  const TensorOfArray4D<DataType> &TraceJacDerivArray,
299  const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
300  const TensorOfArray5D<DataType> &TraceIPSymJacArray)
301 {
302  boost::ignore_unused(flagUpdateDervFlux, qfield, TraceJacDerivArray,
303  TraceJacDerivSign, FwdFluxDeriv, BwdFluxDeriv,
304  TraceIPSymJacArray);
305 
306  int nTracePts = pFields[0]->GetTrace()->GetNpoints();
307  int npoints = pFields[0]->GetNpoints();
308  int nDim = m_spacedim;
309 
310  Array<OneD, Array<OneD, NekDouble>> outpnts(nvariables);
311  for (int i = 0; i < nvariables; i++)
312  {
313  outpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
314  pFields[i]->BwdTrans(outarray[i], outpnts[i]);
315  }
316 
317  // Store forwards/backwards space along trace space
322  TensorOfArray3D<NekDouble> numDerivBwd(nDim);
323  TensorOfArray3D<NekDouble> numDerivFwd(nDim);
324  int indexwspTrace = 0;
325  Fwd = wspTrace[indexwspTrace], indexwspTrace++;
326  Bwd = wspTrace[indexwspTrace], indexwspTrace++;
327  FwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
328  BwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
329 
330  for (int i = 0; i < nvariables; ++i)
331  {
332  pFields[i]->GetFwdBwdTracePhys(outpnts[i], Fwd[i], Bwd[i]);
333  }
334 
335  int indexwspTraceDataType = 0;
336  Array<OneD, Array<OneD, DataType>> Fwdarray(nvariables);
337  for (int m = 0; m < nvariables; ++m)
338  {
339  Fwdarray[m] = wspTraceDataType[indexwspTraceDataType],
340  indexwspTraceDataType++;
341  }
342  Array<OneD, DataType> Fwdreslt;
343  Fwdreslt = wspTraceDataType[indexwspTraceDataType], indexwspTraceDataType++;
344 
345  for (int m = 0; m < nvariables; ++m)
346  {
347  for (int i = 0; i < nTracePts; ++i)
348  {
349  Fwdarray[m][i] = DataType(Fwd[m][i]);
350  }
351  }
352  for (int m = 0; m < nvariables; ++m)
353  {
354  Vmath::Zero(nTracePts, Fwdreslt, 1);
355  for (int n = 0; n < nvariables; ++n)
356  {
357  Vmath::Vvtvp(nTracePts, TraceJacArray[0][m][n], 1, Fwdarray[n], 1,
358  Fwdreslt, 1, Fwdreslt, 1);
359  }
360 
361  for (int i = 0; i < nTracePts; ++i)
362  {
363  FwdFlux[m][i] = NekDouble(Fwdreslt[i]);
364  }
365  }
366 
367  for (int m = 0; m < nvariables; ++m)
368  {
369  for (int i = 0; i < nTracePts; ++i)
370  {
371  Fwdarray[m][i] = DataType(Bwd[m][i]);
372  }
373  }
374  for (int m = 0; m < nvariables; ++m)
375  {
376  Vmath::Zero(nTracePts, Fwdreslt, 1);
377  for (int n = 0; n < nvariables; ++n)
378  {
379  Vmath::Vvtvp(nTracePts, TraceJacArray[1][m][n], 1, Fwdarray[n], 1,
380  Fwdreslt, 1, Fwdreslt, 1);
381  }
382  for (int i = 0; i < nTracePts; ++i)
383  {
384  BwdFlux[m][i] = NekDouble(Fwdreslt[i]);
385  }
386  }
387 
388  for (int i = 0; i < nvariables; ++i)
389  {
390  Vmath::Fill(nCoeffs, 0.0, outarray[i], 1);
391  pFields[i]->AddTraceIntegralToOffDiag(FwdFlux[i], BwdFlux[i],
392  outarray[i]);
393  }
394 
395  for (int i = 0; i < nvariables; ++i)
396  {
397  Vmath::Svtvp(nCoeffs, -m_DtLambdaPreconMat, outarray[i], 1, inarray[i],
398  1, outarray[i], 1);
399  }
400 }
401 
402 template <typename TypeNekBlkMatSharedPtr>
406  const int &nscale)
407 {
408 
409  int nvars = pFields.size();
410  int nelmts = pFields[0]->GetNumElmts();
411  int nelmtcoef;
412  Array<OneD, unsigned int> nelmtmatdim(nelmts);
413  for (int i = 0; i < nelmts; i++)
414  {
415  nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
416  nelmtmatdim[i] = nelmtcoef * nscale;
417  }
418 
419  for (int i = 0; i < nvars; i++)
420  {
421  for (int j = 0; j < nvars; j++)
422  {
423  AllocateNekBlkMatDig(gmtxarray[i][j], nelmtmatdim, nelmtmatdim);
424  }
425  }
426 }
427 
428 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
void DoCalcPreconMatBRJCoeff(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray, SNekBlkMatSharedPtr &gmtVar, Array< OneD, SNekBlkMatSharedPtr > &TraceJac, Array< OneD, SNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, NekSingle >> &TraceJacDerivSign, TensorOfArray4D< NekSingle > &TraceJacArray, TensorOfArray4D< NekSingle > &TraceJacDerivArray, TensorOfArray5D< NekSingle > &TraceIPSymJacArray)
Definition: PreconCfsOp.h:111
virtual void v_InitObject()
virtual void v_DoPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:79
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:80
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:81
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:83
void MinusOffDiag2Rhs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int nvariables, const int nCoeffs, const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, bool flagUpdateDervFlux, Array< OneD, Array< OneD, NekDouble >> &FwdFluxDeriv, Array< OneD, Array< OneD, NekDouble >> &BwdFluxDeriv, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &wspTrace, Array< OneD, Array< OneD, DataType >> &wspTraceDataType, const TensorOfArray4D< DataType > &TraceJacArray, const TensorOfArray4D< DataType > &TraceJacDerivArray, const Array< OneD, const Array< OneD, DataType >> &TraceJacDerivSign, const TensorOfArray5D< DataType > &TraceIPSymJacArray)
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:86
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:82
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const int &nscale=1)
virtual void v_BuildPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble >> &intmp, const NekDouble time, const NekDouble lambda)
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:84
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:130
virtual bool UpdatePreconMatCheck(const Array< OneD, const NekDouble > &res, const NekDouble dtLambda)
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
SNekBlkMatSharedPtr m_PreconMatSingle
Definition: PreconCfsBRJ.h:78
int m_PreconTimesCounter
Definition: PreconCfs.h:85
int m_PreconMatFreezNumb
Definition: PreconCfs.h:84
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
Definition: PreconCfs.cpp:63
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:87
bool m_CalcPreconMatFlag
Definition: PreconCfs.h:90
NekDouble m_BndEvaluateTime
Definition: PreconCfs.h:88
virtual void v_InitObject()
Definition: PreconCfsOp.cpp:56
NekPreconCfsOperators m_operator
Definition: PreconCfsOp.h:196
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
@ eDiagonal
Definition: PreconCfs.h:48
double NekDouble
PreconCfsOpFactory & GetPreconCfsOpFactory()
Declaration of the boundary condition factory singleton.
Definition: PreconCfsOp.cpp:42
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:565
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45