multiphaseMangrovesSource.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 "fvmDdt.H"
34 #include "fvmSup.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace fv
42 {
43  defineTypeNameAndDebug(multiphaseMangrovesSource, 0);
45  (
46  option,
47  multiphaseMangrovesSource,
48  dictionary
49  );
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55 
58 {
59  auto tdragCoeff = volScalarField::New
60  (
61  IOobject::scopedName(typeName, "dragCoeff"),
63  mesh_,
65  );
66  auto& dragCoeff = tdragCoeff.ref();
67 
68  forAll(zoneIDs_, i)
69  {
70  const scalar a = aZone_[i];
71  const scalar N = NZone_[i];
72  const scalar Cd = CdZone_[i];
73 
74  const labelList& zones = zoneIDs_[i];
75 
76  for (label zonei : zones)
77  {
78  const cellZone& cz = mesh_.cellZones()[zonei];
79 
80  for (label celli : cz)
81  {
82  const vector& Uc = U[celli];
83 
84  dragCoeff[celli] = 0.5*Cd*a*N*mag(Uc);
85  }
86  }
87  }
88 
89  dragCoeff.correctBoundaryConditions();
90 
91  return tdragCoeff;
92 }
93 
94 
97 {
98  auto tinertiaCoeff = volScalarField::New
99  (
100  IOobject::scopedName(typeName, "inertiaCoeff"),
102  mesh_,
104  );
105  auto& inertiaCoeff = tinertiaCoeff.ref();
106 
107  const scalar pi = constant::mathematical::pi;
108 
109  forAll(zoneIDs_, i)
110  {
111  const scalar a = aZone_[i];
112  const scalar N = NZone_[i];
113  const scalar Cm = CmZone_[i];
114 
115  const labelList& zones = zoneIDs_[i];
116 
117  for (label zonei : zones)
118  {
119  const cellZone& cz = mesh_.cellZones()[zonei];
120 
121  for (label celli : cz)
122  {
123  inertiaCoeff[celli] = 0.25*(Cm+1)*pi*a*a*N;
124  }
125  }
126  }
127 
128  inertiaCoeff.correctBoundaryConditions();
129 
130  return tinertiaCoeff;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
136 Foam::fv::multiphaseMangrovesSource::multiphaseMangrovesSource
137 (
138  const word& name,
139  const word& modelType,
140  const dictionary& dict,
141  const fvMesh& mesh
142 )
143 :
144  fv::option(name, modelType, dict, mesh),
145  aZone_(),
146  NZone_(),
147  CmZone_(),
148  CdZone_()
149 {
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 (
158  fvMatrix<vector>& eqn,
159  const label fieldi
160 )
161 {
162  const volVectorField& U = eqn.psi();
163 
164  fvMatrix<vector> mangrovesEqn
165  (
166  - fvm::Sp(dragCoeff(U), U)
167  - inertiaCoeff()*fvm::ddt(U)
168  );
170  // Contributions are added to RHS of momentum equation
171  eqn += mangrovesEqn;
172 }
173 
174 
176 (
177  const volScalarField& rho,
178  fvMatrix<vector>& eqn,
179  const label fieldi
180 )
181 {
182  const volVectorField& U = eqn.psi();
183 
184  fvMatrix<vector> mangrovesEqn
185  (
186  - fvm::Sp(rho*dragCoeff(U), U)
187  - rho*inertiaCoeff()*fvm::ddt(U)
188  );
189 
190  // Contributions are added to RHS of momentum equation
191  eqn += mangrovesEqn;
192 }
193 
194 
196 {
197  if (fv::option::read(dict))
198  {
199  if (!coeffs_.readIfPresent("UNames", fieldNames_))
200  {
201  fieldNames_.resize(1);
202  fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
203  }
205 
206  // Create the Mangroves models - 1 per region
207  const dictionary& regionsDict(coeffs_.subDict("regions"));
208  const wordList regionNames(regionsDict.toc());
209  aZone_.setSize(regionNames.size(), 1);
210  NZone_.setSize(regionNames.size(), 1);
211  CmZone_.setSize(regionNames.size(), 1);
212  CdZone_.setSize(regionNames.size(), 1);
213  zoneIDs_.setSize(regionNames.size());
214 
215  forAll(zoneIDs_, i)
216  {
217  const word& regionName = regionNames[i];
218  const dictionary& modelDict = regionsDict.subDict(regionName);
219 
220  const word zoneName(modelDict.get<word>("cellZone"));
221 
222  zoneIDs_[i] = mesh_.cellZones().indices(zoneName);
223  if (zoneIDs_[i].empty())
224  {
226  << "Unable to find cellZone " << zoneName << nl
227  << "Valid cellZones are:" << mesh_.cellZones().names()
228  << exit(FatalError);
229  }
230 
231  modelDict.readEntry("a", aZone_[i]);
232  modelDict.readEntry("N", NZone_[i]);
233  modelDict.readEntry("Cm", CmZone_[i]);
234  modelDict.readEntry("Cd", CdZone_[i]);
235  }
236 
237  return true;
238  }
239 
240  return false;
241 }
242 
243 
244 // ************************************************************************* //
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
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
wordList regionNames
const dimensionSet dimless
Dimensionless.
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
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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...
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)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
scalarList NZone_
Number of plants per unit of area.
constexpr scalar pi(M_PI)
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
tmp< volScalarField > dragCoeff(const volVectorField &U) const
Return the drag force coefficient.
labelListList zoneIDs_
Porosity coefficient.
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)
scalarList aZone_
Width of the vegetation element.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A subset of mesh cells.
Definition: cellZone.H:58
virtual bool read(const dictionary &dict)
Read dictionary.
const Vector< label > N(dict.get< Vector< label >>("N"))
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add implicit contribution to momentum equation.
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
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< volScalarField > inertiaCoeff() const
Return the inertia force coefficient.
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