Loading...
Searching...
No Matches
multiphaseMangrovesTurbulenceModel.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) 2017 IH-Cantabria
9 Copyright (C) 2017-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
31#include "fvMesh.H"
32#include "fvMatrices.H"
33#include "fvmSup.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace fv
41{
44 (
45 option,
48 );
50}
51
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
57(
58 const volVectorField& U
59) const
60{
61 auto tkCoeff = volScalarField::New
62 (
65 mesh_,
67 );
68 auto& kCoeff = tkCoeff.ref();
69
71 {
72 const scalar a = aZone_[i];
73 const scalar N = NZone_[i];
74 const scalar Ckp = CkpZone_[i];
75 const scalar Cd = CdZone_[i];
76
77 for (label zonei : zoneIDs_[i])
78 {
79 const cellZone& cz = mesh_.cellZones()[zonei];
80
81 for (label celli : cz)
82 {
83 kCoeff[celli] = Ckp*Cd*a*N*mag(U[celli]);
84 }
85 }
86 }
87
88 kCoeff.correctBoundaryConditions();
89
90 return tkCoeff;
91}
92
93
96(
97 const volVectorField& U
98) const
99{
100 auto tepsilonCoeff = volScalarField::New
101 (
102 IOobject::scopedName(typeName, "epsilonCoeff"),
103 mesh_,
105 );
106 auto& epsilonCoeff = tepsilonCoeff.ref();
107
108 forAll(zoneIDs_, i)
109 {
110 const scalar a = aZone_[i];
111 const scalar N = NZone_[i];
112 const scalar Cep = CepZone_[i];
113 const scalar Cd = CdZone_[i];
114
115 for (label zonei : zoneIDs_[i])
116 {
117 const cellZone& cz = mesh_.cellZones()[zonei];
118
119 for (label celli : cz)
120 {
121 epsilonCoeff[celli] = Cep*Cd*a*N*mag(U[celli]);
122 }
123 }
124 }
125
126 epsilonCoeff.correctBoundaryConditions();
128 return tepsilonCoeff;
129}
130
131
132// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133
134Foam::fv::multiphaseMangrovesTurbulenceModel::multiphaseMangrovesTurbulenceModel
135(
136 const word& name,
137 const word& modelType,
138 const dictionary& dict,
139 const fvMesh& mesh
140)
141:
142 fv::option(name, modelType, dict, mesh),
143 aZone_(),
144 NZone_(),
145 CkpZone_(),
146 CepZone_(),
147 CdZone_(),
148 UName_("U"),
149 kName_("k"),
150 epsilonName_("epsilon")
152 read(dict);
153}
154
155
156// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157
159(
160 fvMatrix<scalar>& eqn,
161 const label fieldi
162)
163{
164 const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
165
166 if (eqn.psi().name() == epsilonName_)
167 {
168 fvMatrix<scalar> epsilonEqn
169 (
170 - fvm::Sp(epsilonCoeff(U), eqn.psi())
171 );
172 eqn += epsilonEqn;
173 }
174 else if (eqn.psi().name() == kName_)
175 {
177 (
178 - fvm::Sp(kCoeff(U), eqn.psi())
179 );
180 eqn += kEqn;
181 }
182}
183
184
186(
187 const volScalarField& rho,
188 fvMatrix<scalar>& eqn,
189 const label fieldi
190)
191{
192 const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
193
194 if (eqn.psi().name() == epsilonName_)
195 {
196 fvMatrix<scalar> epsilonEqn
197 (
198 - fvm::Sp(rho*epsilonCoeff(U), eqn.psi())
199 );
200 eqn += epsilonEqn;
201 }
202 else if (eqn.psi().name() == kName_)
203 {
204 fvMatrix<scalar> kEqn
205 (
206 - fvm::Sp(rho*kCoeff(U), eqn.psi())
207 );
208 eqn += kEqn;
209 }
210}
211
212
214{
216 {
217 if (!coeffs_.readIfPresent("epsilonNames", fieldNames_))
218 {
219 if (coeffs_.found("epsilon"))
220 {
221 fieldNames_.resize(1);
222 coeffs_.readEntry("epsilon", fieldNames_.first());
223 }
224 else
225 {
226 fieldNames_.resize(2);
227 fieldNames_[0] = "epsilon";
228 fieldNames_[1] = "k";
229 }
230 }
232
233 // Create the Mangroves models - 1 per region
234 const dictionary& regionsDict(coeffs_.subDict("regions"));
235 const wordList regionNames(regionsDict.toc());
236 aZone_.setSize(regionNames.size(), 1);
237 NZone_.setSize(regionNames.size(), 1);
238 CkpZone_.setSize(regionNames.size(), 1);
239 CepZone_.setSize(regionNames.size(), 1);
240 CdZone_.setSize(regionNames.size(), 1);
241 zoneIDs_.setSize(regionNames.size());
242
243 forAll(zoneIDs_, i)
244 {
245 const word& regionName = regionNames[i];
246 const dictionary& modelDict = regionsDict.subDict(regionName);
247
248 const word zoneName(modelDict.get<word>("cellZone"));
249
250 zoneIDs_[i] = mesh_.cellZones().indices(zoneName);
251 if (zoneIDs_[i].empty())
252 {
254 << "Unable to find cellZone " << zoneName << nl
255 << "Valid cellZones are:" << mesh_.cellZones().names()
256 << exit(FatalError);
257 }
258
259 modelDict.readEntry("a", aZone_[i]);
260 modelDict.readEntry("N", NZone_[i]);
261 modelDict.readEntry("Ckp", CkpZone_[i]);
262 modelDict.readEntry("Cep", CepZone_[i]);
263 modelDict.readEntry("Cd", CdZone_[i]);
264 }
265
266 return true;
267 }
268
269 return false;
270}
271
272
273// ************************************************************************* //
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< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A subset of mesh cells.
Definition cellZone.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
wordList toc() const
Return the table of contents.
Definition dictionary.C:587
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
tmp< volScalarField > kCoeff(const volVectorField &U) const
Return the k coefficient.
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add implicit contribution to momentum equation.
virtual bool read(const dictionary &dict)
Read dictionary.
scalarList NZone_
Number of plants per unit of area.
tmp< volScalarField > epsilonCoeff(const volVectorField &U) const
Return the epsilon coefficient.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
option(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition fvOption.C:51
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition fvOptionIO.C:48
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the finiteVolume matrix for implicit and explicit sources.
wordList regionNames
Namespace for finite-volume.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const Vector< label > N(dict.get< Vector< label > >("N"))