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-2021 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 = tmp<volScalarField>::New
60  (
61  IOobject
62  (
63  typeName + ":dragCoeff",
64  mesh_.time().timeName(),
65  mesh_.time(),
68  ),
69  mesh_,
71  );
72  auto& dragCoeff = tdragCoeff.ref();
73 
74  forAll(zoneIDs_, i)
75  {
76  const scalar a = aZone_[i];
77  const scalar N = NZone_[i];
78  const scalar Cd = CdZone_[i];
79 
80  const labelList& zones = zoneIDs_[i];
81 
82  for (label zonei : zones)
83  {
84  const cellZone& cz = mesh_.cellZones()[zonei];
85 
86  for (label celli : cz)
87  {
88  const vector& Uc = U[celli];
89 
90  dragCoeff[celli] = 0.5*Cd*a*N*mag(Uc);
91  }
92  }
93  }
94 
95  dragCoeff.correctBoundaryConditions();
96 
97  return tdragCoeff;
98 }
99 
100 
103 {
104  auto tinertiaCoeff = tmp<volScalarField>::New
105  (
106  IOobject
107  (
108  typeName + ":inertiaCoeff",
109  mesh_.time().timeName(),
110  mesh_.time(),
113  ),
114  mesh_,
116  );
117  auto& inertiaCoeff = tinertiaCoeff.ref();
118 
119  const scalar pi = constant::mathematical::pi;
120 
121  forAll(zoneIDs_, i)
122  {
123  const scalar a = aZone_[i];
124  const scalar N = NZone_[i];
125  const scalar Cm = CmZone_[i];
126 
127  const labelList& zones = zoneIDs_[i];
128 
129  for (label zonei : zones)
130  {
131  const cellZone& cz = mesh_.cellZones()[zonei];
132 
133  for (label celli : cz)
134  {
135  inertiaCoeff[celli] = 0.25*(Cm+1)*pi*a*a*N;
136  }
137  }
138  }
139 
140  inertiaCoeff.correctBoundaryConditions();
141 
142  return tinertiaCoeff;
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
148 Foam::fv::multiphaseMangrovesSource::multiphaseMangrovesSource
149 (
150  const word& name,
151  const word& modelType,
152  const dictionary& dict,
153  const fvMesh& mesh
154 )
155 :
156  fv::option(name, modelType, dict, mesh),
157  aZone_(),
158  NZone_(),
159  CmZone_(),
160  CdZone_()
161 {
163 }
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 (
170  fvMatrix<vector>& eqn,
171  const label fieldi
172 )
173 {
174  const volVectorField& U = eqn.psi();
175 
176  fvMatrix<vector> mangrovesEqn
177  (
178  - fvm::Sp(dragCoeff(U), U)
179  - inertiaCoeff()*fvm::ddt(U)
180  );
182  // Contributions are added to RHS of momentum equation
183  eqn += mangrovesEqn;
184 }
185 
186 
188 (
189  const volScalarField& rho,
190  fvMatrix<vector>& eqn,
191  const label fieldi
192 )
193 {
194  const volVectorField& U = eqn.psi();
195 
196  fvMatrix<vector> mangrovesEqn
197  (
198  - fvm::Sp(rho*dragCoeff(U), U)
199  - rho*inertiaCoeff()*fvm::ddt(U)
200  );
201 
202  // Contributions are added to RHS of momentum equation
203  eqn += mangrovesEqn;
204 }
205 
206 
208 {
209  if (fv::option::read(dict))
210  {
211  if (!coeffs_.readIfPresent("UNames", fieldNames_))
212  {
213  fieldNames_.resize(1);
214  fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
215  }
217 
218  // Create the Mangroves models - 1 per region
219  const dictionary& regionsDict(coeffs_.subDict("regions"));
220  const wordList regionNames(regionsDict.toc());
221  aZone_.setSize(regionNames.size(), 1);
222  NZone_.setSize(regionNames.size(), 1);
223  CmZone_.setSize(regionNames.size(), 1);
224  CdZone_.setSize(regionNames.size(), 1);
225  zoneIDs_.setSize(regionNames.size());
226 
227  forAll(zoneIDs_, i)
228  {
229  const word& regionName = regionNames[i];
230  const dictionary& modelDict = regionsDict.subDict(regionName);
231 
232  const word zoneName(modelDict.get<word>("cellZone"));
233 
234  zoneIDs_[i] = mesh_.cellZones().indices(zoneName);
235  if (zoneIDs_[i].empty())
236  {
238  << "Unable to find cellZone " << zoneName << nl
239  << "Valid cellZones are:" << mesh_.cellZones().names()
240  << exit(FatalError);
241  }
242 
243  modelDict.readEntry("a", aZone_[i]);
244  modelDict.readEntry("N", NZone_[i]);
245  modelDict.readEntry("Cm", CmZone_[i]);
246  modelDict.readEntry("Cd", CdZone_[i]);
247  }
248 
249  return true;
250  }
251 
252  return false;
253 }
254 
255 
256 // ************************************************************************* //
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:598
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
Ignore writing from objectRegistry::writeObject()
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 Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
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...
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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.
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...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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
Nothing to be read.
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:678
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< volScalarField > inertiaCoeff() const
Return the inertia force coefficient.
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