basic.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::PDRDragModels::basic
26 
27 Description
28  Basic sub-grid obstacle drag model.
29  Details supplied by J Puttock 2/7/06.
30 
31  <b> Sub-grid drag term </b>
32 
33  The resistance term (force per unit of volume) is given by:
34 
35  \f[
36  R = -\frac{1}{2} \rho \vert \dwea{\vec{U}} \vert \dwea{\vec{U}}.D
37  \f]
38 
39  where:
40 
41  \f$ D \f$ is the tensor field "CR" in \f$ m^{-1} \f$
42 
43  This is term is treated implicitly in UEqn.H
44 
45  <b> Sub-grid turbulence generation </b>
46 
47  The turbulence source term \f$ G_{R} \f$ occurring in the
48  \f$ \kappa-\epsilon \f$ equations for the generation of turbulence due
49  to interaction with unresolved obstacles :
50 
51  \f$ G_{R} = C_{s}\beta_{\nu}
52  \mu_{eff} A_{w}^{2}(\dwea{\vec{U}}-\dwea{\vec{U}_{s}})^2 + \frac{1}{2}
53  \rho \vert \dwea{\vec{U}} \vert \dwea{\vec{U}}.T.\dwea{\vec{U}} \f$
54 
55  where:
56 
57  \f$ C_{s} \f$ = 1
58 
59  \f$ \beta_{\nu} \f$ is the volume porosity (file "betav").
60 
61  \f$ \mu_{eff} \f$ is the effective viscosity.
62 
63  \f$ A_{w}^{2}\f$ is the obstacle surface area per unit of volume
64  (file "Aw").
65 
66  \f$ \dwea{\vec{U}_{s}} \f$ is the slip velocity and is considered
67  \f$ \frac{1}{2}. \dwea{\vec{U}} \f$.
68 
69  \f$ T \f$ is a tensor in the file CT.
70 
71  The term \f$ G_{R} \f$ is treated explicitly in the \f$ \kappa-\epsilon
72  \f$ Eqs in the \link PDRkEpsilon.C \endlink file.
73 
74 
75 SourceFiles
76  basic.C
77 
78 \*---------------------------------------------------------------------------*/
79 
80 #ifndef basic_H
81 #define basic_H
82 
83 #include "PDRDragModel.H"
84 #include "XiEqModel.H"
85 
86 
87 namespace Foam
88 {
89 namespace PDRDragModels
90 {
91 
92 /*---------------------------------------------------------------------------*\
93  Class basic Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 class basic
97 :
98  public PDRDragModel
99 {
100  // Private Data
101 
102  dimensionedScalar Csu;
103  dimensionedScalar Csk;
104 
105  volScalarField Aw_;
106  volSymmTensorField CR_;
107 
108 
109 public:
110 
111  //- Runtime type information
112  TypeName("basic");
113 
114 
115  // Constructors
116 
117  //- Construct from components
118  basic
119  (
120  const dictionary& PDRProperties,
122  const volScalarField& rho,
123  const volVectorField& U,
124  const surfaceScalarField& phi
125  );
126 
127  //- Disallow default bitwise copy construction
128  basic(const basic&);
129 
130 
131  //- Destructor
132  virtual ~basic();
133 
134 
135  // Member Functions
136 
137  //- Return the momentum drag coefficient
138  virtual tmp<volSymmTensorField> Dcu() const;
139 
140  //- Return the momentum drag turbulence generation rate
141  virtual tmp<volScalarField> Gk() const;
142 
143  //- Update properties from given dictionary
144  virtual bool read(const dictionary& PDRProperties);
145 
146  //- Write fields
147  void writeFields() const;
148 
149 
150  // Member Operators
151 
152  //- Disallow default bitwise assignment
153  void operator=(const basic&) = delete;
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace PDRDragModels
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #endif
165 
166 // ************************************************************************* //
virtual bool read()
Read object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Base-class for sub-grid obstacle drag models. The available drag model is at basic.H.
Definition: PDRDragModel.H:53
virtual ~basic()
Destructor.
virtual tmp< volSymmTensorField > Dcu() const
Return the momentum drag coefficient.
phi
Definition: pEqn.H:104
void writeFields() const
Write fields.
virtual tmp< volScalarField > Gk() const
Return the momentum drag turbulence generation rate.
Basic sub-grid obstacle drag model. Details supplied by J Puttock 2/7/06.
Definition: basic.H:95
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
basic(const dictionary &PDRProperties, const compressible::RASModel &turbulence, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
void operator=(const basic &)=delete
Disallow default bitwise assignment.
TypeName("basic")
Runtime type information.
U
Definition: pEqn.H:72
A class for managing temporary objects.
Definition: PtrList.H:53
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.