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);
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  void PreconCoeff(const Array<OneD, NekDouble> &inarray,
117  Array<OneD, NekDouble> &outarray, const bool &flag);
118 
119  template <typename DataType, typename TypeNekBlkMatSharedPtr>
120  void AddMatNSBlkDiagVol(
121  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
122  const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
124  TensorOfArray4D<DataType> &StdMatDataDBB,
125  TensorOfArray5D<DataType> &StdMatDataDBDB);
126 
127  template <typename DataType>
128  void CalcVolJacStdMat(TensorOfArray4D<DataType> &StdMatDataDBB,
129  TensorOfArray5D<DataType> &StdMatDataDBDB);
130 
131  template <typename DataType, typename TypeNekBlkMatSharedPtr>
132  void AddMatNSBlkDiagBnd(
133  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
138  Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
139  TensorOfArray5D<DataType> &TraceIPSymJacArray);
140 
141  template <typename DataType, typename TypeNekBlkMatSharedPtr>
142  void ElmtVarInvMtrx(
144  TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype);
145 
146  template <typename DataType, typename TypeNekBlkMatSharedPtr>
147  void GetTraceJac(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
151  Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
152  TensorOfArray5D<DataType> &TraceIPSymJacArray);
153 
154  template <typename DataType, typename TypeNekBlkMatSharedPtr>
155  void NumCalcRiemFluxJac(
156  const int nConvectiveFields,
158  const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
159  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
160  TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
161  const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
162  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
163  TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
164  TensorOfArray5D<DataType> &TraceIPSymJacArray);
165 
167  const Array<OneD, NekDouble> &normals,
168  DNekMatSharedPtr &FJac, const NekDouble efix,
169  const NekDouble fsw);
170 
171  template <typename DataType, typename TypeNekBlkMatSharedPtr>
172  void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat,
173  TensorOfArray3D<DataType> &MatArray);
174 
175  template <typename DataType, typename TypeNekBlkMatSharedPtr>
177  const Array<OneD, TypeNekBlkMatSharedPtr> &TraceJac,
178  const Array<OneD, TypeNekBlkMatSharedPtr> &TraceJacDeriv,
179  TensorOfArray4D<DataType> &TraceJacArray,
180  TensorOfArray4D<DataType> &TraceJacDerivArray);
181 
182  template <typename DataType, typename TypeNekBlkMatSharedPtr>
185  const DataType valu);
186 
187  template <typename DataType, typename TypeNekBlkMatSharedPtr>
189  Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu);
190 
192  const Array<OneD, unsigned int> nrow,
193  const Array<OneD, unsigned int> ncol)
194  {
195  mat =
197  SNekMatSharedPtr loc_matNvar;
198  for (int nelm = 0; nelm < nrow.size(); ++nelm)
199  {
200  int nrowsVars = nrow[nelm];
201  int ncolsVars = ncol[nelm];
202 
204  nrowsVars, ncolsVars, 0.0);
205  mat->SetBlock(nelm, nelm, loc_matNvar);
206  }
207  }
208 
210  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
213  Array<OneD, SNekBlkMatSharedPtr> &TraceJacDeriv,
214  Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
215  TensorOfArray4D<NekSingle> &TraceJacArray,
216  TensorOfArray4D<NekSingle> &TraceJacDerivArray,
217  TensorOfArray5D<NekSingle> &TraceIPSymJacArray);
218 
219  template <typename DataType, typename TypeNekBlkMatSharedPtr>
222  const NekDouble dtlamda, const DataType tmpDataType);
223 
225  const int nConvectiveFields, const int nElmtPnt,
226  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
227  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
228  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
229 
230  void GetFluxVectorJacPoint(const int nConvectiveFields,
231  const Array<OneD, NekDouble> &conservVar,
232  const Array<OneD, NekDouble> &normals,
233  DNekMatSharedPtr &fluxJac);
234 
236  const int nConvectiveFields, const int nDim, const int nPts,
237  const int nTracePts, const NekDouble PenaltyFactor2,
239  const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
240  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
241  const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
242  const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
243  const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
244  const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
245  const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
246  const Array<OneD, NekDouble> &MuVarTrace,
247  Array<OneD, int> &nonZeroIndex,
248  Array<OneD, Array<OneD, NekDouble>> &traceflux);
249 
251  const int nConvectiveFields, const int nElmtPnt,
252  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
253  const TensorOfArray3D<NekDouble> &locDerv,
254  const Array<OneD, NekDouble> &locmu,
255  const Array<OneD, NekDouble> &locDmuDT,
256  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
257  Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
258  {
259  v_MinusDiffusionFluxJacPoint(nConvectiveFields, nElmtPnt, locVars,
260  locDerv, locmu, locDmuDT, normals, wspMat,
261  PntJacArray);
262  }
263 
265  const MultiRegions::ExpListSharedPtr &explist,
266  const Array<OneD, const Array<OneD, NekDouble>> &normals,
267  const int nDervDir,
268  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
269  TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir)
270  {
271  v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray,
272  ElmtJacArray, nFluxDir);
273  }
274 
276  const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
277  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
278  const Array<OneD, NekDouble> &locmu,
279  const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
280  DNekMatSharedPtr &wspMat,
281  Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
282  {
283  v_GetFluxDerivJacDirctnElmt(nConvectiveFields, nElmtPnt, nDervDir,
284  locVars, locmu, locnormal, wspMat,
285  PntJacArray);
286  }
287 
289  const MultiRegions::ExpListSharedPtr &explist,
290  const Array<OneD, const Array<OneD, NekDouble>> &normals,
291  const int nDervDir,
292  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
294  {
295  v_GetFluxDerivJacDirctn(explist, normals, nDervDir, inarray, ElmtJac);
296  }
297 
298  void CalcPhysDeriv(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
300  {
301  v_CalcPhysDeriv(inarray, qfield);
302  }
303 
304  void DoDiffusionCoeff(
305  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
306  Array<OneD, Array<OneD, NekDouble>> &outarray,
307  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
308  const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
310  const Array<OneD, const NekDouble> &inarray,
311  Array<OneD, NekDouble> &out, const bool &flag = false);
312 
313  void CalcMuDmuDT(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
315  {
316  v_CalcMuDmuDT(inarray, mu, DmuDT);
317  }
318 
319  virtual void v_DoDiffusionCoeff(
320  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
321  Array<OneD, Array<OneD, NekDouble>> &outarray,
322  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
323  const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
324  {
325  boost::ignore_unused(inarray, outarray, pFwd, pBwd);
326  }
327 
328  virtual void v_CalcMuDmuDT(
329  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
331  {
332  boost::ignore_unused(inarray, mu, DmuDT);
333  }
334 
335  virtual void v_CalcPhysDeriv(
336  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
338  {
339  boost::ignore_unused(inarray, qfield);
340  }
341 
342  virtual void v_MinusDiffusionFluxJacPoint(
343  const int nConvectiveFields, const int nElmtPnt,
344  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
345  const TensorOfArray3D<NekDouble> &locDerv,
346  const Array<OneD, NekDouble> &locmu,
347  const Array<OneD, NekDouble> &locDmuDT,
348  const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
349  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
350 
351  virtual void v_GetFluxDerivJacDirctn(
352  const MultiRegions::ExpListSharedPtr &explist,
353  const Array<OneD, const Array<OneD, NekDouble>> &normals,
354  const int nDervDir,
355  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
356  TensorOfArray5D<NekDouble> &ElmtJacArray, const int nFluxDir);
357 
358  virtual void v_GetFluxDerivJacDirctnElmt(
359  const int nConvectiveFields, const int nElmtPnt, const int nDervDir,
360  const Array<OneD, const Array<OneD, NekDouble>> &locVars,
361  const Array<OneD, NekDouble> &locmu,
362  const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
363  DNekMatSharedPtr &wspMat,
364  Array<OneD, Array<OneD, NekDouble>> &PntJacArray);
365 
367  const MultiRegions::ExpListSharedPtr &explist,
368  const Array<OneD, const Array<OneD, NekDouble>> &normals,
369  const int nDervDir,
370  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
372 
373  virtual bool UpdateTimeStepCheck();
374 };
375 } // namespace Nektar
376 #endif
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)
virtual void v_InitObject(bool DeclareFields=true)
Initialization object for CFSImplicit class.
void Fill1DArrayOfBlkDiagonalMat(Array< OneD, TypeNekBlkMatSharedPtr > &gmtxarray, const DataType valu)
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:1
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