Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
34namespace Foam
35{
36namespace PDRDragModels
37{
40}
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46Foam::PDRDragModels::basic::basic
47(
48 const dictionary& PDRProperties,
49 const compressible::RASModel& turbulence,
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
95Foam::tmp<Foam::volSymmTensorField> Foam::PDRDragModels::basic::Dcu() const
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
119Foam::tmp<Foam::volScalarField> Foam::PDRDragModels::basic::Gk() const
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 =
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
147bool 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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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())
@ NO_REGISTER
Do not request registration (bool: false).
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
virtual bool read()
Inherit read from regIOobject.
Basic sub-grid obstacle drag model. Details supplied by J Puttock 2/7/06.
Definition basic.H:96
virtual ~basic()
Destructor.
virtual tmp< volSymmTensorField > Dcu() const
Return the momentum drag coefficient.
virtual tmp< volScalarField > Gk() const
Return the momentum drag turbulence generation rate.
void writeFields() const
Write fields.
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...
virtual bool read()
Read object.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const Identity< scalar > I
Definition Identity.H:100
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
const volScalarField & betav