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
43
44namespace Nektar
45{
46/**
47 *
48 */
49class CFSImplicit : virtual public CompressibleFlowSystem
50
51{
52public:
55
56 ~CFSImplicit() override = default;
57
58protected:
62
63 int m_nPadding = 1;
64
65 // Flag to control the update of preconditioning matrix.
67 int m_TotLinIts = 0;
69
70 /// Estimate the magnitude of each conserved varibles
72
74
75 /// coefff of spacial derivatives(rhs or m_F in GLM) in calculating the
76 /// residual of the whole equation(used in unsteady time integrations)
80
83
85
87
88 // flag to update shock capturing artificial viscosity. this flag is
89 // switched on/off in DoImplicitSolve() to ensure that the AV is only
90 // updated once every stage in a multi-stage time integration scheme
92
93 void v_InitObject(bool DeclareFields = true) override;
94
96
97 void v_DoSolve() override;
98
99 void v_PrintStatusInformation(const int step,
100 const NekDouble cpuTime) override;
101
102 void v_PrintSummaryStatistics(const NekDouble intTime) override;
103
104 void v_ALEInitObject(
105 int spaceDim,
107
110 const bool &flag);
111
113 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
114 Array<OneD, Array<OneD, NekDouble>> &out, const bool &flag);
115
116 void DoOdeImplicitRhs(
117 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
118 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
119
120 void DoOdeRhsCoeff(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
122 const NekDouble time);
123
124 void DoAdvectionCoeff(
125 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
126 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
127 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
128 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
129
130 void DoDiffusionCoeff(
131 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
133 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
134 const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
135
136 void DoImplicitSolve(
137 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
138 Array<OneD, Array<OneD, NekDouble>> &outpnt, const NekDouble time,
139 const NekDouble lambda);
140
142 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
143 const Array<OneD, const NekDouble> &inarray,
144 Array<OneD, NekDouble> &out, const NekDouble time,
145 const NekDouble lambda);
146
148 const Array<OneD, const NekDouble> &inarray,
149 Array<OneD, NekDouble> &out, const bool &centralDifferenceFlag);
150
151 void CalcRefValues(const Array<OneD, const NekDouble> &inarray);
152
153 void PreconCoeff(const Array<OneD, NekDouble> &inarray,
154 Array<OneD, NekDouble> &outarray, const bool &flag);
155
156 template <typename DataType, typename TypeNekBlkMatSharedPtr>
158 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
159 const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
161 TensorOfArray4D<DataType> &StdMatDataDBB,
162 TensorOfArray5D<DataType> &StdMatDataDBDB);
163
164 template <typename DataType>
166 TensorOfArray5D<DataType> &StdMatDataDBDB);
167
168 template <typename DataType, typename TypeNekBlkMatSharedPtr>
170 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
175 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
176 TensorOfArray5D<DataType> &TraceIPSymJacArray);
177
178 template <typename DataType, typename TypeNekBlkMatSharedPtr>
179 void ElmtVarInvMtrx(
181 TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype);
182
183 template <typename DataType, typename TypeNekBlkMatSharedPtr>
184 void GetTraceJac(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
188 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
189 TensorOfArray5D<DataType> &TraceIPSymJacArray);
190
191 template <typename DataType, typename TypeNekBlkMatSharedPtr>
193 const int nConvectiveFields,
195 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
196 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
197 TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
198 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
199 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
200 TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
201 TensorOfArray5D<DataType> &TraceIPSymJacArray);
202
204 const Array<OneD, NekDouble> &normals,
205 DNekMatSharedPtr &FJac, const NekDouble efix,
206 const NekDouble fsw);
207
208 template <typename DataType, typename TypeNekBlkMatSharedPtr>
209 void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat,
210 TensorOfArray3D<DataType> &MatArray);
211
212 template <typename DataType, typename TypeNekBlkMatSharedPtr>
215 TensorOfArray4D<DataType> &TraceJacDerivArray);
216
217 template <typename DataType, typename TypeNekBlkMatSharedPtr>
220 const DataType valu);
221
222 template <typename DataType, typename TypeNekBlkMatSharedPtr>
224 Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu);
225
227 const Array<OneD, unsigned int> nrow,
228 const Array<OneD, unsigned int> ncol)
229 {
230 mat =
232 SNekMatSharedPtr loc_matNvar;
233 for (int nelm = 0; nelm < nrow.size(); ++nelm)
234 {
235 int nrowsVars = nrow[nelm];
236 int ncolsVars = ncol[nelm];
237
239 nrowsVars, ncolsVars, 0.0);
240 mat->SetBlock(nelm, nelm, loc_matNvar);
241 }
242 }
243
245 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
249 Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
250 TensorOfArray4D<NekSingle> &TraceJacArray,
251 TensorOfArray4D<NekSingle> &TraceJacDerivArray,
252 TensorOfArray5D<NekSingle> &TraceIPSymJacArray);
253
254 template <typename DataType, typename TypeNekBlkMatSharedPtr>
257 const NekDouble dtlamda);
258
260 const int nConvectiveFields, const int nElmtPnt,
261 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
262 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
263 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
264
265 void GetFluxVectorJacPoint(const int nConvectiveFields,
266 const Array<OneD, NekDouble> &conservVar,
267 const Array<OneD, NekDouble> &normals,
268 DNekMatSharedPtr &fluxJac);
269
271 const int nConvectiveFields, const int nDim, const int nPts,
272 const int nTracePts, const NekDouble PenaltyFactor2,
274 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
275 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
276 const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
277 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
278 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
279 const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
280 const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
281 const Array<OneD, NekDouble> &MuVarTrace,
282 Array<OneD, int> &nonZeroIndex,
283 Array<OneD, Array<OneD, NekDouble>> &traceflux);
284
286 const int nConvectiveFields, const int nElmtPnt,
287 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
288 const TensorOfArray3D<NekDouble> &locDerv,
289 const Array<OneD, NekDouble> &locmu,
290 const Array<OneD, NekDouble> &locDmuDT,
291 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
292 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
293 {
294 v_MinusDiffusionFluxJacPoint(nConvectiveFields, nElmtPnt, locVars,
295 locDerv, locmu, locDmuDT, normals, wspMat,
296 PntJacArray);
297 }
298
300 const MultiRegions::ExpListSharedPtr &explist,
301 const Array<OneD, const Array<OneD, NekDouble>> &normals,
302 const int nDervDir,
303 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
304 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir)
305 {
306 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray,
307 ElmtJacArray, nFluxDir);
308 }
309
311 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
312 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
313 const Array<OneD, NekDouble> &locmu,
314 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
315 DNekMatSharedPtr &wspMat,
316 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
317 {
318 v_GetFluxDerivJacDirctnElmt(nConvectiveFields, nElmtPnt, nDervDir,
319 locVars, locmu, locnormal, wspMat,
320 PntJacArray);
321 }
322
324 const MultiRegions::ExpListSharedPtr &explist,
325 const Array<OneD, const Array<OneD, NekDouble>> &normals,
326 const int nDervDir,
327 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
329 {
330 v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray, ElmtJac);
331 }
332
333 void CalcPhysDeriv(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
335 {
336 v_CalcPhysDeriv(inarray, qfield);
337 }
338
339 void CalcMuDmuDT(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
341 {
342 v_CalcMuDmuDT(inarray, mu, DmuDT);
343 }
344
345 virtual void v_DoDiffusionCoeff(
346 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
347 &inarray,
348 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &outarray,
349 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
350 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
351 {
352 }
353
354 virtual void v_CalcMuDmuDT(
355 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
356 &inarray,
357 [[maybe_unused]] Array<OneD, NekDouble> &mu,
358 [[maybe_unused]] Array<OneD, NekDouble> &DmuDT)
359 {
360 }
361
362 virtual void v_CalcPhysDeriv(
363 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>>
364 &inarray,
365 [[maybe_unused]] TensorOfArray3D<NekDouble> &qfield)
366 {
367 }
368
369 virtual void v_MinusDiffusionFluxJacPoint(
370 const int nConvectiveFields, const int nElmtPnt,
371 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
372 const TensorOfArray3D<NekDouble> &locDerv,
373 const Array<OneD, NekDouble> &locmu,
374 const Array<OneD, NekDouble> &locDmuDT,
375 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
376 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
377
378 virtual void v_GetFluxDerivJacDirctn(
379 const MultiRegions::ExpListSharedPtr &explist,
380 const Array<OneD, const Array<OneD, NekDouble>> &normals,
381 const int nDervDir,
382 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
383 TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir);
384
385 virtual void v_GetFluxDerivJacDirctnElmt(
386 const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
387 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
388 const Array<OneD, NekDouble> &locmu,
389 const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
390 DNekMatSharedPtr &wspMat,
391 Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
392
394 const MultiRegions::ExpListSharedPtr &explist,
395 const Array<OneD, const Array<OneD, NekDouble>> &normals,
396 const int nDervDir,
397 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
399
400 bool v_UpdateTimeStepCheck() override;
401};
402} // namespace Nektar
403#endif
void CalcVolJacStdMat(TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void NonlinSysEvaluatorCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &out, const bool &flag)
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)
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 NonlinSysEvaluatorCoeff1D(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag)
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)
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)
void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
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
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)
void MatrixMultiplyMatrixFreeCoeff(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &centralDifferenceFlag)
LibUtilities::NekNonlinSysIterSharedPtr m_nonlinsol
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< NekNonlinSysIter > NekNonlinSysIterSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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