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  size_t nvariables = pFields.size();
64 
67  for (size_t i = 0; i < nvariables; i++)
68  {
70  }
72 
74 }
75 
77 {
79 }
80 
83  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
84  const bool &flag)
85 {
86  boost::ignore_unused(flag);
87 
88  size_t nBRJIterTot = m_PreconItsStep;
89  if (0 == nBRJIterTot)
90  {
91  DoNullPrecon(inarray, outarray, flag);
92  }
93  else
94  {
95  const NekDouble BRJParam = m_BRJRelaxParam;
96  const NekDouble OmBRJParam = 1.0 - BRJParam;
97 
98  size_t nvariables = pFields.size();
99  size_t npoints = pFields[0]->GetNcoeffs();
100  size_t ntotpnt = inarray.size();
101 
102  ASSERTL0(nvariables * npoints == ntotpnt,
103  "nvariables*npoints!=ntotpnt in PreconCoeff");
104 
105  Array<OneD, NekDouble> rhs(ntotpnt);
106 
107  Array<OneD, NekDouble> outN(ntotpnt);
108  Array<OneD, NekDouble> outTmp(ntotpnt);
109  Array<OneD, Array<OneD, NekDouble>> rhs2d(nvariables);
110  Array<OneD, Array<OneD, NekDouble>> out_2d(nvariables);
111  Array<OneD, Array<OneD, NekDouble>> outTmp_2d(nvariables);
112  for (size_t m = 0; m < nvariables; m++)
113  {
114  size_t moffset = m * npoints;
115  rhs2d[m] = rhs + moffset;
116  out_2d[m] = outarray + moffset;
117  outTmp_2d[m] = outTmp + moffset;
118  pFields[m]->MultiplyByMassMatrix(inarray + moffset, rhs2d[m]);
119  }
120 
121  size_t nphysic = pFields[0]->GetNpoints();
122  size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
124  for (int i = 0; i < m_spacedim; i++)
125  {
126  qfield[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
127  for (size_t j = 0; j < nvariables; j++)
128  {
129  qfield[i][j] = Array<OneD, NekDouble>(nphysic);
130  }
131  }
132  size_t ntmpTrace = 4 + 2 * m_spacedim;
133  TensorOfArray3D<NekDouble> tmpTrace(ntmpTrace);
134  for (size_t i = 0; i < ntmpTrace; i++)
135  {
136  tmpTrace[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
137  for (size_t j = 0; j < nvariables; j++)
138  {
139  tmpTrace[i][j] = Array<OneD, NekDouble>(nTracePts);
140  }
141  }
142  Array<OneD, Array<OneD, NekDouble>> FwdFluxDeriv(nvariables);
143  Array<OneD, Array<OneD, NekDouble>> BwdFluxDeriv(nvariables);
144  for (size_t j = 0; j < nvariables; j++)
145  {
146  FwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
147  BwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
148  }
149 
150  bool flagUpdateDervFlux = false;
151 
152  const size_t nwspTraceDataType = nvariables + 1;
153  Array<OneD, Array<OneD, NekSingle>> wspTraceDataType(nwspTraceDataType);
154  for (size_t m = 0; m < nwspTraceDataType; m++)
155  {
156  wspTraceDataType[m] = Array<OneD, NekSingle>(nTracePts);
157  }
158 
159  LibUtilities::Timer timer;
160  timer.Start();
161  PreconBlkDiag(pFields, rhs, outarray);
162  timer.Stop();
163  timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
164 
165  for (size_t nrelax = 0; nrelax < nBRJIterTot - 1; nrelax++)
166  {
167  Vmath::Smul(ntotpnt, OmBRJParam, outarray, 1, outN, 1);
168 
169  timer.Start();
171  pFields, nvariables, npoints, rhs2d, out_2d, flagUpdateDervFlux,
172  FwdFluxDeriv, BwdFluxDeriv, qfield, tmpTrace, wspTraceDataType,
175  timer.Stop();
176  timer.AccumulateRegion("PreconCfsBRJ::MinusOffDiag2Rhs", 2);
177 
178  timer.Start();
179  PreconBlkDiag(pFields, outarray, outTmp);
180  timer.Stop();
181  timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
182 
183  Vmath::Svtvp(ntotpnt, BRJParam, outTmp, 1, outN, 1, outarray, 1);
184  }
185  }
186 }
187 
190  const Array<OneD, const Array<OneD, NekDouble>> &intmp,
191  const NekDouble time, const NekDouble lambda)
192 {
193  boost::ignore_unused(pFields);
194 
195  if (0 < m_PreconItsStep)
196  {
197  SNekBlkMatSharedPtr PreconMatSingle;
198  using vec_t = simd<NekSingle>;
199  int nvariables = pFields.size();
200  int nelmts = pFields[0]->GetNumElmts();
201  Array<OneD, unsigned int> matdim(nelmts);
202  for (int i = 0; i < nelmts; i++)
203  {
204  matdim[i] = pFields[0]->GetExp(i)->GetNcoeffs() * nvariables;
205  }
206  AllocateNekBlkMatDig(PreconMatSingle, matdim, matdim);
207 
209  intmp, m_PreconMatVarsSingle, PreconMatSingle, m_TraceJacSingle,
213 
214  if (m_verbose && m_Comm->GetRank() == 0)
215  {
216  cout << " ## CalcuPreconMat " << endl;
217  }
218 
219  // copy matrix to simd layout
220  // load matrix
221  int cnt = 0;
222  const auto vecwidth = vec_t::width;
223 
224  alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
225 
226  for (int ne = 0; ne < nelmts; ne++)
227  {
228  const auto nElmtDof = matdim[ne];
229  const auto nblocks = nElmtDof / vecwidth;
230 
231  const NekSingle *mmat =
232  PreconMatSingle->GetBlockPtr(ne, ne)->GetRawPtr();
233  /// Copy array into column major blocks of vector width
234  for (int i1 = 0; i1 < nblocks; ++i1)
235  {
236  for (int j = 0; j < nElmtDof; ++j)
237  {
238  for (int i = 0; i < vecwidth; ++i)
239  {
240  tmp[i] = mmat[j + (i1 * vecwidth + i) * nElmtDof];
241  }
242  // store contiguous vec_t array.
243  m_sBlkDiagMat[cnt++].load(tmp.data());
244  }
245  }
246 
247  const auto endwidth = nElmtDof - nblocks * vecwidth;
248 
249  // end rows that do not fit into vector widths
250  if (endwidth)
251  {
252  for (int j = 0; j < nElmtDof; ++j)
253  {
254  for (int i = 0; i < endwidth; ++i)
255  {
256  tmp[i] = mmat[j + (nblocks * vecwidth + i) * nElmtDof];
257  }
258 
259  for (int i = endwidth; i < vecwidth; ++i)
260  {
261  tmp[i] = 0.0;
262  }
263  m_sBlkDiagMat[cnt++].load(tmp.data());
264  }
265  }
266  }
267  }
268 
269  m_BndEvaluateTime = time;
270  m_DtLambdaPreconMat = lambda;
271 
272  m_CalcPreconMatFlag = false;
274 }
275 
277  const Array<OneD, const NekDouble> &res, const NekDouble dtLambda)
278 {
279  boost::ignore_unused(res);
280 
281  bool flag = false;
282 
283  if (m_CalcPreconMatFlag || (m_DtLambdaPreconMat != dtLambda))
284  {
285  flag = true;
286  }
287 
289  {
290  flag = true;
291  }
292 
293  m_CalcPreconMatFlag = flag;
294  return flag;
295 }
296 
299  const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray)
300 {
301  unsigned int nvariables = pFields.size();
302 
303  int nTotElmt = pFields[0]->GetNumElmts();
304 
305  using vec_t = simd<NekSingle>;
306  const auto vecwidth = vec_t::width;
307 
308  // vectorized matrix multiply
309  std::vector<vec_t, tinysimd::allocator<vec_t>> Sinarray(m_max_nblocks);
310  std::vector<vec_t, tinysimd::allocator<vec_t>> Soutarray(m_max_nElmtDof);
311  // std::vector<vec_t, tinysimd::allocator<vec_t>> tmp;
312 
313  alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
314 
315  for (int ne = 0, cnt = 0, icnt = 0, icnt1 = 0; ne < nTotElmt; ne++)
316  {
317  const auto nElmtCoef = pFields[0]->GetNcoeffs(ne);
318  const auto nElmtDof = nElmtCoef * nvariables;
319  const auto nblocks = (nElmtDof % vecwidth) ? nElmtDof / vecwidth + 1
320  : nElmtDof / vecwidth;
321 
322  // gather data into blocks - could probably be done with a
323  // gather call? can be replaced with a gather op when working
324  for (int j = 0; j < nblocks; ++j, icnt += vecwidth)
325  {
326  for (int i = 0; i < vecwidth; ++i)
327  {
328  tmp[i] = inarray[m_inputIdx[icnt + i]];
329  }
330 
331  Sinarray[j].load(tmp.data());
332  }
333 
334  // Do matrix multiply
335  // first block just needs multiplying
336  vec_t in = Sinarray[0];
337  for (int i = 0; i < nElmtDof; ++i)
338  {
339  Soutarray[i] = m_sBlkDiagMat[cnt++] * in;
340  }
341 
342  // rest of blocks are multiply add operations;
343  for (int n = 1; n < nblocks; ++n)
344  {
345  in = Sinarray[n];
346  for (int i = 0; i < nElmtDof; ++i)
347  {
348  Soutarray[i].fma(m_sBlkDiagMat[cnt++], in);
349  }
350  }
351 
352  // get block aligned index for this expansion
353  NekSingle val;
354  for (int i = 0; i < nElmtDof; ++i)
355  {
356  // Get hold of data
357  Soutarray[i].store(tmp.data());
358 
359  // Sum vector width
360  val = tmp[0];
361  for (int j = 1; j < vecwidth; ++j)
362  {
363  val += tmp[j];
364  }
365  // put data into outarray
366  outarray[m_inputIdx[icnt1 + i]] = NekDouble(val);
367  }
368 
369  icnt1 += nblocks * vecwidth;
370  }
371 }
372 
373 template <typename DataType>
376  const size_t nvariables, const size_t nCoeffs,
377  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
378  Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
379  Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
380  Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
382  Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
383  const TensorOfArray4D<DataType> &TraceJacArray,
384  const TensorOfArray4D<DataType> &TraceJacDerivArray,
385  const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
386  const TensorOfArray5D<DataType> &TraceIPSymJacArray)
387 {
388  boost::ignore_unused(flagUpdateDervFlux, qfield, TraceJacDerivArray,
389  TraceJacDerivSign, FwdFluxDeriv, BwdFluxDeriv,
390  TraceIPSymJacArray);
391 
392  size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
393  size_t npoints = pFields[0]->GetNpoints();
394  size_t nDim = m_spacedim;
395 
396  Array<OneD, Array<OneD, NekDouble>> outpnts(nvariables);
397  for (size_t i = 0; i < nvariables; i++)
398  {
399  outpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
400  pFields[i]->BwdTrans(outarray[i], outpnts[i]);
401  }
402 
403  // Store forwards/backwards space along trace space
408  TensorOfArray3D<NekDouble> numDerivBwd(nDim);
409  TensorOfArray3D<NekDouble> numDerivFwd(nDim);
410  size_t indexwspTrace = 0;
411  Fwd = wspTrace[indexwspTrace], indexwspTrace++;
412  Bwd = wspTrace[indexwspTrace], indexwspTrace++;
413  FwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
414  BwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
415 
416  LibUtilities::Timer timer;
417  for (size_t i = 0; i < nvariables; ++i)
418  {
419  timer.Start();
420  pFields[i]->GetFwdBwdTracePhys(outpnts[i], Fwd[i], Bwd[i]);
421  timer.Stop();
422  timer.AccumulateRegion("ExpList::GetFwdBwdTracePhys", 10);
423  }
424 
425  size_t indexwspTraceDataType = 0;
426  Array<OneD, Array<OneD, DataType>> Fwdarray(nvariables);
427  for (size_t m = 0; m < nvariables; ++m)
428  {
429  Fwdarray[m] = wspTraceDataType[indexwspTraceDataType],
430  indexwspTraceDataType++;
431  }
432  Array<OneD, DataType> Fwdreslt;
433  Fwdreslt = wspTraceDataType[indexwspTraceDataType], indexwspTraceDataType++;
434 
435  for (size_t m = 0; m < nvariables; ++m)
436  {
437  for (size_t i = 0; i < nTracePts; ++i)
438  {
439  Fwdarray[m][i] = DataType(Fwd[m][i]);
440  }
441  }
442  for (size_t m = 0; m < nvariables; ++m)
443  {
444  Vmath::Zero(nTracePts, Fwdreslt, 1);
445  for (size_t n = 0; n < nvariables; ++n)
446  {
447  for (size_t p = 0; p < nTracePts; ++p)
448  {
449  Fwdreslt[p] += TraceJacArray[0][m][n][p] * Fwdarray[n][p];
450  }
451  }
452  for (size_t i = 0; i < nTracePts; ++i)
453  {
454  FwdFlux[m][i] = NekDouble(Fwdreslt[i]);
455  }
456  }
457 
458  for (size_t m = 0; m < nvariables; ++m)
459  {
460  for (size_t i = 0; i < nTracePts; ++i)
461  {
462  Fwdarray[m][i] = DataType(Bwd[m][i]);
463  }
464  }
465  for (size_t m = 0; m < nvariables; ++m)
466  {
467  Vmath::Zero(nTracePts, Fwdreslt, 1);
468  for (size_t n = 0; n < nvariables; ++n)
469  {
470  for (size_t p = 0; p < nTracePts; ++p)
471  {
472  Fwdreslt[p] += TraceJacArray[1][m][n][p] * Fwdarray[n][p];
473  }
474  }
475  for (size_t i = 0; i < nTracePts; ++i)
476  {
477  BwdFlux[m][i] = NekDouble(Fwdreslt[i]);
478  }
479  }
480 
481  for (size_t i = 0; i < nvariables; ++i)
482  {
483  Vmath::Fill(nCoeffs, 0.0, outarray[i], 1);
484  timer.Start();
485  pFields[i]->AddTraceIntegralToOffDiag(FwdFlux[i], BwdFlux[i],
486  outarray[i]);
487  timer.Stop();
488  timer.AccumulateRegion("ExpList::AddTraceIntegralToOffDiag", 10);
489  }
490 
491  for (size_t i = 0; i < nvariables; ++i)
492  {
493  for (size_t p = 0; p < nCoeffs; ++p)
494  {
495  outarray[i][p] =
496  -m_DtLambdaPreconMat * outarray[i][p] + inarray[i][p];
497  }
498  }
499 }
500 
501 template <typename TypeNekBlkMatSharedPtr>
505  const int &nscale)
506 {
507 
508  size_t nvars = pFields.size();
509  size_t nelmts = pFields[0]->GetNumElmts();
510  size_t nelmtcoef;
511  Array<OneD, unsigned int> nelmtmatdim(nelmts);
512  for (size_t i = 0; i < nelmts; i++)
513  {
514  nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
515  nelmtmatdim[i] = nelmtcoef * nscale;
516  }
517 
518  for (size_t i = 0; i < nvars; i++)
519  {
520  for (size_t j = 0; j < nvars; j++)
521  {
522  AllocateNekBlkMatDig(gmtxarray[i][j], nelmtmatdim, nelmtmatdim);
523  }
524  }
525 }
526 
527 } // 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
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:88
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:89
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:80
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:90
unsigned int m_max_nblocks
Definition: PreconCfsBRJ.h:82
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:92
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:95
virtual void v_InitObject() override
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:91
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 AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const int &nscale=1)
unsigned int m_max_nElmtDof
Definition: PreconCfsBRJ.h:83
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:93
std::vector< int > m_inputIdx
Definition: PreconCfsBRJ.h:86
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:241
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
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:143
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:85
int m_PreconTimesCounter
Definition: PreconCfs.h:87
int m_PreconMatFreezNumb
Definition: PreconCfs.h:86
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
Definition: PreconCfs.cpp:55
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:89
bool m_CalcPreconMatFlag
Definition: PreconCfs.h:92
NekDouble m_BndEvaluateTime
Definition: PreconCfs.h:90
LibUtilities::CommSharedPtr m_Comm
Definition: PreconCfs.h:82
virtual void v_InitObject() override
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:2
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
@ eDiagonal
Definition: PreconCfs.h:48
tinysimd::simd< NekDouble > vec_t
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 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
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80