Nektar++
DiffusionLDG.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: DiffusionLDG.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: LDG diffusion class.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
36 #define NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
37 
38 #include <boost/core/ignore_unused.hpp>
39 
41 
42 namespace Nektar
43 {
44 namespace SolverUtils
45 {
46 class DiffusionLDG : public Diffusion
47 {
48 public:
49  static DiffusionSharedPtr create(std::string diffType)
50  {
51  boost::ignore_unused(diffType);
52  return DiffusionSharedPtr(new DiffusionLDG());
53  }
54 
55  static std::string type;
56 
57 protected:
58  DiffusionLDG();
59 
60  std::string m_shockCaptureType;
61 
62  /// Coefficient of penalty term
64 
67 
68  virtual void v_InitObject(
71 
72  virtual void v_Diffuse(
73  const std::size_t nConvective,
75  const Array<OneD, Array<OneD, NekDouble>> &inarray,
76  Array<OneD, Array<OneD, NekDouble>> &outarray,
77  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
78  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
79 
80  virtual void v_DiffuseCoeffs(
81  const std::size_t nConvective,
83  const Array<OneD, Array<OneD, NekDouble>> &inarray,
84  Array<OneD, Array<OneD, NekDouble>> &outarray,
85  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
86  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
87 
88  virtual void v_DiffuseCalcDerivative(
90  const Array<OneD, Array<OneD, NekDouble>> &inarray,
92  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
93  const Array<OneD, Array<OneD, NekDouble>> &pBwd) override;
94 
95  virtual void v_DiffuseVolumeFlux(
97  const Array<OneD, Array<OneD, NekDouble>> &inarray,
99  TensorOfArray3D<NekDouble> &VolumeFlux,
100  Array<OneD, int> &nonZeroIndex) override;
101 
102  virtual void v_DiffuseTraceFlux(
104  const Array<OneD, Array<OneD, NekDouble>> &inarray,
106  TensorOfArray3D<NekDouble> &VolumeFlux,
107  Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
108  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
109  const Array<OneD, Array<OneD, NekDouble>> &pBwd,
110  Array<OneD, int> &nonZeroIndex) override;
111 
112  void NumFluxforScalar(
114  const Array<OneD, Array<OneD, NekDouble>> &ufield,
116  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
117  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
118 
119  void ApplyScalarBCs(
121  const std::size_t var, const Array<OneD, const NekDouble> &ufield,
122  const Array<OneD, const NekDouble> &Fwd,
123  const Array<OneD, const NekDouble> &Bwd,
124  Array<OneD, NekDouble> &penaltyflux);
125 
126  void NumFluxforVector(
128  const Array<OneD, Array<OneD, NekDouble>> &ufield,
131 
132  void ApplyVectorBCs(
134  const std::size_t var, const std::size_t dir,
135  const Array<OneD, const NekDouble> &qfield,
136  const Array<OneD, const NekDouble> &qFwd,
137  const Array<OneD, const NekDouble> &qBwd,
138  Array<OneD, NekDouble> &penaltyflux);
139 };
140 } // namespace SolverUtils
141 } // namespace Nektar
142 
143 #endif
virtual void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void ApplyScalarBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionLDG.h:49
virtual void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
void ApplyVectorBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:66
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
virtual void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63
virtual 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
virtual void v_DiffuseCoeffs(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
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble >> &qflux)
Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to p...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:550
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble