Nektar++
DiffusionLFRNS.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DiffusionLFRNS.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: LFRNS diffusion class.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERUTILS_DIFFUSIONLFRNS
36#define NEKTAR_SOLVERUTILS_DIFFUSIONLFRNS
37
39
40namespace Nektar::SolverUtils
41{
43{
44public:
45 static DiffusionSharedPtr create(std::string diffType)
46 {
47 return DiffusionSharedPtr(new DiffusionLFRNS(diffType));
48 }
49
50 static std::string type[];
51
52protected:
53 DiffusionLFRNS(std::string diffType);
54
55 std::string m_diffType;
56
57 void v_InitObject(
60
61 void v_Diffuse(const std::size_t nConvective,
63 const Array<OneD, Array<OneD, NekDouble>> &inarray,
65 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
66 const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
67
69 {
70 m_homoDerivs = deriv;
71 }
72
74 {
75 return m_viscTensor;
76 }
77
78private:
82
85
90
99
103 std::string m_ViscosityType;
108
120
123
125
128
131
134
135 void NumericalFluxO1(
137 const Array<OneD, Array<OneD, NekDouble>> &inarray,
138 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &numericalFluxO1);
139
140 void WeakPenaltyO1(
142 const Array<OneD, Array<OneD, NekDouble>> &inarray,
143 Array<OneD, Array<OneD, NekDouble>> &penaltyfluxO1);
144
145 void NumericalFluxO2(
147 const Array<OneD, Array<OneD, NekDouble>> &ufield,
150
151 void WeakPenaltyO2(
153 const int var, const int dir,
154 const Array<OneD, const NekDouble> &qfield,
155 Array<OneD, NekDouble> &penaltyflux);
156
157 void DerCFlux_1D(const int nConvectiveFields,
160 const Array<OneD, const NekDouble> &iFlux,
161 Array<OneD, NekDouble> &derCFlux);
162
163 void DerCFlux_2D(const int nConvectiveFields, const int direction,
166 const Array<OneD, NekDouble> &iFlux,
167 Array<OneD, NekDouble> &derCFlux);
168
169 void DivCFlux_2D(const int nConvectiveFields,
171 const Array<OneD, const NekDouble> &fluxX1,
172 const Array<OneD, const NekDouble> &fluxX2,
173 const Array<OneD, const NekDouble> &numericalFlux,
174 Array<OneD, NekDouble> &divCFlux);
175
177 const int nConvectiveFields,
179 const Array<OneD, const NekDouble> &fluxX1,
180 const Array<OneD, const NekDouble> &fluxX2,
181 const Array<OneD, const NekDouble> &numericalFlux,
182 Array<OneD, NekDouble> &divCFlux);
183};
184
185typedef std::shared_ptr<DiffusionLFRNS> DiffusionLFRNSSharedPtr;
186} // namespace Nektar::SolverUtils
187
188#endif
void v_Diffuse(const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Calculate FR Diffusion for the Navier-Stokes (NS) equations using an LDG interface flux.
void SetupMetrics(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Setup the metric terms to compute the contravariant fluxes. (i.e. this special metric terms transform...
Array< OneD, Array< OneD, NekDouble > > m_traceVel
Array< OneD, Array< OneD, NekDouble > > m_divFD
Array< OneD, NekDouble > m_jac
void DerCFlux_1D(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &flux, const Array< OneD, const NekDouble > &iFlux, Array< OneD, NekDouble > &derCFlux)
Compute the derivative of the corrective flux for 1D problems.
void NumericalFluxO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &numericalFluxO1)
Builds the numerical flux for the 1st order derivatives.
Array< OneD, Array< OneD, NekDouble > > m_viscFlux
void WeakPenaltyO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &penaltyfluxO1)
Imposes appropriate bcs for the 1st order derivatives.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tmp2
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e1
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_DFC2
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi2
Array< OneD, Array< OneD, NekDouble > > m_divFC
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_DFC1
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi3
Array< OneD, Array< OneD, NekDouble > > m_gmat
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e3
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e2
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi2
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_viscTensor
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > & v_GetFluxTensor() override
void SetupCFunctions(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Setup the derivatives of the correction functions. For more details see J Sci Comput (2011) 47: 50–72...
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tmp1
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi3
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e0
void DivCFlux_2D(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &fluxX2, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 2D problems.
void DerCFlux_2D(const int nConvectiveFields, const int direction, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &flux, const Array< OneD, NekDouble > &iFlux, Array< OneD, NekDouble > &derCFlux)
Compute the derivative of the corrective flux wrt a given coordinate for 2D problems.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_IF1
void DivCFlux_2D_Gauss(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &fluxX2, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 2D problems where POINTSTYPE="GaussGaussLegendre".
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_D1
DiffusionLFRNS(std::string diffType)
DiffusionLFRNS uses the Flux Reconstruction (FR) approach to compute the diffusion term....
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, Array< OneD, NekDouble > > m_homoDerivs
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_BD1
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi1
void NumericalFluxO2(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
Build the numerical flux for the 2nd order derivatives.
static DiffusionSharedPtr create(std::string diffType)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_DU1
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_DD1
void WeakPenaltyO2(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const int var, const int dir, const Array< OneD, const NekDouble > &qfield, Array< OneD, NekDouble > &penaltyflux)
Imposes appropriate bcs for the 2nd order derivatives.
void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv) override
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Initiliase DiffusionLFRNS objects and store them before starting the time-stepping.
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi1
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:55
std::shared_ptr< DiffusionLFRNS > DiffusionLFRNSSharedPtr
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble