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 ()
 
- 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 InitObject ()
 
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)
 
bool UpdatePreconMatCheck (const Array< OneD, const NekDouble > &res, const NekDouble dtLambda)
 
void SetOperators (const NekPreconCfsOperators &in)
 

Static Public Member Functions

static PreconCfsSharedPtr 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 () override
 
virtual void v_DoPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag) override
 
virtual void v_BuildPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble > > &intmp, const NekDouble time, const NekDouble lambda) override
 
virtual bool v_UpdatePreconMatCheck (const Array< OneD, const NekDouble > &res, const NekDouble dtLambda) override
 
virtual void v_InitObject ()=0
 
virtual void v_DoPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)=0
 
virtual void v_BuildPreconCfs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble > > &intmp, const NekDouble time, const NekDouble lambda)=0
 
virtual bool v_UpdatePreconMatCheck (const Array< OneD, const NekDouble > &res, const NekDouble dtLambda)=0
 

Protected Attributes

int m_PreconItsStep
 
int m_BRJRelaxParam
 
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
 
unsigned int m_max_nblocks
 
unsigned int m_max_nElmtDof
 
std::vector< simd< NekSingle >, tinysimd::allocator< simd< NekSingle > > > m_sBlkDiagMat
 
std::vector< int > m_inputIdx
 
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::PreconCfs
LibUtilities::CommSharedPtr m_Comm
 
bool m_verbose
 
int m_spacedim
 
NekPreconCfsOperators m_operator
 
int m_PreconMatFreezNumb
 
int m_PreconTimesCounter
 
NekDouble m_DtLambdaPreconMat = -1.0
 
NekDouble m_BndEvaluateTime
 
bool m_CalcPreconMatFlag = false
 

Private Member Functions

void DoNullPrecon (const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
 
void PreconBlkDiag (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 
template<typename DataType >
void MinusOffDiag2Rhs (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const size_t nvariables, const size_t 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 AllocateSIMDPreconBlkMatDiag (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 This function creates the matrix structure for the block diagonal operator. It organizes the way that the matrix. More...
 
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 49 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 53 of file PreconCfsBRJ.cpp.

57 : PreconCfs(pFields, pSession, vComm)
58{
59 pSession->LoadParameter("PreconItsStep", m_PreconItsStep, 7);
60 pSession->LoadParameter("BRJRelaxParam", m_BRJRelaxParam, 1.0);
61
62 size_t nvariables = pFields.size();
64
66 Array<OneD, Array<OneD, SNekBlkMatSharedPtr>>(nvariables);
67 for (size_t i = 0; i < nvariables; i++)
68 {
69 m_PreconMatVarsSingle[i] = Array<OneD, SNekBlkMatSharedPtr>(nvariables);
70 }
72
74}
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:92
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const int &nscale=1)
void AllocateSIMDPreconBlkMatDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
This function creates the matrix structure for the block diagonal operator. It organizes the way that...
Definition: PreconCfsBRJ.h:146
PreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Definition: PreconCfs.cpp:48
@ eDiagonal
Definition: PreconCfs.h:49

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

◆ ~PreconCfsBRJ()

Nektar::PreconCfsBRJ::~PreconCfsBRJ ( )
inline

Definition at line 71 of file PreconCfsBRJ.h.

71{};

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 244 of file PreconCfsBRJ.h.

247 {
248 mat =
250 SNekMatSharedPtr loc_matNvar;
251 for (size_t nelm = 0; nelm < nrow.size(); ++nelm)
252 {
253 int nrowsVars = nrow[nelm];
254 int ncolsVars = ncol[nelm];
255
257 nrowsVars, ncolsVars, 0.0);
258 mat->SetBlock(nelm, nelm, loc_matNvar);
259 }
260 }
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 v_BuildPreconCfs().

◆ 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 533 of file PreconCfsBRJ.cpp.

537{
538
539 size_t nvars = pFields.size();
540 size_t nelmts = pFields[0]->GetNumElmts();
541 size_t nelmtcoef;
542 Array<OneD, unsigned int> nelmtmatdim(nelmts);
543 for (size_t i = 0; i < nelmts; i++)
544 {
545 nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
546 nelmtmatdim[i] = nelmtcoef * nscale;
547 }
548
549 for (size_t i = 0; i < nvars; i++)
550 {
551 for (size_t j = 0; j < nvars; j++)
552 {
553 AllocateNekBlkMatDig(gmtxarray[i][j], nelmtmatdim, nelmtmatdim);
554 }
555 }
556}
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:244

References AllocateNekBlkMatDig().

Referenced by PreconCfsBRJ().

◆ AllocateSIMDPreconBlkMatDiag()

void Nektar::PreconCfsBRJ::AllocateSIMDPreconBlkMatDiag ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields)
inlineprivate

This function creates the matrix structure for the block diagonal operator. It organizes the way that the matrix.

by SIMD instructions. The degrees of freedom of each element are reorganized, so they are placed in the correct location to perfom SIMD instructions.

Definition at line 146 of file PreconCfsBRJ.h.

148 {
149 using vec_t = simd<NekSingle>;
150
151 int TotMatLen = 0;
152 int TotLen = 0;
153 const auto nTotElmt = pFields[0]->GetNumElmts();
154 const auto nvariables = pFields.size();
155 const auto vecwidth = vec_t::width;
156
157 m_max_nblocks = 0;
158 m_max_nElmtDof = 0;
159
160 for (int ne = 0; ne < nTotElmt; ne++)
161 {
162 const auto nElmtDof = pFields[0]->GetNcoeffs(ne) * nvariables;
163 const auto nblocks = nElmtDof / vecwidth;
164 unsigned int totblocks =
165 (nElmtDof % vecwidth) ? nblocks + 1 : nblocks;
166
168 (m_max_nblocks > totblocks) ? m_max_nblocks : totblocks;
170 (m_max_nElmtDof > nElmtDof) ? m_max_nElmtDof : nElmtDof;
171
172 TotLen += totblocks * vecwidth;
173 TotMatLen += nElmtDof * totblocks;
174 }
175
176 m_sBlkDiagMat.resize(TotMatLen);
177 m_inputIdx.resize(TotLen);
178
179 // generate a index list of vector width aware mapping from
180 // local coeff storage over all variables to elemental storage
181 // over variables
182 unsigned int ncoeffs = pFields[0]->GetNcoeffs();
183 for (int ne = 0, cnt1 = 0; ne < nTotElmt; ne++)
184 {
185 const auto nElmtCoeff = pFields[0]->GetNcoeffs(ne);
186 const auto nElmtDof = nElmtCoeff * nvariables;
187 const auto nblocks = nElmtDof / vecwidth;
188 const auto nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
189 int i = 0;
190 int i0 = 0;
191 int inOffset, j;
192
193 for (int m = 0; m < nvariables; m++)
194 {
195 inOffset = m * ncoeffs + nCoefOffset;
196
197 if (m && (vecwidth - i0 < nElmtCoeff))
198 {
199 // May need to add entries from later variables to
200 // remainder of last variable if the vector width
201 // was not exact multiple of number of elemental
202 // coeffs
203 for (i = 0; i0 < vecwidth; ++i, ++i0)
204 {
205 m_inputIdx[cnt1++] = inOffset + i;
206 }
207 }
208 else
209 {
210 i = 0;
211 }
212
213 // load up other vectors in variable that fit into vector
214 // width
215 for (j = 0; (j + 1) * vecwidth < nElmtCoeff - i; ++j)
216 {
217 for (i0 = 0; i0 < vecwidth; ++i0)
218 {
219 m_inputIdx[cnt1++] = inOffset + i + j * vecwidth + i0;
220 }
221 }
222
223 // load up any residual data for this variable
224 for (i0 = 0, j = j * vecwidth; j < nElmtCoeff - i; ++j, ++i0)
225 {
226 m_inputIdx[cnt1++] = inOffset + i + j;
227 }
228 }
229
230 const auto endwidth = nElmtDof - nblocks * vecwidth;
231
232 // fill out rest of index to match vector width with last entry
233 if (endwidth)
234 {
235 for (i0 = endwidth; i0 < vecwidth; ++i0)
236 {
237 m_inputIdx[cnt1++] = inOffset + i + j - 1;
238 }
239 }
240 ASSERTL1(cnt1 <= TotLen, "m_inputIdx over extended");
241 }
242 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
unsigned int m_max_nblocks
Definition: PreconCfsBRJ.h:79
unsigned int m_max_nElmtDof
Definition: PreconCfsBRJ.h:80
std::vector< int > m_inputIdx
Definition: PreconCfsBRJ.h:83
std::vector< simd< NekSingle >, tinysimd::allocator< simd< NekSingle > > > m_sBlkDiagMat
Definition: PreconCfsBRJ.h:82
tinysimd::simd< NekDouble > vec_t

References ASSERTL1, m_inputIdx, m_max_nblocks, m_max_nElmtDof, and m_sBlkDiagMat.

Referenced by PreconCfsBRJ().

◆ create()

static PreconCfsSharedPtr 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 55 of file PreconCfsBRJ.h.

59 {
61 pFields, pSession, vComm);
62 return p;
63 }
std::shared_ptr< PreconCfs > PreconCfsSharedPtr
Definition: PreconCfs.h:126

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

◆ DoNullPrecon()

void Nektar::PreconCfsBRJ::DoNullPrecon ( const Array< OneD, NekDouble > &  pInput,
Array< OneD, NekDouble > &  pOutput,
const bool &  flag 
)
private

Definition at line 311 of file PreconCfsBRJ.cpp.

314{
315 boost::ignore_unused(flag);
316 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
317}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

References Vmath::Vcopy().

Referenced by v_DoPreconCfs().

◆ MinusOffDiag2Rhs()

template<typename DataType >
void Nektar::PreconCfsBRJ::MinusOffDiag2Rhs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const size_t  nvariables,
const size_t  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 402 of file PreconCfsBRJ.cpp.

415{
416 boost::ignore_unused(flagUpdateDervFlux, qfield, TraceJacDerivArray,
417 TraceJacDerivSign, FwdFluxDeriv, BwdFluxDeriv,
418 TraceIPSymJacArray);
419
420 size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
421 size_t npoints = pFields[0]->GetNpoints();
422 size_t nDim = m_spacedim;
423
424 Array<OneD, Array<OneD, NekDouble>> outpnts(nvariables);
425 for (size_t i = 0; i < nvariables; i++)
426 {
427 outpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
428 pFields[i]->BwdTrans(outarray[i], outpnts[i]);
429 }
430
431 // Store forwards/backwards space along trace space
432 Array<OneD, Array<OneD, NekDouble>> Fwd;
433 Array<OneD, Array<OneD, NekDouble>> Bwd;
434 Array<OneD, Array<OneD, NekDouble>> FwdFlux;
435 Array<OneD, Array<OneD, NekDouble>> BwdFlux;
436 TensorOfArray3D<NekDouble> numDerivBwd(nDim);
437 TensorOfArray3D<NekDouble> numDerivFwd(nDim);
438 size_t indexwspTrace = 0;
439 Fwd = wspTrace[indexwspTrace], indexwspTrace++;
440 Bwd = wspTrace[indexwspTrace], indexwspTrace++;
441 FwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
442 BwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
443
444 LibUtilities::Timer timer;
445 for (size_t i = 0; i < nvariables; ++i)
446 {
447 timer.Start();
448 pFields[i]->GetFwdBwdTracePhys(outpnts[i], Fwd[i], Bwd[i]);
449 timer.Stop();
450 timer.AccumulateRegion("ExpList::GetFwdBwdTracePhys", 10);
451 }
452
453 size_t indexwspTraceDataType = 0;
454 Array<OneD, Array<OneD, DataType>> Fwdarray(nvariables);
455 for (size_t m = 0; m < nvariables; ++m)
456 {
457 Fwdarray[m] = wspTraceDataType[indexwspTraceDataType],
458 indexwspTraceDataType++;
459 }
460 Array<OneD, DataType> Fwdreslt;
461 Fwdreslt = wspTraceDataType[indexwspTraceDataType], indexwspTraceDataType++;
462
463 for (size_t m = 0; m < nvariables; ++m)
464 {
465 for (size_t i = 0; i < nTracePts; ++i)
466 {
467 Fwdarray[m][i] = DataType(Fwd[m][i]);
468 }
469 }
470 for (size_t m = 0; m < nvariables; ++m)
471 {
472 Vmath::Zero(nTracePts, Fwdreslt, 1);
473 for (size_t n = 0; n < nvariables; ++n)
474 {
475 for (size_t p = 0; p < nTracePts; ++p)
476 {
477 Fwdreslt[p] += TraceJacArray[0][m][n][p] * Fwdarray[n][p];
478 }
479 }
480 for (size_t i = 0; i < nTracePts; ++i)
481 {
482 FwdFlux[m][i] = NekDouble(Fwdreslt[i]);
483 }
484 }
485
486 for (size_t m = 0; m < nvariables; ++m)
487 {
488 for (size_t i = 0; i < nTracePts; ++i)
489 {
490 Fwdarray[m][i] = DataType(Bwd[m][i]);
491 }
492 }
493 for (size_t m = 0; m < nvariables; ++m)
494 {
495 Vmath::Zero(nTracePts, Fwdreslt, 1);
496 for (size_t n = 0; n < nvariables; ++n)
497 {
498 for (size_t p = 0; p < nTracePts; ++p)
499 {
500 Fwdreslt[p] += TraceJacArray[1][m][n][p] * Fwdarray[n][p];
501 }
502 }
503 for (size_t i = 0; i < nTracePts; ++i)
504 {
505 BwdFlux[m][i] = NekDouble(Fwdreslt[i]);
506 }
507 }
508
509 for (size_t i = 0; i < nvariables; ++i)
510 {
511 Vmath::Fill(nCoeffs, 0.0, outarray[i], 1);
512 timer.Start();
513 pFields[i]->AddTraceIntegralToOffDiag(FwdFlux[i], BwdFlux[i],
514 outarray[i]);
515 timer.Stop();
516 timer.AccumulateRegion("ExpList::AddTraceIntegralToOffDiag", 10);
517 }
518
519 for (size_t i = 0; i < nvariables; ++i)
520 {
521 for (size_t p = 0; p < nCoeffs; ++p)
522 {
523 outarray[i][p] =
524 -m_DtLambdaPreconMat * outarray[i][p] + inarray[i][p];
525 }
526 }
527}
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:106
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:43

References Nektar::LibUtilities::Timer::AccumulateRegion(), Vmath::Fill(), Nektar::PreconCfs::m_DtLambdaPreconMat, Nektar::PreconCfs::m_spacedim, CellMLToNektar.cellml_metadata::p, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Vmath::Zero().

Referenced by v_DoPreconCfs().

◆ PreconBlkDiag()

void Nektar::PreconCfsBRJ::PreconBlkDiag ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
private

Definition at line 322 of file PreconCfsBRJ.cpp.

325{
326 unsigned int nvariables = pFields.size();
327
328 int nTotElmt = pFields[0]->GetNumElmts();
329
330 using vec_t = simd<NekSingle>;
331 const auto vecwidth = vec_t::width;
332
333 // vectorized matrix multiply
334 std::vector<vec_t, tinysimd::allocator<vec_t>> Sinarray(m_max_nblocks);
335 std::vector<vec_t, tinysimd::allocator<vec_t>> Soutarray(m_max_nElmtDof);
336 // std::vector<vec_t, tinysimd::allocator<vec_t>> tmp;
337
338 alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
339
340 for (int ne = 0, cnt = 0, icnt = 0, icnt1 = 0; ne < nTotElmt; ne++)
341 {
342 const auto nElmtCoef = pFields[0]->GetNcoeffs(ne);
343 const auto nElmtDof = nElmtCoef * nvariables;
344 const auto nblocks = (nElmtDof % vecwidth) ? nElmtDof / vecwidth + 1
345 : nElmtDof / vecwidth;
346
347 // gather data into blocks - could probably be done with a
348 // gather call? can be replaced with a gather op when working
349 for (int j = 0; j < nblocks; ++j, icnt += vecwidth)
350 {
351 for (int i = 0; i < vecwidth; ++i)
352 {
353 tmp[i] = inarray[m_inputIdx[icnt + i]];
354 }
355
356 Sinarray[j].load(tmp.data());
357 }
358
359 // Do matrix multiply
360 // first block just needs multiplying
361 vec_t in = Sinarray[0];
362 for (int i = 0; i < nElmtDof; ++i)
363 {
364 Soutarray[i] = m_sBlkDiagMat[cnt++] * in;
365 }
366
367 // rest of blocks are multiply add operations;
368 for (int n = 1; n < nblocks; ++n)
369 {
370 in = Sinarray[n];
371 for (int i = 0; i < nElmtDof; ++i)
372 {
373 Soutarray[i].fma(m_sBlkDiagMat[cnt++], in);
374 }
375 }
376
377 // get block aligned index for this expansion
378 NekSingle val;
379 for (int i = 0; i < nElmtDof; ++i)
380 {
381 // Get hold of data
382 Soutarray[i].store(tmp.data());
383
384 // Sum vector width
385 val = tmp[0];
386 for (int j = 1; j < vecwidth; ++j)
387 {
388 val += tmp[j];
389 }
390 // put data into outarray
391 outarray[m_inputIdx[icnt1 + i]] = NekDouble(val);
392 }
393
394 icnt1 += nblocks * vecwidth;
395 }
396}

References m_inputIdx, m_max_nblocks, m_max_nElmtDof, and m_sBlkDiagMat.

Referenced by v_DoPreconCfs().

◆ 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 
)
overrideprotectedvirtual

Copy array into column major blocks of vector width

Implements Nektar::PreconCfs.

Definition at line 196 of file PreconCfsBRJ.cpp.

200{
201 boost::ignore_unused(pFields);
202
203 if (0 < m_PreconItsStep)
204 {
205 SNekBlkMatSharedPtr PreconMatSingle;
206 using vec_t = simd<NekSingle>;
207 int nvariables = pFields.size();
208 int nelmts = pFields[0]->GetNumElmts();
209 Array<OneD, unsigned int> matdim(nelmts);
210 for (int i = 0; i < nelmts; i++)
211 {
212 matdim[i] = pFields[0]->GetExp(i)->GetNcoeffs() * nvariables;
213 }
214 AllocateNekBlkMatDig(PreconMatSingle, matdim, matdim);
215
217 intmp, m_PreconMatVarsSingle, PreconMatSingle, m_TraceJacSingle,
221
222 if (m_verbose && m_Comm->GetRank() == 0)
223 {
224 cout << " ## CalcuPreconMat " << endl;
225 }
226
227 // copy matrix to simd layout
228 // load matrix
229 int cnt = 0;
230 const auto vecwidth = vec_t::width;
231
232 alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
233
234 for (int ne = 0; ne < nelmts; ne++)
235 {
236 const auto nElmtDof = matdim[ne];
237 const auto nblocks = nElmtDof / vecwidth;
238
239 const NekSingle *mmat =
240 PreconMatSingle->GetBlockPtr(ne, ne)->GetRawPtr();
241 /// Copy array into column major blocks of vector width
242 for (int i1 = 0; i1 < nblocks; ++i1)
243 {
244 for (int j = 0; j < nElmtDof; ++j)
245 {
246 for (int i = 0; i < vecwidth; ++i)
247 {
248 tmp[i] = mmat[j + (i1 * vecwidth + i) * nElmtDof];
249 }
250 // store contiguous vec_t array.
251 m_sBlkDiagMat[cnt++].load(tmp.data());
252 }
253 }
254
255 const auto endwidth = nElmtDof - nblocks * vecwidth;
256
257 // end rows that do not fit into vector widths
258 if (endwidth)
259 {
260 for (int j = 0; j < nElmtDof; ++j)
261 {
262 for (int i = 0; i < endwidth; ++i)
263 {
264 tmp[i] = mmat[j + (nblocks * vecwidth + i) * nElmtDof];
265 }
266
267 for (int i = endwidth; i < vecwidth; ++i)
268 {
269 tmp[i] = 0.0;
270 }
271 m_sBlkDiagMat[cnt++].load(tmp.data());
272 }
273 }
274 }
275 }
276
277 m_BndEvaluateTime = time;
278 m_DtLambdaPreconMat = lambda;
279
280 m_CalcPreconMatFlag = false;
282}
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:97
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:85
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:86
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:87
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:89
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:88
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:90
NekPreconCfsOperators m_operator
Definition: PreconCfs.h:101
int m_PreconTimesCounter
Definition: PreconCfs.h:104
bool m_CalcPreconMatFlag
Definition: PreconCfs.h:109
NekDouble m_BndEvaluateTime
Definition: PreconCfs.h:107
LibUtilities::CommSharedPtr m_Comm
Definition: PreconCfs.h:97
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr

References AllocateNekBlkMatDig(), Nektar::NekPreconCfsOperators::DoCalcPreconMatBRJCoeff(), Nektar::PreconCfs::m_BndEvaluateTime, Nektar::PreconCfs::m_CalcPreconMatFlag, Nektar::PreconCfs::m_Comm, Nektar::PreconCfs::m_DtLambdaPreconMat, Nektar::PreconCfs::m_operator, m_PreconItsStep, m_PreconMatVarsSingle, Nektar::PreconCfs::m_PreconTimesCounter, m_sBlkDiagMat, 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 
)
overrideprotectedvirtual

Implements Nektar::PreconCfs.

Definition at line 86 of file PreconCfsBRJ.cpp.

90{
91 boost::ignore_unused(flag);
92
93 size_t nBRJIterTot = m_PreconItsStep;
94 if (0 == nBRJIterTot)
95 {
96 DoNullPrecon(inarray, outarray, flag);
97 }
98 else
99 {
100 const NekDouble BRJParam = m_BRJRelaxParam;
101 const NekDouble OmBRJParam = 1.0 - BRJParam;
102
103 size_t nvariables = pFields.size();
104 size_t npoints = pFields[0]->GetNcoeffs();
105 size_t ntotpnt = inarray.size();
106
107 ASSERTL0(nvariables * npoints == ntotpnt,
108 "nvariables*npoints!=ntotpnt in PreconCoeff");
109
110 Array<OneD, NekDouble> rhs(ntotpnt);
111
112 Array<OneD, NekDouble> outN(ntotpnt);
113 Array<OneD, NekDouble> outTmp(ntotpnt);
114 Array<OneD, Array<OneD, NekDouble>> rhs2d(nvariables);
115 Array<OneD, Array<OneD, NekDouble>> out_2d(nvariables);
116 Array<OneD, Array<OneD, NekDouble>> outTmp_2d(nvariables);
117 for (size_t m = 0; m < nvariables; m++)
118 {
119 size_t moffset = m * npoints;
120 rhs2d[m] = rhs + moffset;
121 out_2d[m] = outarray + moffset;
122 outTmp_2d[m] = outTmp + moffset;
123 pFields[m]->MultiplyByMassMatrix(inarray + moffset, rhs2d[m]);
124 }
125
126 size_t nphysic = pFields[0]->GetNpoints();
127 size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
128 TensorOfArray3D<NekDouble> qfield(m_spacedim);
129 for (int i = 0; i < m_spacedim; i++)
130 {
131 qfield[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
132 for (size_t j = 0; j < nvariables; j++)
133 {
134 qfield[i][j] = Array<OneD, NekDouble>(nphysic);
135 }
136 }
137 size_t ntmpTrace = 4 + 2 * m_spacedim;
138 TensorOfArray3D<NekDouble> tmpTrace(ntmpTrace);
139 for (size_t i = 0; i < ntmpTrace; i++)
140 {
141 tmpTrace[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
142 for (size_t j = 0; j < nvariables; j++)
143 {
144 tmpTrace[i][j] = Array<OneD, NekDouble>(nTracePts);
145 }
146 }
147 Array<OneD, Array<OneD, NekDouble>> FwdFluxDeriv(nvariables);
148 Array<OneD, Array<OneD, NekDouble>> BwdFluxDeriv(nvariables);
149 for (size_t j = 0; j < nvariables; j++)
150 {
151 FwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
152 BwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
153 }
154
155 bool flagUpdateDervFlux = false;
156
157 const size_t nwspTraceDataType = nvariables + 1;
158 Array<OneD, Array<OneD, NekSingle>> wspTraceDataType(nwspTraceDataType);
159 for (size_t m = 0; m < nwspTraceDataType; m++)
160 {
161 wspTraceDataType[m] = Array<OneD, NekSingle>(nTracePts);
162 }
163
164 LibUtilities::Timer timer;
165 timer.Start();
166 PreconBlkDiag(pFields, rhs, outarray);
167 timer.Stop();
168 timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
169
170 for (size_t nrelax = 0; nrelax < nBRJIterTot - 1; nrelax++)
171 {
172 Vmath::Smul(ntotpnt, OmBRJParam, outarray, 1, outN, 1);
173
174 timer.Start();
176 pFields, nvariables, npoints, rhs2d, out_2d, flagUpdateDervFlux,
177 FwdFluxDeriv, BwdFluxDeriv, qfield, tmpTrace, wspTraceDataType,
180 timer.Stop();
181 timer.AccumulateRegion("PreconCfsBRJ::MinusOffDiag2Rhs", 2);
182
183 timer.Start();
184 PreconBlkDiag(pFields, outarray, outTmp);
185 timer.Stop();
186 timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
187
188 Vmath::Svtvp(ntotpnt, BRJParam, outTmp, 1, outN, 1, outarray, 1);
189 }
190 }
191}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void MinusOffDiag2Rhs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const size_t nvariables, const size_t 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)
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
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:617
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:245

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL0, DoNullPrecon(), m_BRJRelaxParam, m_PreconItsStep, Nektar::PreconCfs::m_spacedim, m_TraceIPSymJacArraySingle, m_TraceJacArraySingle, m_TraceJacDerivArraySingle, m_TraceJacDerivSignSingle, MinusOffDiag2Rhs(), PreconBlkDiag(), Vmath::Smul(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Vmath::Svtvp().

◆ v_InitObject()

void Nektar::PreconCfsBRJ::v_InitObject ( )
overrideprotectedvirtual

Implements Nektar::PreconCfs.

Definition at line 79 of file PreconCfsBRJ.cpp.

80{
81}

◆ v_UpdatePreconMatCheck()

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

Implements Nektar::PreconCfs.

Definition at line 287 of file PreconCfsBRJ.cpp.

289{
290 boost::ignore_unused(res);
291
292 bool flag = false;
293
294 if (m_CalcPreconMatFlag || (m_DtLambdaPreconMat != dtLambda))
295 {
296 flag = true;
297 }
298
300 {
301 flag = true;
302 }
303
304 m_CalcPreconMatFlag = flag;
305 return flag;
306}
int m_PreconMatFreezNumb
Definition: PreconCfs.h:103

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

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:198
static PreconCfsSharedPtr create(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Creates an instance of this class.
Definition: PreconCfsBRJ.h:55
PreconCfsFactory & GetPreconCfsFactory()
Declaration of the boundary condition factory singleton.
Definition: PreconCfs.cpp:42

Name of the class.

Definition at line 66 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_inputIdx

std::vector<int> Nektar::PreconCfsBRJ::m_inputIdx
protected

Definition at line 83 of file PreconCfsBRJ.h.

Referenced by AllocateSIMDPreconBlkMatDiag(), and PreconBlkDiag().

◆ m_max_nblocks

unsigned int Nektar::PreconCfsBRJ::m_max_nblocks
protected

Definition at line 79 of file PreconCfsBRJ.h.

Referenced by AllocateSIMDPreconBlkMatDiag(), and PreconBlkDiag().

◆ m_max_nElmtDof

unsigned int Nektar::PreconCfsBRJ::m_max_nElmtDof
protected

Definition at line 80 of file PreconCfsBRJ.h.

Referenced by AllocateSIMDPreconBlkMatDiag(), and PreconBlkDiag().

◆ 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_PreconMatStorage

PrecType Nektar::PreconCfsBRJ::m_PreconMatStorage
protected

Definition at line 92 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_sBlkDiagMat

std::vector<simd<NekSingle>, tinysimd::allocator<simd<NekSingle> > > Nektar::PreconCfsBRJ::m_sBlkDiagMat
protected

Definition at line 82 of file PreconCfsBRJ.h.

Referenced by AllocateSIMDPreconBlkMatDiag(), PreconBlkDiag(), and v_BuildPreconCfs().

◆ m_TraceIPSymJacArraySingle

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

Definition at line 90 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacArraySingle

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

Definition at line 86 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacDerivArraySingle

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

Definition at line 88 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 89 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs(), and v_DoPreconCfs().

◆ m_TraceJacDerivSingle

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

Definition at line 87 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs().

◆ m_TraceJacSingle

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

Definition at line 85 of file PreconCfsBRJ.h.

Referenced by v_BuildPreconCfs().