DarcyForchheimer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "DarcyForchheimer.H"
28 #include "geometricOneField.H"
29 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace porosityModels
36  {
39  }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& name,
48  const word& modelType,
49  const fvMesh& mesh,
50  const dictionary& dict,
51  const word& cellZoneName
52 )
53 :
54  porosityModel(name, modelType, mesh, dict, cellZoneName),
55  dXYZ_("d", dimless/sqr(dimLength), coeffs_),
56  fXYZ_("f", dimless/dimLength, coeffs_),
57  rhoName_(coeffs_.lookupOrDefault<word>("rho", "rho")),
58  muName_(coeffs_.lookupOrDefault<word>("mu", "mu")),
59  nuName_(coeffs_.lookupOrDefault<word>("nu", "nu"))
60 {
63 
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  if (coordSys_.R().uniform())
79  {
80  D_.setSize(1);
81  F_.setSize(1);
82 
83  D_[0] = Zero;
84  D_[0].xx() = dXYZ_.value().x();
85  D_[0].yy() = dXYZ_.value().y();
86  D_[0].zz() = dXYZ_.value().z();
87 
88  D_[0] = coordSys_.R().transform(Zero, D_[0]);
89 
90  // leading 0.5 is from 1/2*rho
91  F_[0] = Zero;
92  F_[0].xx() = 0.5*fXYZ_.value().x();
93  F_[0].yy() = 0.5*fXYZ_.value().y();
94  F_[0].zz() = 0.5*fXYZ_.value().z();
95 
96  F_[0] = coordSys_.R().transform(Zero, F_[0]);
97  }
98  else
99  {
100  const labelList& cells = mesh_.cellZones()[zoneName_];
101 
102  D_.setSize(cells.size());
103  F_.setSize(cells.size());
104 
105  forAll(cells, i)
106  {
107  D_[i] = Zero;
108  D_[i].xx() = dXYZ_.value().x();
109  D_[i].yy() = dXYZ_.value().y();
110  D_[i].zz() = dXYZ_.value().z();
111 
112  // leading 0.5 is from 1/2*rho
113  F_[i] = Zero;
114  F_[i].xx() = 0.5*fXYZ_.value().x();
115  F_[i].yy() = 0.5*fXYZ_.value().y();
116  F_[i].zz() = 0.5*fXYZ_.value().z();
117  }
118 
119  const coordinateRotation& R = coordSys_.R
120  (
121  UIndirectList<vector>(mesh_.C(), cells)()
122  );
123 
124  D_ = R.transform(D_);
125  F_ = R.transform(F_);
126  }
127 
128  if (debug && (mesh_.time().writeTime() || mesh_.time().timeIndex() == 0))
129  {
130  volTensorField Dout
131  (
132  IOobject
133  (
134  typedName("D"),
135  mesh_.time().name(),
136  mesh_,
139  ),
140  mesh_,
141  dimensionedTensor(dXYZ_.dimensions(), Zero)
142  );
143  volTensorField Fout
144  (
145  IOobject
146  (
147  typedName("F"),
148  mesh_.time().name(),
149  mesh_,
152  ),
153  mesh_,
154  dimensionedTensor(fXYZ_.dimensions(), Zero)
155  );
156 
157  UIndirectList<tensor>(Dout, mesh_.cellZones()[zoneName_]) = D_;
158  UIndirectList<tensor>(Fout, mesh_.cellZones()[zoneName_]) = F_;
159 
160  Dout.write();
161  Fout.write();
162  }
163 }
164 
165 
167 (
168  const volVectorField& U,
169  const volScalarField& rho,
170  const volScalarField& mu,
171  vectorField& force
172 ) const
173 {
174  scalarField Udiag(U.size(), 0.0);
175  vectorField Usource(U.size(), Zero);
176  const scalarField& V = mesh_.V();
177 
178  apply(Udiag, Usource, V, rho, mu, U);
179 
180  force = Udiag*U - Usource;
181 }
182 
183 
185 (
187 ) const
188 {
189  const volVectorField& U = UEqn.psi();
190  const scalarField& V = mesh_.V();
191  scalarField& Udiag = UEqn.diag();
192  vectorField& Usource = UEqn.source();
193 
194  word rhoName(IOobject::groupName(rhoName_, U.group()));
195  word muName(IOobject::groupName(muName_, U.group()));
196  word nuName(IOobject::groupName(nuName_, U.group()));
197 
198  if (UEqn.dimensions() == dimForce)
199  {
200  const volScalarField& rho = mesh_.lookupObject<volScalarField>(rhoName);
201 
202  if (mesh_.foundObject<volScalarField>(muName))
203  {
204  const volScalarField& mu =
205  mesh_.lookupObject<volScalarField>(muName);
206 
207  apply(Udiag, Usource, V, rho, mu, U);
208  }
209  else
210  {
211  const volScalarField& nu =
212  mesh_.lookupObject<volScalarField>(nuName);
213 
214  apply(Udiag, Usource, V, rho, rho*nu, U);
215  }
216  }
217  else
218  {
219  if (mesh_.foundObject<volScalarField>(nuName))
220  {
221  const volScalarField& nu =
222  mesh_.lookupObject<volScalarField>(nuName);
223 
224  apply(Udiag, Usource, V, geometricOneField(), nu, U);
225  }
226  else
227  {
228  const volScalarField& rho =
229  mesh_.lookupObject<volScalarField>(rhoName);
230  const volScalarField& mu =
231  mesh_.lookupObject<volScalarField>(muName);
232 
233  apply(Udiag, Usource, V, geometricOneField(), mu/rho, U);
234  }
235  }
236 }
237 
238 
240 (
241  const fvVectorMatrix& UEqn,
242  volTensorField& AU
243 ) const
244 {
245  const volVectorField& U = UEqn.psi();
246 
247  word rhoName(IOobject::groupName(rhoName_, U.group()));
248  word muName(IOobject::groupName(muName_, U.group()));
249  word nuName(IOobject::groupName(nuName_, U.group()));
250 
251  if (UEqn.dimensions() == dimForce)
252  {
253  const volScalarField& rho = mesh_.lookupObject<volScalarField>(rhoName);
254  const volScalarField& mu = mesh_.lookupObject<volScalarField>(muName);
255 
256  apply(AU, rho, mu, U);
257  }
258  else
259  {
260  if (mesh_.foundObject<volScalarField>(nuName))
261  {
262  const volScalarField& nu =
263  mesh_.lookupObject<volScalarField>(nuName);
264 
265  apply(AU, geometricOneField(), nu, U);
266  }
267  else
268  {
269  const volScalarField& rho =
270  mesh_.lookupObject<volScalarField>(rhoName);
271  const volScalarField& mu =
272  mesh_.lookupObject<volScalarField>(muName);
273 
274  apply(AU, geometricOneField(), mu/rho, U);
275  }
276  }
277 }
278 
279 
281 {
282  os << indent << name_ << endl;
283  dict_.write(os);
284 
285  return true;
286 }
287 
288 
289 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A List with indirect addressing.
Definition: UIndirectList.H:60
Abstract base class for coordinate rotation.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Top level model for porosity models.
Definition: porosityModel.H:57
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:40
Darcy-Forchheimer law porosity model, given by:
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
bool writeData(Ostream &os) const
Write.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
DarcyForchheimer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName)
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for handling words, derived from string.
Definition: word.H:62
fvVectorMatrix & UEqn
Definition: UEqn.H:11
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const cellShapeList & cells
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimForce
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
dictionary dict