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 -------------------------------------------------------------------------------
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 
30 #include "mathematicalConstants.H"
31 #include "fvMesh.H"
32 #include "fvMatrices.H"
33 #include "fvmSup.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(multiphaseMangrovesTurbulenceModel, 0);
44  (
45  option,
46  multiphaseMangrovesTurbulenceModel,
47  dictionary
48  );
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
57 (
58  const volVectorField& U
59 ) const
60 {
61  auto tkCoeff = volScalarField::New
62  (
63  IOobject::scopedName(typeName, "kCoeff"),
65  mesh_,
67  );
68  auto& kCoeff = tkCoeff.ref();
69 
70  forAll(zoneIDs_, i)
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();
127 
128  return tepsilonCoeff;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
134 Foam::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")
151 {
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  {
176  fvMatrix<scalar> kEqn
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 {
215  if (fv::option::read(dict))
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 // ************************************************************************* //
virtual bool read(const dictionary &dict)
Read dictionary.
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
scalarList NZone_
Number of plants per unit of area.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
wordList regionNames
const dimensionSet dimless
Dimensionless.
tmp< volScalarField > kCoeff(const volVectorField &U) const
Return the k coefficient.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add implicit contribution to momentum equation.
A subset of mesh cells.
Definition: cellZone.H:58
const Vector< label > N(dict.get< Vector< label >>("N"))
List< word > wordList
List of word.
Definition: fileName.H:59
U
Definition: pEqn.H:72
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
tmp< volScalarField > epsilonCoeff(const volVectorField &U) const
Return the epsilon coefficient.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
scalarList aZone_
Width of the vegetation element.
Do not request registration (bool: false)
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127