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
39using namespace std;
40
41namespace Nektar
42{
43/**
44 * @class PreconCfsBRJ
45 *
46 * Solves a linear system using iterative methods.
47 */
48std::string PreconCfsBRJ::className =
50 "PreconCfsBRJ", PreconCfsBRJ::create,
51 "Block Relaxed Jacobi Preconditioner for CFS.");
52
56 const LibUtilities::CommSharedPtr &vComm)
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
67 for (size_t i = 0; i < nvariables; i++)
68 {
70 }
72
74}
75
76/**
77 *
78 */
80{
81}
82
83/**
84 *
85 */
88 const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
89 const bool &flag)
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();
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
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}
192
193/**
194 *
195 */
198 const Array<OneD, const Array<OneD, NekDouble>> &intmp,
199 const NekDouble time, const NekDouble lambda)
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}
283
284/**
285 *
286 */
288 const Array<OneD, const NekDouble> &res, const NekDouble dtLambda)
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}
307
308/**
309 *
310 */
312 Array<OneD, NekDouble> &pOutput,
313 const bool &flag)
314{
315 boost::ignore_unused(flag);
316 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
317}
318
319/**
320 *
321 */
324 const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray)
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}
397
398/**
399 *
400 */
401template <typename DataType>
404 const size_t nvariables, const size_t nCoeffs,
405 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
406 Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
407 Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
408 Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
410 Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
411 const TensorOfArray4D<DataType> &TraceJacArray,
412 const TensorOfArray4D<DataType> &TraceJacDerivArray,
413 const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
414 const TensorOfArray5D<DataType> &TraceIPSymJacArray)
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
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
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}
528
529/**
530 *
531 */
532template <typename TypeNekBlkMatSharedPtr>
536 const int &nscale)
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}
557
558} // 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:72
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
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)
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 void v_DoPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag) override
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:86
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
static std::string className
Name of the class.
Definition: PreconCfsBRJ.h:66
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:87
unsigned int m_max_nblocks
Definition: PreconCfsBRJ.h:79
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:89
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
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:92
virtual void v_InitObject() override
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:88
PreconCfsBRJ(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
unsigned int m_max_nElmtDof
Definition: PreconCfsBRJ.h:80
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const int &nscale=1)
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:90
std::vector< int > m_inputIdx
Definition: PreconCfsBRJ.h:83
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:244
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
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual bool v_UpdatePreconMatCheck(const Array< OneD, const NekDouble > &res, const NekDouble dtLambda) override
std::vector< simd< NekSingle >, tinysimd::allocator< simd< NekSingle > > > m_sBlkDiagMat
Definition: PreconCfsBRJ.h:82
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
NekPreconCfsOperators m_operator
Definition: PreconCfs.h:101
int m_PreconTimesCounter
Definition: PreconCfs.h:104
int m_PreconMatFreezNumb
Definition: PreconCfs.h:103
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:106
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< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
PreconCfsFactory & GetPreconCfsFactory()
Declaration of the boundary condition factory singleton.
Definition: PreconCfs.cpp:42
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
@ eDiagonal
Definition: PreconCfs.h:49
tinysimd::simd< NekDouble > vec_t
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: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
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
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80