basic.C
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-2023 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 \*---------------------------------------------------------------------------*/
25 
26 #include "basic.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace PDRDragModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& PDRProperties,
47  const volScalarField& rho,
48  const volVectorField& U,
49  const surfaceScalarField& phi
50 )
51 :
52  PDRDragModel(PDRProperties, turbulence, rho, U, phi),
53  Csu("Csu", dimless, PDRDragModelCoeffs_),
54  Csk("Csk", dimless, PDRDragModelCoeffs_),
55 
56  Aw_
57  (
58  IOobject
59  (
60  "Aw",
61  U_.mesh().facesInstance(),
62  U_.mesh(),
63  IOobject::MUST_READ,
64  IOobject::NO_WRITE
65  ),
66  U_.mesh()
67  ),
68 
69  CR_
70  (
71  IOobject
72  (
73  "CR",
74  U_.mesh().facesInstance(),
75  U_.mesh(),
76  IOobject::MUST_READ,
77  IOobject::NO_WRITE
78  ),
79  U_.mesh()
80  )
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
85 
87 {}
88 
89 
90 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
91 
93 {
95  (
97  (
98  "tDragDcu",
99  U_.mesh(),
101  (
103  Zero
104  )
105  )
106  );
107 
108  volSymmTensorField& DragDcu = tDragDcu.ref();
109 
110  if (on_)
111  {
112  const volScalarField& betav =
113  U_.db().lookupObject<volScalarField>("betav");
114 
115  DragDcu =
116  rho_*(0.5*CR_*mag(U_) + (Csu*I)*betav*turbulence_.nuEff()*sqr(Aw_));
117  }
118 
119  return tDragDcu;
120 }
121 
122 
124 {
126  (
128  (
129  "tGk",
130  U_.mesh(),
132  )
133  );
134 
135  volScalarField& Gk = tGk.ref();
136 
137  if (on_)
138  {
139  const volScalarField& betav =
140  U_.db().lookupObject<volScalarField>("betav");
141 
142  const volSymmTensorField& CT =
143  U_.db().lookupObject<volSymmTensorField>("CT");
144 
145  Gk =
146  rho_
147  *(
148  0.5*mag(U_)*(U_ & CT & U_)
149  + Csk*betav*turbulence_.nuEff()*sqr(Aw_)*magSqr(U_)
150  );
151  }
152 
153  return tGk;
154 }
155 
156 
158 {
159  PDRDragModel::read(PDRProperties);
160 
161  PDRDragModelCoeffs_.lookup("Csu") >> Csu.value();
162  PDRDragModelCoeffs_.lookup("Csk") >> Csk.value();
163 
164  return true;
165 }
166 
167 
169 {
170  Aw_.write();
171  CR_.write();
172 }
173 
174 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Base-class for sub-grid obstacle drag models. The available drag model is at basic....
Definition: PDRDragModel.H:56
virtual bool read()
Inherit read from regIOobject.
Basic sub-grid obstacle drag model. Details supplied by J Puttock 2/7/06.
Definition: basic.H:98
virtual tmp< volSymmTensorField > Dcu() const
Return the momentum drag coefficient.
Definition: basic.C:92
virtual ~basic()
Destructor.
Definition: basic.C:86
basic(const dictionary &PDRProperties, const compressible::RASModel &turbulence, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Definition: basic.C:44
void writeFields() const
Write fields.
Definition: basic.C:168
virtual tmp< volScalarField > Gk() const
Return the momentum drag turbulence generation rate.
Definition: basic.C:123
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Generic dimensioned Type class.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(PDRDragModel, basic, dictionary)
defineTypeNameAndDebug(basic, 0)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
const dimensionSet dimLength
static const Identity< scalar > I
Definition: Identity.H:93
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMass
dimensioned< scalar > magSqr(const dimensioned< Type > &)
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))