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