basic.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 {
35  defineTypeNameAndDebug(basic, 0);
36  addToRunTimeSelectionTable(PDRDragModel, basic, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::PDRDragModels::basic::basic
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 {
94  tmp<volSymmTensorField> tDragDcu
95  (
97  (
98  IOobject
99  (
100  "tDragDcu",
101  U_.mesh().time().constant(),
102  U_.mesh(),
105  ),
106  U_.mesh(),
108  (
109  "zero",
111  Zero
112  )
113  )
114  );
115 
116  volSymmTensorField& DragDcu = tDragDcu.ref();
117 
118  if (on_)
119  {
120  const volScalarField& betav =
121  U_.db().lookupObject<volScalarField>("betav");
122 
123  DragDcu =
124  (0.5*rho_)*CR_*mag(U_) + (Csu*I)*betav*turbulence_.muEff()*sqr(Aw_);
125  }
126 
127  return tDragDcu;
128 }
129 
130 
132 {
133  tmp<volScalarField> tGk
134  (
135  new volScalarField
136  (
137  IOobject
138  (
139  "tGk",
140  U_.mesh().time().constant(),
141  U_.mesh(),
144  ),
145  U_.mesh(),
146  dimensionedScalar("zero", dimMass/dimLength/pow(dimTime, 3), 0.0)
147  )
148  );
149 
150  volScalarField& Gk = tGk.ref();
151 
152  if (on_)
153  {
154  const volScalarField& betav =
155  U_.db().lookupObject<volScalarField>("betav");
156 
157  const volSymmTensorField& CT =
158  U_.db().lookupObject<volSymmTensorField>("CT");
159 
160  Gk =
161  (0.5*rho_)*mag(U_)*(U_ & CT & U_)
162  + Csk*betav*turbulence_.muEff()*sqr(Aw_)*magSqr(U_);
163  }
164 
165  return tGk;
166 }
167 
168 
169 bool Foam::PDRDragModels::basic::read(const dictionary& PDRProperties)
170 {
171  PDRDragModel::read(PDRProperties);
172 
173  PDRDragModelCoeffs_.lookup("Csu") >> Csu.value();
174  PDRDragModelCoeffs_.lookup("Csk") >> Csk.value();
175 
176  return true;
177 }
178 
179 
181 {
182  Aw_.write();
183  CR_.write();
184 }
185 
186 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
U
Definition: pEqn.H:83
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
virtual bool read()
Read object.
void writeFields() const
Write fields.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual ~basic()
Destructor.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual tmp< volScalarField > Gk() const
Return the momentum drag turbulence generation rate.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dynamicFvMesh & mesh
static const Identity< scalar > I
Definition: Identity.H:93
static const zero Zero
Definition: zero.H:91
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual tmp< volSymmTensorField > Dcu() const
Return the momentum drag coefficient.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:54
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.