DESModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2015 OpenFOAM Foundation
9  Copyright (C) 2022 Upstream CFD GmbH
10  Copyright (C) 2016, 2022 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Class
29  Foam::LESModels::DESModel
30 
31 Description
32  Templated abstract base class for DES models
33 
34 SourceFiles
35  DESModel.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef LESModels_DESModel_H
40 #define LESModels_DESModel_H
41 
42 #include "DESModelBase.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace LESModels
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class DESModel Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class BasicTurbulenceModel>
56 class DESModel
57 :
58  public DESModelBase,
59  public LESeddyViscosity<BasicTurbulenceModel>
60 {
61  // Private Member Functions
62 
63  //- No copy construct
64  DESModel(const DESModel&) = delete;
65 
66  //- No copy assignment
67  void operator=(const DESModel&) = delete;
68 
69 
70 protected:
71 
72  // Protected Data
73 
74  //- Model-specific transition constant
76 
77 
78 public:
79 
80  typedef typename BasicTurbulenceModel::alphaField alphaField;
81  typedef typename BasicTurbulenceModel::rhoField rhoField;
82  typedef typename BasicTurbulenceModel::transportModel transportModel;
83 
84  // Constructors
85 
86  //- Construct from components
87  DESModel
88  (
89  const word& type,
90  const alphaField& alpha,
91  const rhoField& rho,
92  const volVectorField& U,
93  const surfaceScalarField& alphaRhoPhi,
94  const surfaceScalarField& phi,
95  const transportModel& transport,
96  const word& propertiesName = turbulenceModel::propertiesName
97  );
98 
99 
100  //- Destructor
101  virtual ~DESModel() = default;
102 
103 
104  // Public Member Functions
105 
106  //- Re-read model coefficients if they have changed
107  virtual bool read();
108 
109  //- Return the LES field indicator
110  virtual tmp<volScalarField> LESRegion() const = 0;
111 
112  //- Return modified strain rate
114  (
115  const volTensorField& gradU,
116  const dimensionedScalar& coeff
117  );
118 
119  //- Return modified strain rate
120  // Note: uses Ctrans_ coefficient
121  virtual tmp<volScalarField> Ssigma(const volTensorField& gradU) const;
122 
123  //- Return the shielding function
124  virtual tmp<volScalarField> fd() const;
125 };
126 
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 } // End namespace LESModels
131 } // End namespace Foam
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 #ifdef NoRepository
136  #include "DESModel.C"
137 #endif
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 #endif
142 
143 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual tmp< volScalarField > LESRegion() const =0
Return the LES field indicator.
BasicTurbulenceModel::rhoField rhoField
Definition: DESModel.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
Base class for DES models providing an interfaces to DES fields.
Definition: DESModelBase.H:49
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< volScalarField > Ssigma(const volTensorField &gradU, const dimensionedScalar &coeff)
Return modified strain rate.
Definition: DESModel.C:82
dimensionedScalar Ctrans_
Model-specific transition constant.
Definition: DESModel.H:76
BasicTurbulenceModel::transportModel transportModel
Definition: DESModel.H:83
virtual ~DESModel()=default
Destructor.
U
Definition: pEqn.H:72
virtual tmp< volScalarField > fd() const
Return the shielding function.
Definition: DESModel.C:162
BasicTurbulenceModel::alphaField alphaField
Definition: DESModel.H:81
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool read()
Re-read model coefficients if they have changed.
Definition: DESModel.C:67
Templated abstract base class for DES models.
Definition: DESModel.H:51
Namespace for OpenFOAM.