DarcyForchheimer.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2018-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 "DarcyForchheimer.H"
31 #include "geometricOneField.H"
32 #include "fvMatrices.H"
33 #include "pointIndList.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace porosityModels
40  {
41  defineTypeNameAndDebug(DarcyForchheimer, 0);
42  addToRunTimeSelectionTable(porosityModel, DarcyForchheimer, mesh);
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::porosityModels::DarcyForchheimer::DarcyForchheimer
50 (
51  const word& name,
52  const word& modelType,
53  const fvMesh& mesh,
54  const dictionary& dict,
55  const wordRe& cellZoneName
56 )
57 :
58  porosityModel(name, modelType, mesh, dict, cellZoneName),
59  dXYZ_("d", dimless/sqr(dimLength), coeffs_),
60  fXYZ_("f", dimless/dimLength, coeffs_),
61  D_(cellZoneIDs_.size()),
62  F_(cellZoneIDs_.size()),
63  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
64  muName_(coeffs_.getOrDefault<word>("mu", "thermo:mu")),
65  nuName_(coeffs_.getOrDefault<word>("nu", "nu"))
66 {
69 
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  // The Darcy coefficient as a tensor
79  tensor darcyCoeff(Zero);
80  darcyCoeff.xx() = dXYZ_.value().x();
81  darcyCoeff.yy() = dXYZ_.value().y();
82  darcyCoeff.zz() = dXYZ_.value().z();
83 
84  // The Forchheimer coefficient as a tensor
85  // - the leading 0.5 is from 1/2*rho
86  tensor forchCoeff(Zero);
87  forchCoeff.xx() = 0.5*fXYZ_.value().x();
88  forchCoeff.yy() = 0.5*fXYZ_.value().y();
89  forchCoeff.zz() = 0.5*fXYZ_.value().z();
90 
91  if (csys().uniform())
92  {
93  forAll(cellZoneIDs_, zonei)
94  {
95  D_[zonei].resize(1);
96  F_[zonei].resize(1);
97 
98  D_[zonei] = csys().transform(darcyCoeff);
99  F_[zonei] = csys().transform(forchCoeff);
100  }
101  }
102  else
103  {
104  forAll(cellZoneIDs_, zonei)
105  {
106  const pointUIndList cc
107  (
108  mesh_.cellCentres(),
109  mesh_.cellZones()[cellZoneIDs_[zonei]]
110  );
111 
112  D_[zonei] = csys().transform(cc, darcyCoeff);
113  F_[zonei] = csys().transform(cc, forchCoeff);
114  }
115  }
116 
117 
118  if
119  (
120  debug
121  && (
122  mesh_.time().writeTime()
123  || mesh_.time().timeIndex() == mesh_.time().startTimeIndex()
124  )
125  )
126  {
127  volTensorField Dout
128  (
129  mesh_.newIOobject(IOobject::scopedName(typeName, "D")),
130  mesh_,
131  dimensionedTensor(dXYZ_.dimensions(), Zero)
132  );
133 
134  volTensorField Fout
135  (
136  mesh_.newIOobject(IOobject::scopedName(typeName, "F")),
137  mesh_,
138  dimensionedTensor(fXYZ_.dimensions(), Zero)
139  );
140 
141  forAll(cellZoneIDs_, zonei)
142  {
143  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zonei]];
144 
145  if (csys().uniform())
146  {
147  UIndirectList<tensor>(Dout, cells) = D_[zonei].first();
148  UIndirectList<tensor>(Fout, cells) = F_[zonei].first();
149  }
150  else
151  {
152  UIndirectList<tensor>(Dout, cells) = D_[zonei];
153  UIndirectList<tensor>(Fout, cells) = F_[zonei];
154  }
155  }
156 
157  Dout.write();
158  Fout.write();
159  }
160 }
161 
162 
164 (
165  const volVectorField& U,
166  const volScalarField& rho,
167  const volScalarField& mu,
168  vectorField& force
169 ) const
170 {
171  scalarField Udiag(U.size(), Zero);
172  vectorField Usource(U.size(), Zero);
173  const scalarField& V = mesh_.V();
174 
175  apply(Udiag, Usource, V, rho, mu, U);
176 
177  force = Udiag*U - Usource;
178 }
179 
180 
182 (
184 ) const
185 {
186  const volVectorField& U = UEqn.psi();
187  const scalarField& V = mesh_.V();
188  scalarField& Udiag = UEqn.diag();
189  vectorField& Usource = UEqn.source();
190 
191  word rhoName(IOobject::groupName(rhoName_, U.group()));
192  word muName(IOobject::groupName(muName_, U.group()));
193  word nuName(IOobject::groupName(nuName_, U.group()));
194 
195  if (UEqn.dimensions() == dimForce)
196  {
197  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
198 
199  if (mesh_.foundObject<volScalarField>(muName))
200  {
201  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
202 
203  apply(Udiag, Usource, V, rho, mu, U);
204  }
205  else
206  {
207  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
208 
209  apply(Udiag, Usource, V, rho, rho*nu, U);
210  }
211  }
212  else
213  {
214  if (mesh_.foundObject<volScalarField>(nuName))
215  {
216  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
217 
218  apply(Udiag, Usource, V, geometricOneField(), nu, U);
219  }
220  else
221  {
222  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
223  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
224 
225  apply(Udiag, Usource, V, geometricOneField(), mu/rho, U);
226  }
227  }
228 }
229 
230 
232 (
234  const volScalarField& rho,
235  const volScalarField& mu
236 ) const
237 {
238  const vectorField& U = UEqn.psi();
239  const scalarField& V = mesh_.V();
240  scalarField& Udiag = UEqn.diag();
241  vectorField& Usource = UEqn.source();
242 
243  apply(Udiag, Usource, V, rho, mu, U);
244 }
245 
246 
248 (
249  const fvVectorMatrix& UEqn,
250  volTensorField& AU
251 ) const
252 {
253  const volVectorField& U = UEqn.psi();
254 
255  word rhoName(IOobject::groupName(rhoName_, U.group()));
256  word muName(IOobject::groupName(muName_, U.group()));
257  word nuName(IOobject::groupName(nuName_, U.group()));
258 
259  if (UEqn.dimensions() == dimForce)
260  {
261  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
262  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
263 
264  apply(AU, rho, mu, U);
265  }
266  else
267  {
268  if (mesh_.foundObject<volScalarField>(nuName))
269  {
270  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
271 
272  apply(AU, geometricOneField(), nu, U);
273  }
274  else
275  {
276  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
277  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
279  apply(AU, geometricOneField(), mu/rho, U);
280  }
281  }
282 }
283 
284 
286 {
287  dict_.writeEntry(name_, os);
288 
289  return true;
290 }
291 
292 
293 // ************************************************************************* //
dictionary dict
bool writeData(Ostream &os) const
Write.
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:88
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:36
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
const dimensionSet dimForce
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
UIndirectList< point > pointUIndList
UIndirectList of point.
Definition: pointIndList.H:36
Tensor of scalars, i.e. Tensor<scalar>.
volScalarField & nu
Top level model for porosity models.
Definition: porosityModel.H:53
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127