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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "basic.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace PDRDragModels
37 {
39  addToRunTimeSelectionTable(PDRDragModel, basic, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 Foam::PDRDragModels::basic::basic
47 (
48  const dictionary& PDRProperties,
50  const volScalarField& rho,
51  const volVectorField& U,
52  const surfaceScalarField& phi
53 )
54 :
55  PDRDragModel(PDRProperties, turbulence, rho, U, phi),
56  Csu("Csu", dimless, PDRDragModelCoeffs_),
57  Csk("Csk", dimless, PDRDragModelCoeffs_),
58 
59  Aw_
60  (
61  IOobject
62  (
63  "Aw",
64  U_.mesh().facesInstance(),
65  U_.mesh(),
66  IOobject::MUST_READ,
67  IOobject::NO_WRITE
68  ),
69  U_.mesh()
70  ),
71 
72  CR_
73  (
74  IOobject
75  (
76  "CR",
77  U_.mesh().facesInstance(),
78  U_.mesh(),
79  IOobject::MUST_READ,
80  IOobject::NO_WRITE
81  ),
82  U_.mesh()
83  )
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
90 {}
91 
92 
93 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
94 
96 {
97  auto tDragDcu = volSymmTensorField::New
98  (
99  "tDragDcu",
101  U_.mesh(),
103  );
104  auto& DragDcu = tDragDcu.ref();
105 
106  if (on_)
107  {
108  const volScalarField& betav =
109  U_.db().lookupObject<volScalarField>("betav");
110 
111  DragDcu =
112  (0.5*rho_)*CR_*mag(U_) + (Csu*I)*betav*turbulence_.muEff()*sqr(Aw_);
113  }
114 
115  return tDragDcu;
116 }
117 
118 
120 {
121  auto tGk = volScalarField::New
122  (
123  "tGk",
125  U_.mesh(),
127  );
128  auto& Gk = tGk.ref();
129 
130  if (on_)
131  {
132  const volScalarField& betav =
133  U_.db().lookupObject<volScalarField>("betav");
134 
135  const volSymmTensorField& CT =
136  U_.db().lookupObject<volSymmTensorField>("CT");
137 
138  Gk =
139  (0.5*rho_)*mag(U_)*(U_ & CT & U_)
140  + Csk*betav*turbulence_.muEff()*sqr(Aw_)*magSqr(U_);
141  }
142 
143  return tGk;
144 }
145 
146 
147 bool Foam::PDRDragModels::basic::read(const dictionary& PDRProperties)
148 {
149  PDRDragModel::read(PDRProperties);
150 
151  PDRDragModelCoeffs_.readEntry("Csu", Csu.value());
152  PDRDragModelCoeffs_.readEntry("Csk", Csk.value());
153 
154  return true;
155 }
156 
157 
159 {
160  Aw_.write();
161  CR_.write();
162 }
163 
164 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:84
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)
const dimensionSet dimless
Dimensionless.
virtual ~basic()
Destructor.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
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:72
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:450
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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
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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
Do not request registration (bool: false)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127