Nektar++
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 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 /**
44  * @class PreconCfsBRJ
45  *
46  * Solves a linear system using iterative methods.
47  */
48 std::string PreconCfsBRJ::className =
50  "PreconCfsBRJ", PreconCfsBRJ::create,
51  "Block Relaxed Jacobi Preconditioner for CFS.");
52 
53 PreconCfsBRJ::PreconCfsBRJ(
56  const LibUtilities::CommSharedPtr &vComm)
57  : PreconCfsOp(pFields, pSession, vComm)
58 {
59  pSession->LoadParameter("PreconItsStep", m_PreconItsStep, 7);
60  pSession->LoadParameter("BRJRelaxParam", m_BRJRelaxParam, 1.0);
61 
62  int nvariables = pFields.size();
64 
66  for (int i = 0; i < nvariables; i++)
67  {
69  }
71 
72  int nelmts = pFields[0]->GetNumElmts();
73  int nelmtcoef;
74  Array<OneD, unsigned int> nelmtmatdim(nelmts);
75  for (int i = 0; i < nelmts; i++)
76  {
77  nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
78  nelmtmatdim[i] = nelmtcoef * nvariables;
79  }
80  AllocateNekBlkMatDig(m_PreconMatSingle, nelmtmatdim, nelmtmatdim);
81 }
82 
84 {
86 }
87 
90  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
91  const bool &flag)
92 {
93  boost::ignore_unused(flag);
94 
95  int nBRJIterTot = m_PreconItsStep;
96  if (0 == nBRJIterTot)
97  {
98  DoNullPrecon(inarray, outarray, flag);
99  }
100  else
101  {
102  const NekDouble BRJParam = m_BRJRelaxParam;
103  const NekDouble OmBRJParam = 1.0 - BRJParam;
104 
105  unsigned int nvariables = pFields.size();
106  unsigned int npoints = pFields[0]->GetNcoeffs();
107  unsigned int ntotpnt = inarray.size();
108 
109  ASSERTL0(nvariables * npoints == ntotpnt,
110  "nvariables*npoints!=ntotpnt in PreconCoeff");
111 
112  Array<OneD, NekDouble> rhs(ntotpnt);
113 
114  Array<OneD, NekDouble> outN(ntotpnt);
115  Array<OneD, NekDouble> outTmp(ntotpnt);
116  Array<OneD, Array<OneD, NekDouble>> rhs2d(nvariables);
117  Array<OneD, Array<OneD, NekDouble>> out_2d(nvariables);
118  Array<OneD, Array<OneD, NekDouble>> outTmp_2d(nvariables);
119  for (int m = 0; m < nvariables; m++)
120  {
121  int moffset = m * npoints;
122  rhs2d[m] = rhs + moffset;
123  out_2d[m] = outarray + moffset;
124  outTmp_2d[m] = outTmp + moffset;
125  pFields[m]->MultiplyByMassMatrix(inarray + moffset, rhs2d[m]);
126  }
127 
128  int nphysic = pFields[0]->GetNpoints();
129  int nTracePts = pFields[0]->GetTrace()->GetNpoints();
131  for (int i = 0; i < m_spacedim; i++)
132  {
133  qfield[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
134  for (int j = 0; j < nvariables; j++)
135  {
136  qfield[i][j] = Array<OneD, NekDouble>(nphysic);
137  }
138  }
139  int ntmpTrace = 4 + 2 * m_spacedim;
140  TensorOfArray3D<NekDouble> tmpTrace(ntmpTrace);
141  for (int i = 0; i < ntmpTrace; i++)
142  {
143  tmpTrace[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
144  for (int j = 0; j < nvariables; j++)
145  {
146  tmpTrace[i][j] = Array<OneD, NekDouble>(nTracePts);
147  }
148  }
149  Array<OneD, Array<OneD, NekDouble>> FwdFluxDeriv(nvariables);
150  Array<OneD, Array<OneD, NekDouble>> BwdFluxDeriv(nvariables);
151  for (int j = 0; j < nvariables; j++)
152  {
153  FwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
154  BwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
155  }
156 
157  bool flagUpdateDervFlux = false;
158 
159  const int nwspTraceDataType = nvariables + 1;
160  NekSingle tmpSingle;
161  Array<OneD, Array<OneD, NekSingle>> wspTraceDataType(nwspTraceDataType);
162  for (int m = 0; m < nwspTraceDataType; m++)
163  {
164  wspTraceDataType[m] = Array<OneD, NekSingle>(nTracePts);
165  }
166 
167  LibUtilities::Timer timer;
168  timer.Start();
169  PreconBlkDiag(pFields, rhs, outarray, m_PreconMatSingle, tmpSingle);
170  timer.Stop();
171  timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
172 
173  for (int nrelax = 0; nrelax < nBRJIterTot - 1; nrelax++)
174  {
175  Vmath::Smul(ntotpnt, OmBRJParam, outarray, 1, outN, 1);
176 
177  timer.Start();
179  pFields, nvariables, npoints, rhs2d, out_2d, flagUpdateDervFlux,
180  FwdFluxDeriv, BwdFluxDeriv, qfield, tmpTrace, wspTraceDataType,
183  timer.Stop();
184  timer.AccumulateRegion("PreconCfsBRJ::MinusOffDiag2Rhs", 2);
185 
186  timer.Start();
187  PreconBlkDiag(pFields, outarray, outTmp, m_PreconMatSingle,
188  tmpSingle);
189  timer.Stop();
190  timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
191 
192  Vmath::Svtvp(ntotpnt, BRJParam, outTmp, 1, outN, 1, outarray, 1);
193  }
194  }
195 }
196 
199  const Array<OneD, const Array<OneD, NekDouble>> &intmp,
200  const NekDouble time, const NekDouble lambda)
201 {
202  if (0 < m_PreconItsStep)
203  {
209 
210  if (m_verbose && m_root)
211  {
212  cout << " ## CalcuPreconMat " << endl;
213  }
214  }
215 
216  m_BndEvaluateTime = time;
217  m_DtLambdaPreconMat = lambda;
218 
219  m_CalcPreconMatFlag = false;
221 }
222 
224  const NekDouble dtLambda)
225 {
226  boost::ignore_unused(res);
227 
228  bool flag = false;
229 
230  if (m_CalcPreconMatFlag || (m_DtLambdaPreconMat != dtLambda))
231  {
232  flag = true;
233  }
234 
236  {
237  flag = true;
238  }
239 
240  m_CalcPreconMatFlag = flag;
241  return flag;
242 }
243 
244 template <typename DataType, typename TypeNekBlkMatSharedPtr>
247  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
248  const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
249 {
250  boost::ignore_unused(tmpDataType);
251 
252  unsigned int nvariables = pFields.size();
253  unsigned int npoints = pFields[0]->GetNcoeffs();
254  unsigned int npointsVar = nvariables * npoints;
255  Array<OneD, DataType> Sinarray(npointsVar);
256  Array<OneD, DataType> Soutarray(npointsVar);
257  NekVector<DataType> tmpVect(npointsVar, Sinarray, eWrapper);
258  NekVector<DataType> outVect(npointsVar, Soutarray, eWrapper);
259 
260  std::shared_ptr<LocalRegions::ExpansionVector> expvect =
261  pFields[0]->GetExp();
262  int nTotElmt = (*expvect).size();
263 
264  for (int m = 0; m < nvariables; m++)
265  {
266  int nVarOffset = m * npoints;
267  for (int ne = 0; ne < nTotElmt; ne++)
268  {
269  int nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
270  int nElmtCoef = pFields[0]->GetNcoeffs(ne);
271  int inOffset = nVarOffset + nCoefOffset;
272  int outOffset = nCoefOffset * nvariables + m * nElmtCoef;
273  for (int i = 0; i < nElmtCoef; i++)
274  {
275  Sinarray[outOffset + i] = DataType(inarray[inOffset + i]);
276  }
277  }
278  }
279 
280  outVect = (*PreconMatVars) * tmpVect;
281 
282  for (int m = 0; m < nvariables; m++)
283  {
284  int nVarOffset = m * npoints;
285  for (int ne = 0; ne < nTotElmt; ne++)
286  {
287  int nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
288  int nElmtCoef = pFields[0]->GetNcoeffs(ne);
289  int inOffset = nVarOffset + nCoefOffset;
290  int outOffset = nCoefOffset * nvariables + m * nElmtCoef;
291  for (int i = 0; i < nElmtCoef; i++)
292  {
293  outarray[inOffset + i] = NekDouble(Soutarray[outOffset + i]);
294  }
295  }
296  }
297 }
298 
299 template <typename DataType>
302  const int nvariables, const int nCoeffs,
303  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
304  Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
305  Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
306  Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
308  Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
309  const TensorOfArray4D<DataType> &TraceJacArray,
310  const TensorOfArray4D<DataType> &TraceJacDerivArray,
311  const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
312  const TensorOfArray5D<DataType> &TraceIPSymJacArray)
313 {
314  boost::ignore_unused(flagUpdateDervFlux, qfield, TraceJacDerivArray,
315  TraceJacDerivSign, FwdFluxDeriv, BwdFluxDeriv,
316  TraceIPSymJacArray);
317 
318  int nTracePts = pFields[0]->GetTrace()->GetNpoints();
319  int npoints = pFields[0]->GetNpoints();
320  int nDim = m_spacedim;
321 
322  Array<OneD, Array<OneD, NekDouble>> outpnts(nvariables);
323  for (int i = 0; i < nvariables; i++)
324  {
325  outpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
326  pFields[i]->BwdTrans(outarray[i], outpnts[i]);
327  }
328 
329  // Store forwards/backwards space along trace space
334  TensorOfArray3D<NekDouble> numDerivBwd(nDim);
335  TensorOfArray3D<NekDouble> numDerivFwd(nDim);
336  int indexwspTrace = 0;
337  Fwd = wspTrace[indexwspTrace], indexwspTrace++;
338  Bwd = wspTrace[indexwspTrace], indexwspTrace++;
339  FwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
340  BwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
341 
342  LibUtilities::Timer timer;
343  for (int i = 0; i < nvariables; ++i)
344  {
345  timer.Start();
346  pFields[i]->GetFwdBwdTracePhys(outpnts[i], Fwd[i], Bwd[i]);
347  timer.Stop();
348  timer.AccumulateRegion("ExpList::GetFwdBwdTracePhys", 10);
349  }
350 
351  int indexwspTraceDataType = 0;
352  Array<OneD, Array<OneD, DataType>> Fwdarray(nvariables);
353  for (int m = 0; m < nvariables; ++m)
354  {
355  Fwdarray[m] = wspTraceDataType[indexwspTraceDataType],
356  indexwspTraceDataType++;
357  }
358  Array<OneD, DataType> Fwdreslt;
359  Fwdreslt = wspTraceDataType[indexwspTraceDataType], indexwspTraceDataType++;
360 
361  for (int m = 0; m < nvariables; ++m)
362  {
363  for (int i = 0; i < nTracePts; ++i)
364  {
365  Fwdarray[m][i] = DataType(Fwd[m][i]);
366  }
367  }
368  for (int m = 0; m < nvariables; ++m)
369  {
370  Vmath::Zero(nTracePts, Fwdreslt, 1);
371  for (int n = 0; n < nvariables; ++n)
372  {
373  Vmath::Vvtvp(nTracePts, TraceJacArray[0][m][n], 1, Fwdarray[n], 1,
374  Fwdreslt, 1, Fwdreslt, 1);
375  }
376 
377  for (int i = 0; i < nTracePts; ++i)
378  {
379  FwdFlux[m][i] = NekDouble(Fwdreslt[i]);
380  }
381  }
382 
383  for (int m = 0; m < nvariables; ++m)
384  {
385  for (int i = 0; i < nTracePts; ++i)
386  {
387  Fwdarray[m][i] = DataType(Bwd[m][i]);
388  }
389  }
390  for (int m = 0; m < nvariables; ++m)
391  {
392  Vmath::Zero(nTracePts, Fwdreslt, 1);
393  for (int n = 0; n < nvariables; ++n)
394  {
395  Vmath::Vvtvp(nTracePts, TraceJacArray[1][m][n], 1, Fwdarray[n], 1,
396  Fwdreslt, 1, Fwdreslt, 1);
397  }
398  for (int i = 0; i < nTracePts; ++i)
399  {
400  BwdFlux[m][i] = NekDouble(Fwdreslt[i]);
401  }
402  }
403 
404  for (int i = 0; i < nvariables; ++i)
405  {
406  Vmath::Fill(nCoeffs, 0.0, outarray[i], 1);
407  timer.Start();
408  pFields[i]->AddTraceIntegralToOffDiag(FwdFlux[i], BwdFlux[i],
409  outarray[i]);
410  timer.Stop();
411  timer.AccumulateRegion("ExpList::AddTraceIntegralToOffDiag", 10);
412  }
413 
414  for (int i = 0; i < nvariables; ++i)
415  {
416  Vmath::Svtvp(nCoeffs, -m_DtLambdaPreconMat, outarray[i], 1, inarray[i],
417  1, outarray[i], 1);
418  }
419 }
420 
421 template <typename TypeNekBlkMatSharedPtr>
425  const int &nscale)
426 {
427 
428  int nvars = pFields.size();
429  int nelmts = pFields[0]->GetNumElmts();
430  int nelmtcoef;
431  Array<OneD, unsigned int> nelmtmatdim(nelmts);
432  for (int i = 0; i < nelmts; i++)
433  {
434  nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
435  nelmtmatdim[i] = nelmtcoef * nscale;
436  }
437 
438  for (int i = 0; i < nvars; i++)
439  {
440  for (int j = 0; j < nvars; j++)
441  {
442  AllocateNekBlkMatDig(gmtxarray[i][j], nelmtmatdim, nelmtmatdim);
443  }
444  }
445 }
446 
447 } // 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
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:73
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:622
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:574
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:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45