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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 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 \*---------------------------------------------------------------------------*/
27 
28 #include "basic.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace PDRDragModels
36 {
38  addToRunTimeSelectionTable(PDRDragModel, basic, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::PDRDragModels::basic::basic
46 (
47  const dictionary& PDRProperties,
49  const volScalarField& rho,
50  const volVectorField& U,
51  const surfaceScalarField& phi
52 )
53 :
54  PDRDragModel(PDRProperties, turbulence, rho, U, phi),
55  Csu("Csu", dimless, PDRDragModelCoeffs_),
56  Csk("Csk", dimless, PDRDragModelCoeffs_),
57 
58  Aw_
59  (
60  IOobject
61  (
62  "Aw",
63  U_.mesh().facesInstance(),
64  U_.mesh(),
65  IOobject::MUST_READ,
66  IOobject::NO_WRITE
67  ),
68  U_.mesh()
69  ),
70 
71  CR_
72  (
73  IOobject
74  (
75  "CR",
76  U_.mesh().facesInstance(),
77  U_.mesh(),
78  IOobject::MUST_READ,
79  IOobject::NO_WRITE
80  ),
81  U_.mesh()
82  )
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
93 
95 {
96  tmp<volSymmTensorField> tDragDcu
97  (
99  (
100  IOobject
101  (
102  "tDragDcu",
103  U_.mesh().time().constant(),
104  U_.mesh(),
107  ),
108  U_.mesh(),
110  )
111  );
112 
113  volSymmTensorField& DragDcu = tDragDcu.ref();
114 
115  if (on_)
116  {
117  const volScalarField& betav =
118  U_.db().lookupObject<volScalarField>("betav");
119 
120  DragDcu =
121  (0.5*rho_)*CR_*mag(U_) + (Csu*I)*betav*turbulence_.muEff()*sqr(Aw_);
122  }
123 
124  return tDragDcu;
125 }
126 
127 
129 {
130  tmp<volScalarField> tGk
131  (
132  new volScalarField
133  (
134  IOobject
135  (
136  "tGk",
137  U_.mesh().time().constant(),
138  U_.mesh(),
141  ),
142  U_.mesh(),
144  )
145  );
146 
147  volScalarField& Gk = tGk.ref();
148 
149  if (on_)
150  {
151  const volScalarField& betav =
152  U_.db().lookupObject<volScalarField>("betav");
153 
154  const volSymmTensorField& CT =
155  U_.db().lookupObject<volSymmTensorField>("CT");
156 
157  Gk =
158  (0.5*rho_)*mag(U_)*(U_ & CT & U_)
159  + Csk*betav*turbulence_.muEff()*sqr(Aw_)*magSqr(U_);
160  }
161 
162  return tGk;
163 }
164 
165 
166 bool Foam::PDRDragModels::basic::read(const dictionary& PDRProperties)
167 {
168  PDRDragModel::read(PDRProperties);
169 
170  PDRDragModelCoeffs_.readEntry("Csu", Csu.value());
171  PDRDragModelCoeffs_.readEntry("Csk", Csk.value());
172 
173  return true;
174 }
175 
176 
178 {
179  Aw_.write();
180  CR_.write();
181 }
182 
183 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual bool read()
Read object.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
virtual ~basic()
Destructor.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual tmp< volSymmTensorField > Dcu() const
Return the momentum drag coefficient.
Macros for easy insertion into run-time selection tables.
void writeFields() const
Write fields.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual tmp< volScalarField > Gk() const
Return the momentum drag turbulence generation rate.
dynamicFvMesh & mesh
static const Identity< scalar > I
Definition: Identity.H:100
const volScalarField & betav
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:425
const wordList basic
Standard basic field types (label, scalar, vector, tensor, etc)
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("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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133