Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::PreconCfsBRJ Class Reference

#include <PreconCfsBRJ.h>

Inheritance diagram for Nektar::PreconCfsBRJ:
[legend]

Public Member Functions

 PreconCfsBRJ (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
 
 ~PreconCfsBRJ ()
 
virtual bool UpdatePreconMatCheck (const Array< OneD, const NekDouble > &res, const NekDouble dtLambda)
 
- Public Member Functions inherited from Nektar::PreconCfsOp
 PreconCfsOp (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
 
virtual ~PreconCfsOp ()
 
void SetOperators (const NekPreconCfsOperators &in)
 
- Public Member Functions inherited from Nektar::PreconCfs
 PreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
 
virtual ~PreconCfs ()
 
void DoPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
 
void BuildPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble >> &intmp, const NekDouble time, const NekDouble lambda)
 
void InitObject ()
 

Static Public Member Functions

static PreconCfsOpSharedPtr create (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

virtual void v_InitObject ()
 
- Protected Member Functions inherited from Nektar::PreconCfs
void DoNullPrecon (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
 

Protected Attributes

int m_PreconItsStep
 
int m_BRJRelaxParam
 
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
 
SNekBlkMatSharedPtr m_PreconMatSingle
 
Array< OneD, SNekBlkMatSharedPtrm_TraceJacSingle
 
TensorOfArray4D< NekSinglem_TraceJacArraySingle
 
Array< OneD, SNekBlkMatSharedPtrm_TraceJacDerivSingle
 
TensorOfArray4D< NekSinglem_TraceJacDerivArraySingle
 
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
 
TensorOfArray5D< NekSinglem_TraceIPSymJacArraySingle
 
PrecType m_PreconMatStorage
 
- Protected Attributes inherited from Nektar::PreconCfsOp
NekPreconCfsOperators m_operator
 
- Protected Attributes inherited from Nektar::PreconCfs
LibUtilities::CommSharedPtr m_Comm
 
bool m_root
 
bool m_verbose
 
int m_spacedim
 
int m_PreconMatFreezNumb
 
int m_PreconTimesCounter
 
NekDouble m_DtLambdaPreconMat = -1.0
 
NekDouble m_BndEvaluateTime
 
bool m_CalcPreconMatFlag = false
 

Private Member Functions

virtual void v_DoPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
 
virtual void v_BuildPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble >> &intmp, const NekDouble time, const NekDouble lambda)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void PreconBlkDiag (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
 
template<typename DataType >
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)
 
template<typename TypeNekBlkMatSharedPtr >
void AllocatePreconBlkDiagCoeff (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const int &nscale=1)
 
void AllocateNekBlkMatDig (SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
 

Friends

class MemoryManager< PreconCfsBRJ >
 

Detailed Description

Block Relaxed(weighted) Jacobi iterative (BRJ) Preconditioner for CFS

Solves a linear system using iterative methods.

Definition at line 46 of file PreconCfsBRJ.h.

Constructor & Destructor Documentation

◆ PreconCfsBRJ()

Nektar::PreconCfsBRJ::PreconCfsBRJ ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vComm 
)

Definition at line 52 of file PreconCfsBRJ.cpp.

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 
64  m_PreconMatVarsSingle = TensorOfArray2D<SNekBlkMatSharedPtr>(nvariables);
65  for (int i = 0; i < nvariables; i++)
66  {
67  m_PreconMatVarsSingle[i] = Array<OneD, SNekBlkMatSharedPtr>(nvariables);
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 }
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:86
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const int &nscale=1)
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:130
SNekBlkMatSharedPtr m_PreconMatSingle
Definition: PreconCfsBRJ.h:78
PreconCfsOp(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Definition: PreconCfsOp.cpp:48
@ eDiagonal
Definition: PreconCfs.h:48

References AllocateNekBlkMatDig(), AllocatePreconBlkDiagCoeff(), Nektar::eDiagonal, m_BRJRelaxParam, m_PreconItsStep, m_PreconMatSingle, m_PreconMatStorage, and m_PreconMatVarsSingle.

◆ ~PreconCfsBRJ()

Nektar::PreconCfsBRJ::~PreconCfsBRJ ( )
inline

Definition at line 68 of file PreconCfsBRJ.h.

68 {};

Member Function Documentation

◆ AllocateNekBlkMatDig()

void Nektar::PreconCfsBRJ::AllocateNekBlkMatDig ( SNekBlkMatSharedPtr mat,
const Array< OneD, unsigned int >  nrow,
const Array< OneD, unsigned int >  ncol 
)
inlineprivate

Definition at line 130 of file PreconCfsBRJ.h.

133  {
134  mat =
136  SNekMatSharedPtr loc_matNvar;
137  for (int nelm = 0; nelm < nrow.size(); ++nelm)
138  {
139  int nrowsVars = nrow[nelm];
140  int ncolsVars = ncol[nelm];
141 
143  nrowsVars, ncolsVars, 0.0);
144  mat->SetBlock(nelm, nelm, loc_matNvar);
145  }
146  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< SNekMat > SNekMatSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::eDIAGONAL.

Referenced by AllocatePreconBlkDiagCoeff(), and PreconCfsBRJ().

◆ AllocatePreconBlkDiagCoeff()

template<typename TypeNekBlkMatSharedPtr >
void Nektar::PreconCfsBRJ::AllocatePreconBlkDiagCoeff ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &  gmtxarray,
const int &  nscale = 1 
)
private

Definition at line 403 of file PreconCfsBRJ.cpp.

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 }

References AllocateNekBlkMatDig().

Referenced by PreconCfsBRJ().

◆ create()

static PreconCfsOpSharedPtr Nektar::PreconCfsBRJ::create ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const LibUtilities::SessionReaderSharedPtr pSession,
const LibUtilities::CommSharedPtr vComm 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file PreconCfsBRJ.h.

56  {
58  pFields, pSession, vComm);
59  return p;
60  }
std::shared_ptr< PreconCfsOp > PreconCfsOpSharedPtr
Definition: PreconCfsOp.h:152

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ MinusOffDiag2Rhs()

template<typename DataType >
void Nektar::PreconCfsBRJ::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 
)
private

Definition at line 288 of file PreconCfsBRJ.cpp.

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
318  Array<OneD, Array<OneD, NekDouble>> Fwd;
319  Array<OneD, Array<OneD, NekDouble>> Bwd;
320  Array<OneD, Array<OneD, NekDouble>> FwdFlux;
321  Array<OneD, Array<OneD, NekDouble>> BwdFlux;
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 }
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:87
double NekDouble
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 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

References Vmath::Fill(), Nektar::PreconCfs::m_DtLambdaPreconMat, Nektar::PreconCfs::m_spacedim, Vmath::Svtvp(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_DoPreconCfs().

◆ PreconBlkDiag()

template<typename DataType , typename TypeNekBlkMatSharedPtr >
void Nektar::PreconCfsBRJ::PreconBlkDiag ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const TypeNekBlkMatSharedPtr &  PreconMatVars,
const DataType &  tmpDataType 
)
private

Definition at line 233 of file PreconCfsBRJ.cpp.

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 }

References Nektar::eWrapper.

Referenced by v_DoPreconCfs().

◆ UpdatePreconMatCheck()

bool Nektar::PreconCfsBRJ::UpdatePreconMatCheck ( const Array< OneD, const NekDouble > &  res,
const NekDouble  dtLambda 
)
virtual

Reimplemented from Nektar::PreconCfs.

Definition at line 211 of file PreconCfsBRJ.cpp.

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 }
int m_PreconTimesCounter
Definition: PreconCfs.h:85
int m_PreconMatFreezNumb
Definition: PreconCfs.h:84
bool m_CalcPreconMatFlag
Definition: PreconCfs.h:90

References Nektar::PreconCfs::m_CalcPreconMatFlag, Nektar::PreconCfs::m_DtLambdaPreconMat, Nektar::PreconCfs::m_PreconMatFreezNumb, and Nektar::PreconCfs::m_PreconTimesCounter.

◆ v_BuildPreconCfs()

void Nektar::PreconCfsBRJ::v_BuildPreconCfs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, const Array< OneD, NekDouble >> &  intmp,
const NekDouble  time,
const NekDouble  lambda 
)
privatevirtual

Reimplemented from Nektar::PreconCfsOp.

Definition at line 185 of file PreconCfsBRJ.cpp.

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 }
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
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:79
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:80
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:81
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:83
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:82
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:84
NekDouble m_BndEvaluateTime
Definition: PreconCfs.h:88
NekPreconCfsOperators m_operator
Definition: PreconCfsOp.h:196

References Nektar::NekPreconCfsOperators::DoCalcPreconMatBRJCoeff(), Nektar::PreconCfs::m_BndEvaluateTime, Nektar::PreconCfs::m_CalcPreconMatFlag, Nektar::PreconCfs::m_DtLambdaPreconMat, Nektar::PreconCfsOp::m_operator, m_PreconItsStep, m_PreconMatSingle, m_PreconMatVarsSingle, Nektar::PreconCfs::m_PreconTimesCounter, Nektar::PreconCfs::m_root, m_TraceIPSymJacArraySingle, m_TraceJacArraySingle, m_TraceJacDerivArraySingle, m_TraceJacDerivSignSingle, m_TraceJacDerivSingle, m_TraceJacSingle, and Nektar::PreconCfs::m_verbose.

◆ v_DoPreconCfs()

void Nektar::PreconCfsBRJ::v_DoPreconCfs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  flag 
)
privatevirtual

Reimplemented from Nektar::PreconCfsOp.

Definition at line 87 of file PreconCfsBRJ.cpp.

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();
129  TensorOfArray3D<NekDouble> qfield(m_spacedim);
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 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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)
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const TypeNekBlkMatSharedPtr &PreconMatVars, const DataType &tmpDataType)
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
Definition: PreconCfs.cpp:63
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

References ASSERTL0, Nektar::PreconCfs::DoNullPrecon(), m_BRJRelaxParam, m_PreconItsStep, m_PreconMatSingle, Nektar::PreconCfs::m_spacedim, m_TraceIPSymJacArraySingle, m_TraceJacArraySingle, m_TraceJacDerivArraySingle, m_TraceJacDerivSignSingle, MinusOffDiag2Rhs(), PreconBlkDiag(), Vmath::Smul(), and Vmath::Svtvp().

◆ v_InitObject()

void Nektar::PreconCfsBRJ::v_InitObject ( )
protectedvirtual

Reimplemented from Nektar::PreconCfsOp.

Definition at line 82 of file PreconCfsBRJ.cpp.

83 {
85 }
virtual void v_InitObject()
Definition: PreconCfsOp.cpp:56

References Nektar::PreconCfsOp::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< PreconCfsBRJ >

friend class MemoryManager< PreconCfsBRJ >
friend

Definition at line 1 of file PreconCfsBRJ.h.

Member Data Documentation

◆ className

std::string Nektar::PreconCfsBRJ::className
static
Initial value:
=
"PreconCfsBRJ", PreconCfsBRJ::create,
"Block Relaxed Jacobi Preconditioner for CFS.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static PreconCfsOpSharedPtr create(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Creates an instance of this class.
Definition: PreconCfsBRJ.h:52
PreconCfsOpFactory & GetPreconCfsOpFactory()
Declaration of the boundary condition factory singleton.
Definition: PreconCfsOp.cpp:42

Name of the class.

Definition at line 63 of file PreconCfsBRJ.h.

◆ m_BRJRelaxParam

int Nektar::PreconCfsBRJ::m_BRJRelaxParam
protected

Definition at line 75 of file PreconCfsBRJ.h.

Referenced by PreconCfsBRJ(), and v_DoPreconCfs().

◆ m_PreconItsStep

int Nektar::PreconCfsBRJ::m_PreconItsStep
protected

Definition at line 74 of file PreconCfsBRJ.h.

Referenced by PreconCfsBRJ(), v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_PreconMatSingle

SNekBlkMatSharedPtr Nektar::PreconCfsBRJ::m_PreconMatSingle
protected

Definition at line 78 of file PreconCfsBRJ.h.

Referenced by PreconCfsBRJ(), v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_PreconMatStorage

PrecType Nektar::PreconCfsBRJ::m_PreconMatStorage
protected

Definition at line 86 of file PreconCfsBRJ.h.

Referenced by PreconCfsBRJ().

◆ m_PreconMatVarsSingle

Array<OneD, Array<OneD, SNekBlkMatSharedPtr> > Nektar::PreconCfsBRJ::m_PreconMatVarsSingle
protected

Definition at line 77 of file PreconCfsBRJ.h.

Referenced by PreconCfsBRJ(), and v_BuildPreconCfs().

◆ m_TraceIPSymJacArraySingle

TensorOfArray5D<NekSingle> Nektar::PreconCfsBRJ::m_TraceIPSymJacArraySingle
protected

Definition at line 84 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacArraySingle

TensorOfArray4D<NekSingle> Nektar::PreconCfsBRJ::m_TraceJacArraySingle
protected

Definition at line 80 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacDerivArraySingle

TensorOfArray4D<NekSingle> Nektar::PreconCfsBRJ::m_TraceJacDerivArraySingle
protected

Definition at line 82 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacDerivSignSingle

Array<OneD, Array<OneD, NekSingle> > Nektar::PreconCfsBRJ::m_TraceJacDerivSignSingle
protected

Definition at line 83 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacDerivSingle

Array<OneD, SNekBlkMatSharedPtr> Nektar::PreconCfsBRJ::m_TraceJacDerivSingle
protected

Definition at line 81 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs().

◆ m_TraceJacSingle

Array<OneD, SNekBlkMatSharedPtr> Nektar::PreconCfsBRJ::m_TraceJacSingle
protected

Definition at line 79 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs().