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-2022 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 (debug && mesh_.time().writeTime())
119  {
120  volTensorField Dout
121  (
122  IOobject
123  (
124  typeName + ":D",
125  mesh_.time().timeName(),
126  mesh_,
129  ),
130  mesh_,
131  dimensionedTensor(dXYZ_.dimensions(), Zero)
132  );
133  volTensorField Fout
134  (
135  IOobject
136  (
137  typeName + ":F",
138  mesh_.time().timeName(),
139  mesh_,
142  ),
143  mesh_,
144  dimensionedTensor(fXYZ_.dimensions(), Zero)
145  );
146 
147 
148  forAll(cellZoneIDs_, zonei)
149  {
150  const labelList& cells = mesh_.cellZones()[cellZoneIDs_[zonei]];
151 
152  if (csys().uniform())
153  {
154  UIndirectList<tensor>(Dout, cells) = D_[zonei].first();
155  UIndirectList<tensor>(Fout, cells) = F_[zonei].first();
156  }
157  else
158  {
159  UIndirectList<tensor>(Dout, cells) = D_[zonei];
160  UIndirectList<tensor>(Fout, cells) = F_[zonei];
161  }
162  }
163 
164  Dout.write();
165  Fout.write();
166  }
167 }
168 
169 
171 (
172  const volVectorField& U,
173  const volScalarField& rho,
174  const volScalarField& mu,
175  vectorField& force
176 ) const
177 {
178  scalarField Udiag(U.size(), Zero);
179  vectorField Usource(U.size(), Zero);
180  const scalarField& V = mesh_.V();
181 
182  apply(Udiag, Usource, V, rho, mu, U);
183 
184  force = Udiag*U - Usource;
185 }
186 
187 
189 (
191 ) const
192 {
193  const volVectorField& U = UEqn.psi();
194  const scalarField& V = mesh_.V();
195  scalarField& Udiag = UEqn.diag();
196  vectorField& Usource = UEqn.source();
197 
198  word rhoName(IOobject::groupName(rhoName_, U.group()));
199  word muName(IOobject::groupName(muName_, U.group()));
200  word nuName(IOobject::groupName(nuName_, U.group()));
201 
202  if (UEqn.dimensions() == dimForce)
203  {
204  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
205 
206  if (mesh_.foundObject<volScalarField>(muName))
207  {
208  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
209 
210  apply(Udiag, Usource, V, rho, mu, U);
211  }
212  else
213  {
214  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
215 
216  apply(Udiag, Usource, V, rho, rho*nu, U);
217  }
218  }
219  else
220  {
221  if (mesh_.foundObject<volScalarField>(nuName))
222  {
223  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
224 
225  apply(Udiag, Usource, V, geometricOneField(), nu, U);
226  }
227  else
228  {
229  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
230  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
231 
232  apply(Udiag, Usource, V, geometricOneField(), mu/rho, U);
233  }
234  }
235 }
236 
237 
239 (
241  const volScalarField& rho,
242  const volScalarField& mu
243 ) const
244 {
245  const vectorField& U = UEqn.psi();
246  const scalarField& V = mesh_.V();
247  scalarField& Udiag = UEqn.diag();
248  vectorField& Usource = UEqn.source();
249 
250  apply(Udiag, Usource, V, rho, mu, U);
251 }
252 
253 
255 (
256  const fvVectorMatrix& UEqn,
257  volTensorField& AU
258 ) const
259 {
260  const volVectorField& U = UEqn.psi();
261 
262  word rhoName(IOobject::groupName(rhoName_, U.group()));
263  word muName(IOobject::groupName(muName_, U.group()));
264  word nuName(IOobject::groupName(nuName_, U.group()));
265 
266  if (UEqn.dimensions() == dimForce)
267  {
268  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
269  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
270 
271  apply(AU, rho, mu, U);
272  }
273  else
274  {
275  if (mesh_.foundObject<volScalarField>(nuName))
276  {
277  const auto& nu = mesh_.lookupObject<volScalarField>(nuName);
278 
279  apply(AU, geometricOneField(), nu, U);
280  }
281  else
282  {
283  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName);
284  const auto& mu = mesh_.lookupObject<volScalarField>(muName);
286  apply(AU, geometricOneField(), mu/rho, U);
287  }
288  }
289 }
290 
291 
293 {
294  dict_.writeEntry(name_, os);
295 
296  return true;
297 }
298 
299 
300 // ************************************************************************* //
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:120
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:89
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)
Ignore writing from objectRegistry::writeObject()
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:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
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 INVALID.
Definition: exprTraits.C:52
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:79
UIndirectList< point > pointUIndList
A UIndirectList of points.
Definition: pointIndList.H:47
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
Nothing to be read.
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
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:157