Nektar++
CompressibleFlowSystemImplicit.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: CompressibleFlowSystemImplicit.h
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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Auxiliary functions for the incompressible
32// compressible flow system
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_COMPRESSIBLEFLOWSYSTEMIMPLICIT_H
37#define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_COMPRESSIBLEFLOWSYSTEMIMPLICIT_H
38
40
41namespace Nektar
42{
43/**
44 *
45 */
46class CFSImplicit : virtual public CompressibleFlowSystem
47
48{
49public:
52
53 ~CFSImplicit() override = default;
54
55protected:
59
60 int m_nPadding = 1;
61
62 // Flag to control the update of preconditioning matrix.
64 int m_TotLinIts = 0;
66
67 /// Estimate the magnitude of each conserved varibles
69
71
72 /// coefff of spacial derivatives(rhs or m_F in GLM) in calculating the
73 /// residual of the whole equation(used in unsteady time integrations)
78
81
83
85
86 // flag to update shock capturing artificial viscosity. this flag is
87 // switched on/off in DoImplicitSolve() to ensure that the AV is only
88 // updated once every stage in a multi-stage time integration scheme
90
91 void v_InitObject(bool DeclareFields = true) override;
92
94
95 void v_DoSolve() override;
96
97 void v_PrintStatusInformation(const int step,
98 const NekDouble cpuTime) override;
99
100 void v_PrintSummaryStatistics(const NekDouble intTime) override;
101
104 const bool &flag = true);
105
107 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
109
110 void DoOdeImplicitRhs(
111 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
112 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
113
114 void DoOdeRhsCoeff(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
116 const NekDouble time);
117
118 void DoAdvectionCoeff(
119 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
120 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
121 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
122 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
123
124 void DoDiffusionCoeff(
125 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
127 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
128 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
129
130 void DoImplicitSolve(
131 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
132 Array<OneD, Array<OneD, NekDouble>> &outpnt, const NekDouble time,
133 const NekDouble lambda);
134
136 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
137 const Array<OneD, const NekDouble> &inarray,
138 Array<OneD, NekDouble> &out, const NekDouble time,
139 const NekDouble lambda);
140
142 const Array<OneD, const NekDouble> &inarray,
143 Array<OneD, NekDouble> &out, const bool &flag = false);
144
145 void CalcRefValues(const Array<OneD, const NekDouble> &inarray);
146
147 void PreconCoeff(const Array<OneD, NekDouble> &inarray,
148 Array<OneD, NekDouble> &outarray, const bool &flag);
149
150 template <typename DataType, typename TypeNekBlkMatSharedPtr>
152 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
153 const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
155 TensorOfArray4D<DataType> &StdMatDataDBB,
156 TensorOfArray5D<DataType> &StdMatDataDBDB);
157
158 template <typename DataType>
160 TensorOfArray5D<DataType> &StdMatDataDBDB);
161
162 template <typename DataType, typename TypeNekBlkMatSharedPtr>
164 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
169 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
170 TensorOfArray5D<DataType> &TraceIPSymJacArray);
171
172 template <typename DataType, typename TypeNekBlkMatSharedPtr>
173 void ElmtVarInvMtrx(
175 TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype);
176
177 template <typename DataType, typename TypeNekBlkMatSharedPtr>
178 void GetTraceJac(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
182 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
183 TensorOfArray5D<DataType> &TraceIPSymJacArray);
184
185 template <typename DataType, typename TypeNekBlkMatSharedPtr>
187 const int nConvectiveFields,
189 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
190 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
191 TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
192 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
193 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
194 TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
195 TensorOfArray5D<DataType> &TraceIPSymJacArray);
196
198 const Array<OneD, NekDouble> &normals,
199 DNekMatSharedPtr &FJac, const NekDouble efix,
200 const NekDouble fsw);
201
202 template <typename DataType, typename TypeNekBlkMatSharedPtr>
203 void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat,
204 TensorOfArray3D<DataType> &MatArray);
205
206 template <typename DataType, typename TypeNekBlkMatSharedPtr>
209 TensorOfArray4D<DataType> &TraceJacDerivArray);
210
211 template <typename DataType, typename TypeNekBlkMatSharedPtr>
214 const DataType valu);
215
216 template <typename DataType, typename TypeNekBlkMatSharedPtr>
218 Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu);
219
221 const Array<OneD, unsigned int> nrow,
222 const Array<OneD, unsigned int> ncol)
223 {
224 mat =
226 SNekMatSharedPtr loc_matNvar;
227 for (int nelm = 0; nelm < nrow.size(); ++nelm)
228 {
229 int nrowsVars = nrow[nelm];
230 int ncolsVars = ncol[nelm];
231
233 nrowsVars, ncolsVars, 0.0);
234 mat->SetBlock(nelm, nelm, loc_matNvar);
235 }
236 }
237
239 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
243 Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
244 TensorOfArray4D<NekSingle> &TraceJacArray,
245 TensorOfArray4D<NekSingle> &TraceJacDerivArray,
246 TensorOfArray5D<NekSingle> &TraceIPSymJacArray);
247
248 template <typename DataType, typename TypeNekBlkMatSharedPtr>
251 const NekDouble dtlamda);
252
254 const int nConvectiveFields, const int nElmtPnt,
255 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
256 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
257 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
258
259 void GetFluxVectorJacPoint(const int nConvectiveFields,
260 const Array<OneD, NekDouble> &conservVar,
261 const Array<OneD, NekDouble> &normals,
262 DNekMatSharedPtr &fluxJac);
263
265 const int nConvectiveFields, const int nDim, const int nPts,
266 const int nTracePts, const NekDouble PenaltyFactor2,
268 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
269 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
270 const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
271 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
272 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
273 const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
274 const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
275 const Array<OneD, NekDouble> &MuVarTrace,
276 Array<OneD, int> &nonZeroIndex,
277 Array<OneD, Array<OneD, NekDouble>> &traceflux);
278
280 const int nConvectiveFields, const int nElmtPnt,
281 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
282 const TensorOfArray3D<NekDouble> &locDerv,
283 const Array<OneD, NekDouble> &locmu,
284 const Array<OneD, NekDouble> &locDmuDT,
285 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
286 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
287 {
288 v_MinusDiffusionFluxJacPoint(nConvectiveFields, nElmtPnt, locVars,
289 locDerv, locmu, locDmuDT, normals, wspMat,
290 PntJacArray);
291 }
292
294 const MultiRegions::ExpListSharedPtr &explist,
295 const Array<OneD, const Array<OneD, NekDouble>> &normals,
296 const int nDervDir,
297 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
298 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir)
299 {
300 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray,
301 ElmtJacArray, nFluxDir);
302 }
303
305 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
306 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
307 const Array<OneD, NekDouble> &locmu,
308 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
309 DNekMatSharedPtr &wspMat,
310 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
311 {
312 v_GetFluxDerivJacDirctnElmt(nConvectiveFields, nElmtPnt, nDervDir,
313 locVars, locmu, locnormal, wspMat,
314 PntJacArray);
315 }
316
318 const MultiRegions::ExpListSharedPtr &explist,
319 const Array<OneD, const Array<OneD, NekDouble>> &normals,
320 const int nDervDir,
321 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
323 {
324 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray, ElmtJac);
325 }
326
327 void CalcPhysDeriv(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
329 {
330 v_CalcPhysDeriv(inarray, qfield);
331 }
332
333 void CalcMuDmuDT(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
335 {
336 v_CalcMuDmuDT(inarray, mu, DmuDT);
337 }
338
339 virtual void v_DoDiffusionCoeff(
340 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
341 &inarray,
342 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &outarray,
343 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
344 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
345 {
346 }
347
348 virtual void v_CalcMuDmuDT(
349 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
350 &inarray,
351 [[maybe_unused]] Array<OneD, NekDouble> &mu,
352 [[maybe_unused]] Array<OneD, NekDouble> &DmuDT)
353 {
354 }
355
356 virtual void v_CalcPhysDeriv(
357 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
358 &inarray,
359 [[maybe_unused]] TensorOfArray3D<NekDouble> &qfield)
360 {
361 }
362
363 virtual void v_MinusDiffusionFluxJacPoint(
364 const int nConvectiveFields, const int nElmtPnt,
365 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
366 const TensorOfArray3D<NekDouble> &locDerv,
367 const Array<OneD, NekDouble> &locmu,
368 const Array<OneD, NekDouble> &locDmuDT,
369 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
370 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
371
372 virtual void v_GetFluxDerivJacDirctn(
373 const MultiRegions::ExpListSharedPtr &explist,
374 const Array<OneD, const Array<OneD, NekDouble>> &normals,
375 const int nDervDir,
376 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
377 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir);
378
379 virtual void v_GetFluxDerivJacDirctnElmt(
380 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
381 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
382 const Array<OneD, NekDouble> &locmu,
383 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
384 DNekMatSharedPtr &wspMat,
385 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
386
388 const MultiRegions::ExpListSharedPtr &explist,
389 const Array<OneD, const Array<OneD, NekDouble>> &normals,
390 const int nDervDir,
391 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
393
394 bool v_UpdateTimeStepCheck() override;
395};
396} // namespace Nektar
397#endif
void CalcVolJacStdMat(TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void NumCalcRiemFluxJac(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &AdvVel, const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, const NekDouble &time, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac, TensorOfArray5D< DataType > &TraceIPSymJacArray)
void DoOdeImplicitRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void AddMatNSBlkDiagVol(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const TensorOfArray2D< NekDouble > > &qfield, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void GetFluxDerivJacDirctnElmt(const int nConvectiveFields, const int nElmtPnt, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &locmu, const Array< OneD, const Array< OneD, NekDouble > > &locnormal, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void CalcTraceNumericalFlux(const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &AdvVel, const Array< OneD, const Array< OneD, NekDouble > > &inarray, const NekDouble time, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, const Array< OneD, const TensorOfArray2D< NekDouble > > &qFwd, const Array< OneD, const TensorOfArray2D< NekDouble > > &qBwd, const Array< OneD, NekDouble > &MuVarTrace, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &traceflux)
LibUtilities::NekNonlinSysSharedPtr m_nonlinsol
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, Array< OneD, Array< OneD, NekDouble > > &outpnt, const NekDouble time, const NekDouble lambda)
void DoImplicitSolveCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const NekDouble time, const NekDouble lambda)
void GetFluxVectorJacDirElmt(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void ElmtVarInvMtrx(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype)
virtual void v_CalcMuDmuDT(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
void CalcPreconMatBRJCoeff(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)
virtual void v_MinusDiffusionFluxJacPoint(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const TensorOfArray3D< NekDouble > &locDerv, const Array< OneD, NekDouble > &locmu, const Array< OneD, NekDouble > &locDmuDT, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void PointFluxJacobianPoint(const Array< OneD, NekDouble > &Fwd, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &FJac, const NekDouble efix, const NekDouble fsw)
void v_DoSolve() override
Solves an unsteady problem.
~CFSImplicit() override=default
CFSImplicit(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void GetTraceJac(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, DataType > > &TraceJacDerivSign, TensorOfArray5D< DataType > &TraceIPSymJacArray)
void PreconCoeff(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const bool &flag)
void MultiplyElmtInvMassPlusSource(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const NekDouble dtlamda)
void GetFluxDerivJacDirctn(const MultiRegions::ExpListSharedPtr &explist, const Array< OneD, const Array< OneD, NekDouble > > &normals, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray5D< NekDouble > &ElmtJacArray, const int nFluxDir)
void CalcRefValues(const Array< OneD, const NekDouble > &inarray)
void TransTraceJacMatToArray(const Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, TensorOfArray4D< DataType > &TraceJacDerivArray)
void v_PrintStatusInformation(const int step, const NekDouble cpuTime) override
Print Status Information.
void CalcPhysDeriv(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield)
void MinusDiffusionFluxJacPoint(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const TensorOfArray3D< NekDouble > &locDerv, const Array< OneD, NekDouble > &locmu, const Array< OneD, NekDouble > &locDmuDT, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void CalcMuDmuDT(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
void NonlinSysEvaluatorCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &out)
Array< OneD, Array< OneD, NekDouble > > m_solutionPhys
void DoOdeRhsCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
void DoDiffusionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
Add the diffusions terms to the right-hand side Similar to DoDiffusion() but with outarray in coeffic...
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
virtual void v_GetFluxDerivJacDirctnElmt(const int nConvectiveFields, const int nElmtPnt, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &locmu, const Array< OneD, const Array< OneD, NekDouble > > &locnormal, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
TensorOfArray4D< NekSingle > m_stdSMatDataDBB
void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat, TensorOfArray3D< DataType > &MatArray)
virtual void v_GetFluxDerivJacDirctn(const MultiRegions::ExpListSharedPtr &explist, const Array< OneD, const Array< OneD, NekDouble > > &normals, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray5D< NekDouble > &ElmtJacArray, const int nFluxDir)
void v_PrintSummaryStatistics(const NekDouble intTime) override
Print Summary Statistics.
virtual void v_DoDiffusionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
void GetFluxDerivJacDirctn(const MultiRegions::ExpListSharedPtr &explist, const Array< OneD, const Array< OneD, NekDouble > > &normals, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, DNekMatSharedPtr > > &ElmtJac)
void GetFluxVectorJacPoint(const int nConvectiveFields, const Array< OneD, NekDouble > &conservVar, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &fluxJac)
TensorOfArray5D< NekSingle > m_stdSMatDataDBDB
void NonlinSysEvaluatorCoeff1D(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag=true)
void MatrixMultiplyMatrixFreeCoeff(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag=false)
NekDouble m_TimeIntegLambda
coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(us...
void Fill2DArrayOfBlkDiagonalMat(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const DataType valu)
Array< OneD, NekDouble > m_magnitdEstimat
Estimate the magnitude of each conserved varibles.
void AddMatNSBlkDiagBnd(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray2D< TypeNekBlkMatSharedPtr > &gmtxarray, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, DataType > > &TraceJacDerivSign, TensorOfArray5D< DataType > &TraceIPSymJacArray)
virtual void v_CalcPhysDeriv(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield)
void Fill1DArrayOfBlkDiagonalMat(Array< OneD, TypeNekBlkMatSharedPtr > &gmtxarray, const DataType valu)
virtual void v_GetFluxDerivJacDirctn(const MultiRegions::ExpListSharedPtr &explist, const Array< OneD, const Array< OneD, NekDouble > > &normals, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, DNekMatSharedPtr > > &ElmtJac)
void v_InitObject(bool DeclareFields=true) override
Initialization object for CFSImplicit class.
void DoAdvectionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
Compute the advection terms for the right-hand side.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< NekNonlinSys > NekNonlinSysSharedPtr
Definition: NekNonlinSys.h:46
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
std::shared_ptr< SNekMat > SNekMatSharedPtr
std::shared_ptr< PreconCfs > PreconCfsSharedPtr
Definition: PreconCfs.h:112
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble