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#include <boost/core/ignore_unused.hpp>
41
42namespace Nektar
43{
44/**
45 *
46 */
47class CFSImplicit : virtual public CompressibleFlowSystem
48
49{
50public:
53
54 virtual ~CFSImplicit();
55
56 virtual void v_InitObject(bool DeclareFields = true) override;
57
58 virtual void v_DoSolve() override;
59
60 virtual void v_PrintStatusInformation(const int step,
61 const NekDouble cpuTime) override;
62
63 virtual void v_PrintSummaryStatistics(const NekDouble intTime) override;
64
66
69 const bool &flag);
70
72 const Array<OneD, const NekDouble> &inarray,
73 Array<OneD, NekDouble> &out, const bool &flag = true,
75
77 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
79 const Array<OneD, const Array<OneD, NekDouble>> &source =
81
82 void DoOdeRhsCoeff(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
84 const NekDouble time);
85
87 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
88 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
89 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
90 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
91
92 void DoImplicitSolve(
93 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
94 Array<OneD, Array<OneD, NekDouble>> &outpnt, const NekDouble time,
95 const NekDouble lambda);
96
98 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
99 const Array<OneD, const NekDouble> &inarray,
100 Array<OneD, NekDouble> &out, const NekDouble time,
101 const NekDouble lambda);
102
103 void CalcRefValues(const Array<OneD, const NekDouble> &inarray);
104
105protected:
109
110 int m_nPadding = 1;
111
112 // Flag to control the update of preconditioning matrix.
114 int m_TotLinIts = 0;
116
117 /// Estimate the magnitude of each conserved varibles
119
121
122 /// coefff of spacial derivatives(rhs or m_F in GLM) in calculating the
123 /// residual of the whole equation(used in unsteady time integrations)
128
131
133
135
136 // flag to update shock capturing artificial viscosity. this flag is
137 // switched on/off in DoImplicitSolve() to ensure that the AV is only
138 // updated once every stage in a multi-stage time integration scheme
140
141 void PreconCoeff(const Array<OneD, NekDouble> &inarray,
142 Array<OneD, NekDouble> &outarray, const bool &flag);
143
144 template <typename DataType, typename TypeNekBlkMatSharedPtr>
146 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
147 const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
149 TensorOfArray4D<DataType> &StdMatDataDBB,
150 TensorOfArray5D<DataType> &StdMatDataDBDB);
151
152 template <typename DataType>
154 TensorOfArray5D<DataType> &StdMatDataDBDB);
155
156 template <typename DataType, typename TypeNekBlkMatSharedPtr>
158 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
163 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
164 TensorOfArray5D<DataType> &TraceIPSymJacArray);
165
166 template <typename DataType, typename TypeNekBlkMatSharedPtr>
167 void ElmtVarInvMtrx(
169 TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype);
170
171 template <typename DataType, typename TypeNekBlkMatSharedPtr>
172 void GetTraceJac(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
176 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
177 TensorOfArray5D<DataType> &TraceIPSymJacArray);
178
179 template <typename DataType, typename TypeNekBlkMatSharedPtr>
181 const int nConvectiveFields,
183 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
184 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
185 TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
186 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
187 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
188 TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
189 TensorOfArray5D<DataType> &TraceIPSymJacArray);
190
192 const Array<OneD, NekDouble> &normals,
193 DNekMatSharedPtr &FJac, const NekDouble efix,
194 const NekDouble fsw);
195
196 template <typename DataType, typename TypeNekBlkMatSharedPtr>
197 void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat,
198 TensorOfArray3D<DataType> &MatArray);
199
200 template <typename DataType, typename TypeNekBlkMatSharedPtr>
203 const Array<OneD, TypeNekBlkMatSharedPtr> &TraceJacDeriv,
204 TensorOfArray4D<DataType> &TraceJacArray,
205 TensorOfArray4D<DataType> &TraceJacDerivArray);
206
207 template <typename DataType, typename TypeNekBlkMatSharedPtr>
210 const DataType valu);
211
212 template <typename DataType, typename TypeNekBlkMatSharedPtr>
214 Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu);
215
217 const Array<OneD, unsigned int> nrow,
218 const Array<OneD, unsigned int> ncol)
219 {
220 mat =
222 SNekMatSharedPtr loc_matNvar;
223 for (int nelm = 0; nelm < nrow.size(); ++nelm)
224 {
225 int nrowsVars = nrow[nelm];
226 int ncolsVars = ncol[nelm];
227
229 nrowsVars, ncolsVars, 0.0);
230 mat->SetBlock(nelm, nelm, loc_matNvar);
231 }
232 }
233
235 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
239 Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
240 TensorOfArray4D<NekSingle> &TraceJacArray,
241 TensorOfArray4D<NekSingle> &TraceJacDerivArray,
242 TensorOfArray5D<NekSingle> &TraceIPSymJacArray);
243
244 template <typename DataType, typename TypeNekBlkMatSharedPtr>
247 const NekDouble dtlamda, const DataType tmpDataType);
248
250 const int nConvectiveFields, const int nElmtPnt,
251 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
252 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
253 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
254
255 void GetFluxVectorJacPoint(const int nConvectiveFields,
256 const Array<OneD, NekDouble> &conservVar,
257 const Array<OneD, NekDouble> &normals,
258 DNekMatSharedPtr &fluxJac);
259
261 const int nConvectiveFields, const int nDim, const int nPts,
262 const int nTracePts, const NekDouble PenaltyFactor2,
264 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
265 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
266 const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
267 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
268 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
269 const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
270 const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
271 const Array<OneD, NekDouble> &MuVarTrace,
272 Array<OneD, int> &nonZeroIndex,
273 Array<OneD, Array<OneD, NekDouble>> &traceflux);
274
276 const int nConvectiveFields, const int nElmtPnt,
277 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
278 const TensorOfArray3D<NekDouble> &locDerv,
279 const Array<OneD, NekDouble> &locmu,
280 const Array<OneD, NekDouble> &locDmuDT,
281 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
282 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
283 {
284 v_MinusDiffusionFluxJacPoint(nConvectiveFields, nElmtPnt, locVars,
285 locDerv, locmu, locDmuDT, normals, wspMat,
286 PntJacArray);
287 }
288
290 const MultiRegions::ExpListSharedPtr &explist,
291 const Array<OneD, const Array<OneD, NekDouble>> &normals,
292 const int nDervDir,
293 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
294 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir)
295 {
296 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray,
297 ElmtJacArray, nFluxDir);
298 }
299
301 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
302 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
303 const Array<OneD, NekDouble> &locmu,
304 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
305 DNekMatSharedPtr &wspMat,
306 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
307 {
308 v_GetFluxDerivJacDirctnElmt(nConvectiveFields, nElmtPnt, nDervDir,
309 locVars, locmu, locnormal, wspMat,
310 PntJacArray);
311 }
312
314 const MultiRegions::ExpListSharedPtr &explist,
315 const Array<OneD, const Array<OneD, NekDouble>> &normals,
316 const int nDervDir,
317 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
319 {
320 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray, ElmtJac);
321 }
322
323 void CalcPhysDeriv(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
325 {
326 v_CalcPhysDeriv(inarray, qfield);
327 }
328
329 void DoDiffusionCoeff(
330 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
332 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
333 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
335 const Array<OneD, const NekDouble> &inarray,
336 Array<OneD, NekDouble> &out, const bool &flag = false);
337
338 void CalcMuDmuDT(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
340 {
341 v_CalcMuDmuDT(inarray, mu, DmuDT);
342 }
343
344 virtual void v_DoDiffusionCoeff(
345 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
347 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
348 const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
349 {
350 boost::ignore_unused(inarray, outarray, pFwd, pBwd);
351 }
352
353 virtual void v_CalcMuDmuDT(
354 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
356 {
357 boost::ignore_unused(inarray, mu, DmuDT);
358 }
359
360 virtual void v_CalcPhysDeriv(
361 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
363 {
364 boost::ignore_unused(inarray, qfield);
365 }
366
367 virtual void v_MinusDiffusionFluxJacPoint(
368 const int nConvectiveFields, const int nElmtPnt,
369 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
370 const TensorOfArray3D<NekDouble> &locDerv,
371 const Array<OneD, NekDouble> &locmu,
372 const Array<OneD, NekDouble> &locDmuDT,
373 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
374 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
375
376 virtual void v_GetFluxDerivJacDirctn(
377 const MultiRegions::ExpListSharedPtr &explist,
378 const Array<OneD, const Array<OneD, NekDouble>> &normals,
379 const int nDervDir,
380 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
381 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir);
382
383 virtual void v_GetFluxDerivJacDirctnElmt(
384 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
385 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
386 const Array<OneD, NekDouble> &locmu,
387 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
388 DNekMatSharedPtr &wspMat,
389 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
390
392 const MultiRegions::ExpListSharedPtr &explist,
393 const Array<OneD, const Array<OneD, NekDouble>> &normals,
394 const int nDervDir,
395 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
397
398 virtual bool v_UpdateTimeStepCheck() override;
399};
400} // namespace Nektar
401#endif
virtual bool v_UpdateTimeStepCheck() override
void MultiplyElmtInvMassPlusSource(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const NekDouble dtlamda, const DataType tmpDataType)
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 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)
virtual void v_DoSolve() override
Solves an unsteady problem.
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 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 NonlinSysEvaluatorCoeff1D(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag)
void CalcRefValues(const Array< OneD, const NekDouble > &inarray)
virtual 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)
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 NonlinSysEvaluatorCoeff(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag=true, const Array< OneD, const NekDouble > &source=NullNekDouble1DArray)
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)
virtual void v_PrintSummaryStatistics(const NekDouble intTime) override
Print Summary Statistics.
virtual ~CFSImplicit()
Destructor for CFSImplicit class.
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 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 TransTraceJacMatToArray(const Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, const Array< OneD, TypeNekBlkMatSharedPtr > &TraceJacDeriv, TensorOfArray4D< DataType > &TraceJacArray, TensorOfArray4D< DataType > &TraceJacDerivArray)
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)
virtual 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:48
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
std::shared_ptr< SNekMat > SNekMatSharedPtr
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
std::shared_ptr< PreconCfs > PreconCfsSharedPtr
Definition: PreconCfs.h:126
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble