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 
42 namespace Nektar
43 {
44 /**
45  *
46  */
47 class CFSImplicit : virtual public CompressibleFlowSystem
48 
49 {
50 public:
53 
54  virtual ~CFSImplicit();
55 
56  virtual void v_InitObject(bool DeclareFields = true) override;
57 
59 
62  const bool &flag);
63 
65  const Array<OneD, const NekDouble> &inarray,
66  Array<OneD, NekDouble> &out, const bool &flag = true,
68 
70  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
72  const Array<OneD, const Array<OneD, NekDouble>> &source =
74 
75  void DoOdeRhsCoeff(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
76  Array<OneD, Array<OneD, NekDouble>> &outarray,
77  const NekDouble time);
78 
79  void DoAdvectionCoeff(
80  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
81  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
82  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
83  const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
84 
85  void DoImplicitSolve(
86  const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
87  Array<OneD, Array<OneD, NekDouble>> &outpnt, const NekDouble time,
88  const NekDouble lambda);
89 
91  const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
92  const Array<OneD, const NekDouble> &inarray,
93  Array<OneD, NekDouble> &out, const NekDouble time,
94  const NekDouble lambda);
95 
96  void CalcRefValues(const Array<OneD, const NekDouble> &inarray);
97 
98 protected:
101 
102  int m_nPadding = 1;
103 
105 
108 
111 
113 
115 
116  // flag to update shock capturing artificial viscosity. this flag is
117  // switched on/off in DoImplicitSolve() to ensure that the AV is only
118  // updated once every stage in a multi-stage time integration scheme
120 
121  void PreconCoeff(const Array<OneD, NekDouble> &inarray,
122  Array<OneD, NekDouble> &outarray, const bool &flag);
123 
124  template <typename DataType, typename TypeNekBlkMatSharedPtr>
125  void AddMatNSBlkDiagVol(
126  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
127  const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
129  TensorOfArray4D<DataType> &StdMatDataDBB,
130  TensorOfArray5D<DataType> &StdMatDataDBDB);
131 
132  template <typename DataType>
133  void CalcVolJacStdMat(TensorOfArray4D<DataType> &StdMatDataDBB,
134  TensorOfArray5D<DataType> &StdMatDataDBDB);
135 
136  template <typename DataType, typename TypeNekBlkMatSharedPtr>
137  void AddMatNSBlkDiagBnd(
138  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
143  Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
144  TensorOfArray5D<DataType> &TraceIPSymJacArray);
145 
146  template <typename DataType, typename TypeNekBlkMatSharedPtr>
147  void ElmtVarInvMtrx(
149  TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype);
150 
151  template <typename DataType, typename TypeNekBlkMatSharedPtr>
152  void GetTraceJac(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
156  Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
157  TensorOfArray5D<DataType> &TraceIPSymJacArray);
158 
159  template <typename DataType, typename TypeNekBlkMatSharedPtr>
160  void NumCalcRiemFluxJac(
161  const int nConvectiveFields,
163  const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
164  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
165  TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
166  const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
167  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
168  TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
169  TensorOfArray5D<DataType> &TraceIPSymJacArray);
170 
172  const Array<OneD, NekDouble> &normals,
173  DNekMatSharedPtr &FJac, const NekDouble efix,
174  const NekDouble fsw);
175 
176  template <typename DataType, typename TypeNekBlkMatSharedPtr>
177  void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat,
178  TensorOfArray3D<DataType> &MatArray);
179 
180  template <typename DataType, typename TypeNekBlkMatSharedPtr>
182  const Array<OneD, TypeNekBlkMatSharedPtr> &TraceJac,
183  const Array<OneD, TypeNekBlkMatSharedPtr> &TraceJacDeriv,
184  TensorOfArray4D<DataType> &TraceJacArray,
185  TensorOfArray4D<DataType> &TraceJacDerivArray);
186 
187  template <typename DataType, typename TypeNekBlkMatSharedPtr>
190  const DataType valu);
191 
192  template <typename DataType, typename TypeNekBlkMatSharedPtr>
194  Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu);
195 
197  const Array<OneD, unsigned int> nrow,
198  const Array<OneD, unsigned int> ncol)
199  {
200  mat =
202  SNekMatSharedPtr loc_matNvar;
203  for (int nelm = 0; nelm < nrow.size(); ++nelm)
204  {
205  int nrowsVars = nrow[nelm];
206  int ncolsVars = ncol[nelm];
207 
209  nrowsVars, ncolsVars, 0.0);
210  mat->SetBlock(nelm, nelm, loc_matNvar);
211  }
212  }
213 
215  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
218  Array<OneD, SNekBlkMatSharedPtr> &TraceJacDeriv,
219  Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
220  TensorOfArray4D<NekSingle> &TraceJacArray,
221  TensorOfArray4D<NekSingle> &TraceJacDerivArray,
222  TensorOfArray5D<NekSingle> &TraceIPSymJacArray);
223 
224  template <typename DataType, typename TypeNekBlkMatSharedPtr>
227  const NekDouble dtlamda, const DataType tmpDataType);
228 
230  const int nConvectiveFields, const int nElmtPnt,
231  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
232  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
233  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
234 
235  void GetFluxVectorJacPoint(const int nConvectiveFields,
236  const Array<OneD, NekDouble> &conservVar,
237  const Array<OneD, NekDouble> &normals,
238  DNekMatSharedPtr &fluxJac);
239 
241  const int nConvectiveFields, const int nDim, const int nPts,
242  const int nTracePts, const NekDouble PenaltyFactor2,
244  const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
245  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
246  const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
247  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
248  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
249  const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
250  const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
251  const Array<OneD, NekDouble> &MuVarTrace,
252  Array<OneD, int> &nonZeroIndex,
253  Array<OneD, Array<OneD, NekDouble>> &traceflux);
254 
256  const int nConvectiveFields, const int nElmtPnt,
257  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
258  const TensorOfArray3D<NekDouble> &locDerv,
259  const Array<OneD, NekDouble> &locmu,
260  const Array<OneD, NekDouble> &locDmuDT,
261  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
262  Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
263  {
264  v_MinusDiffusionFluxJacPoint(nConvectiveFields, nElmtPnt, locVars,
265  locDerv, locmu, locDmuDT, normals, wspMat,
266  PntJacArray);
267  }
268 
270  const MultiRegions::ExpListSharedPtr &explist,
271  const Array<OneD, const Array<OneD, NekDouble>> &normals,
272  const int nDervDir,
273  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
274  TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir)
275  {
276  v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray,
277  ElmtJacArray, nFluxDir);
278  }
279 
281  const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
282  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
283  const Array<OneD, NekDouble> &locmu,
284  const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
285  DNekMatSharedPtr &wspMat,
286  Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
287  {
288  v_GetFluxDerivJacDirctnElmt(nConvectiveFields, nElmtPnt, nDervDir,
289  locVars, locmu, locnormal, 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,
299  {
300  v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray, ElmtJac);
301  }
302 
303  void CalcPhysDeriv(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
305  {
306  v_CalcPhysDeriv(inarray, qfield);
307  }
308 
309  void DoDiffusionCoeff(
310  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
311  Array<OneD, Array<OneD, NekDouble>> &outarray,
312  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
313  const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
315  const Array<OneD, const NekDouble> &inarray,
316  Array<OneD, NekDouble> &out, const bool &flag = false);
317 
318  void CalcMuDmuDT(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
320  {
321  v_CalcMuDmuDT(inarray, mu, DmuDT);
322  }
323 
324  virtual void v_DoDiffusionCoeff(
325  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
326  Array<OneD, Array<OneD, NekDouble>> &outarray,
327  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
328  const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
329  {
330  boost::ignore_unused(inarray, outarray, pFwd, pBwd);
331  }
332 
333  virtual void v_CalcMuDmuDT(
334  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
336  {
337  boost::ignore_unused(inarray, mu, DmuDT);
338  }
339 
340  virtual void v_CalcPhysDeriv(
341  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
343  {
344  boost::ignore_unused(inarray, qfield);
345  }
346 
347  virtual void v_MinusDiffusionFluxJacPoint(
348  const int nConvectiveFields, const int nElmtPnt,
349  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
350  const TensorOfArray3D<NekDouble> &locDerv,
351  const Array<OneD, NekDouble> &locmu,
352  const Array<OneD, NekDouble> &locDmuDT,
353  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
354  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
355 
356  virtual void v_GetFluxDerivJacDirctn(
357  const MultiRegions::ExpListSharedPtr &explist,
358  const Array<OneD, const Array<OneD, NekDouble>> &normals,
359  const int nDervDir,
360  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
361  TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir);
362 
363  virtual void v_GetFluxDerivJacDirctnElmt(
364  const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
365  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
366  const Array<OneD, NekDouble> &locmu,
367  const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
368  DNekMatSharedPtr &wspMat,
369  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
370 
372  const MultiRegions::ExpListSharedPtr &explist,
373  const Array<OneD, const Array<OneD, NekDouble>> &normals,
374  const int nDervDir,
375  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
377 
378  virtual bool v_UpdateTimeStepCheck() override;
379 };
380 } // namespace Nektar
381 #endif
virtual bool v_UpdateTimeStepCheck() override
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 CalcVolJacStdMat(TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void ElmtVarInvMtrx(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype)
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.
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)
LibUtilities::NekNonlinSysSharedPtr m_nonlinsol
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 PointFluxJacobianPoint(const Array< OneD, NekDouble > &Fwd, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &FJac, const NekDouble efix, const NekDouble fsw)
CFSImplicit(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void PreconCoeff(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const bool &flag)
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 MultiplyElmtInvMassPlusSource(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const NekDouble dtlamda, const DataType tmpDataType)
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_CalcPhysDeriv(const Array< OneD, const Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield)
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 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 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 Fill2DArrayOfBlkDiagonalMat(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const DataType valu)
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 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)
Array< OneD, Array< OneD, NekDouble > > m_solutionPhys
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
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)
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)
virtual ~CFSImplicit()
Destructor for CFSImplicit class.
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 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 GetFluxVectorJacPoint(const int nConvectiveFields, const Array< OneD, NekDouble > &conservVar, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &fluxJac)
TensorOfArray5D< NekSingle > m_stdSMatDataDBDB
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)
virtual void v_CalcMuDmuDT(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
void MatrixMultiplyMatrixFreeCoeff(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag=false)
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 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)
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 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)
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 DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble >> &inpnts, Array< OneD, Array< OneD, NekDouble >> &outpnt, const NekDouble time, const NekDouble lambda)
void Fill1DArrayOfBlkDiagonalMat(Array< OneD, TypeNekBlkMatSharedPtr > &gmtxarray, const DataType valu)
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for CFSImplicit class.
void CalcMuDmuDT(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
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)
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:172
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
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< PreconCfsOp > PreconCfsOpSharedPtr
Definition: PreconCfsOp.h:152
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble