fixedCoeff.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-2024 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 "fixedCoeff.H"
28 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace porosityModels
35  {
38  }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::porosityModels::fixedCoeff::apply
45 (
46  scalarField& Udiag,
47  vectorField& Usource,
48  const scalarField& V,
49  const vectorField& U,
50  const scalar rho
51 ) const
52 {
54 
55  forAll(cells, i)
56  {
57  const label celli = cells[i];
58  const label j = fieldIndex(i);
59  const tensor Cd = rho*(alpha_[j] + beta_[j]*mag(U[celli]));
60  const scalar isoCd = tr(Cd);
61 
62  Udiag[celli] += V[celli]*isoCd;
63  Usource[celli] -= V[celli]*((Cd - I*isoCd) & U[celli]);
64  }
65 }
66 
67 
68 void Foam::porosityModels::fixedCoeff::apply
69 (
70  tensorField& AU,
71  const vectorField& U,
72  const scalar rho
73 ) const
74 {
75  const labelList& cells = mesh_.cellZones()[zoneName_];
76 
77  forAll(cells, i)
78  {
79  const label celli = cells[i];
80  const label j = fieldIndex(i);
81  const tensor alpha = alpha_[j];
82  const tensor beta = beta_[j];
83 
84  AU[celli] += rho*(alpha + beta*mag(U[celli]));
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const word& name,
94  const fvMesh& mesh,
95  const dictionary& dict,
96  const dictionary& coeffDict,
97  const word& cellZoneName
98 )
99 :
100  porosityModel(name, mesh, dict, coeffDict, cellZoneName),
101  alphaXYZ_("alpha", dimless/dimTime, coeffDict),
102  betaXYZ_("beta", dimless/dimLength, coeffDict),
103  rhoRefFound_(coeffDict.found("rhoRef")),
104  rhoRef_(coeffDict.lookupOrDefault<scalar>("rhoRef", 1.0))
105 {
106  adjustNegativeResistance(alphaXYZ_);
107  adjustNegativeResistance(betaXYZ_);
108 
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
114 
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  if (coordSys_.R().uniform())
124  {
125  alpha_.setSize(1);
126  beta_.setSize(1);
127 
128  alpha_[0] = Zero;
129  alpha_[0].xx() = alphaXYZ_.value().x();
130  alpha_[0].yy() = alphaXYZ_.value().y();
131  alpha_[0].zz() = alphaXYZ_.value().z();
132  alpha_[0] = coordSys_.R().transform(Zero, alpha_[0]);
133 
134  beta_[0] = Zero;
135  beta_[0].xx() = betaXYZ_.value().x();
136  beta_[0].yy() = betaXYZ_.value().y();
137  beta_[0].zz() = betaXYZ_.value().z();
138  beta_[0] = coordSys_.R().transform(Zero, beta_[0]);
139  }
140  else
141  {
142  const labelList& cells = mesh_.cellZones()[zoneName_];
143 
144  alpha_.setSize(cells.size());
145  beta_.setSize(cells.size());
146 
147  forAll(cells, i)
148  {
149  alpha_[i] = Zero;
150  alpha_[i].xx() = alphaXYZ_.value().x();
151  alpha_[i].yy() = alphaXYZ_.value().y();
152  alpha_[i].zz() = alphaXYZ_.value().z();
153 
154  beta_[i] = Zero;
155  beta_[i].xx() = betaXYZ_.value().x();
156  beta_[i].yy() = betaXYZ_.value().y();
157  beta_[i].zz() = betaXYZ_.value().z();
158  }
159 
160  const coordinateRotation& R = coordSys_.R
161  (
162  UIndirectList<vector>(mesh_.C(), cells)()
163  );
164 
165  alpha_ = R.transform(alpha_);
166  beta_ = R.transform(beta_);
167  }
168 }
169 
170 
172 (
173  const volVectorField& U,
174  const volScalarField& rho,
175  const volScalarField& mu,
176  vectorField& force
177 ) const
178 {
179  scalarField Udiag(U.size(), 0.0);
180  vectorField Usource(U.size(), Zero);
181  const scalarField& V = mesh_.V();
182 
183  if (!rhoRefFound_)
184  {
186  << "rhoRef not specified" << exit(FatalError);
187  }
188 
189  apply(Udiag, Usource, V, U, rhoRef_);
190 
191  force = Udiag*U - Usource;
192 }
193 
194 
196 (
198 ) const
199 {
200  const vectorField& U = UEqn.psi();
201  const scalarField& V = mesh_.V();
202  scalarField& Udiag = UEqn.diag();
203  vectorField& Usource = UEqn.source();
204 
205  if (UEqn.dimensions() == dimForce && !rhoRefFound_)
206  {
208  << "rhoRef not specified" << exit(FatalError);
209  }
210 
211  apply(Udiag, Usource, V, U, rhoRef_);
212 }
213 
214 
216 (
217  const fvVectorMatrix& UEqn,
218  volTensorField& AU
219 ) const
220 {
221  const vectorField& U = UEqn.psi();
222 
223  if (UEqn.dimensions() == dimForce && !rhoRefFound_)
224  {
226  << "rhoRef not specified" << exit(FatalError);
227  }
228 
229  apply(AU, U, rhoRef_);
230 }
231 
232 
233 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A List with indirect addressing.
Definition: UIndirectList.H:60
Abstract base class for coordinate rotation.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:96
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
Top level model for porosity models.
Definition: porosityModel.H:57
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:75
word zoneName_
Name of cellZone.
Definition: porosityModel.H:78
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:40
label fieldIndex(const label index) const
Return label index.
Definition: porosityModel.C:66
Fixed coefficient form of porosity model.
Definition: fixedCoeff.H:61
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
Definition: fixedCoeff.C:172
fixedCoeff(const word &name, const fvMesh &mesh, const dictionary &dict, const dictionary &coeffDict, const word &cellZoneName)
Definition: fixedCoeff.C:92
virtual ~fixedCoeff()
Destructor.
Definition: fixedCoeff.C:115
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
Definition: fixedCoeff.C:196
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
Definition: fixedCoeff.C:121
A class for handling words, derived from string.
Definition: word.H:62
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
A special matrix type and solver, designed for finite volume solutions of scalar equations.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar mu
Atomic mass unit.
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
const dimensionSet dimless
const dimensionSet dimLength
static const Identity< scalar > I
Definition: Identity.H:93
const dimensionSet dimForce
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
Field< vector > vectorField
Specialisation of Field<T> for vector.
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
dictionary dict