porosityModel.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 
26 #include "porosityModel.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
39 
41 {
42  scalar maxCmpt = cmptMax(resist.value());
43 
44  if (maxCmpt < 0)
45  {
47  << "Cannot have all resistances set to negative, resistance = "
48  << resist
49  << exit(FatalError);
50  }
51  else
52  {
53  maxCmpt = max(0, maxCmpt);
54  vector& val = resist.value();
55  for (label cmpt = 0; cmpt < vector::nComponents; cmpt++)
56  {
57  if (val[cmpt] < 0)
58  {
59  val[cmpt] *= -maxCmpt;
60  }
61  }
62  }
63 }
64 
65 
67 {
68  label index = 0;
69  if (!coordSys_.R().uniform())
70  {
71  index = i;
72  }
73  return index;
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 Foam::porosityModel::porosityModel
80 (
81  const word& name,
82  const fvMesh& mesh,
83  const dictionary& dict,
84  const dictionary& coeffDict,
85  const word& cellZoneName
86 )
87 :
89  (
90  IOobject
91  (
92  name,
93  mesh.time().name(),
94  mesh,
95  IOobject::NO_READ,
96  IOobject::NO_WRITE
97  )
98  ),
99  name_(name),
100  mesh_(mesh),
101  zoneName_
102  (
103  cellZoneName != word::null
104  ? cellZoneName
105  : dict.lookup<word>("cellZone")
106  ),
107  coordSys_(coordinateSystem::New(mesh, coeffDict))
108 {
109  Info<< " creating porous zone: " << zoneName_ << endl;
110 
111  if (mesh_.cellZones().findIndex(zoneName_) == -1)
112  {
114  << "cannot find porous cellZone " << zoneName_
115  << exit(FatalError);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 Foam::tmp<Foam::vectorField> Foam::porosityModel::porosityModel::force
129 (
130  const volVectorField& U,
131  const volScalarField& rho,
132  const volScalarField& mu
133 ) const
134 {
135  tmp<vectorField> tforce(new vectorField(U.size(), Zero));
136  this->calcForce(U, rho, mu, tforce.ref());
137 
138  return tforce;
139 }
140 
141 
143 {
144  this->correct(UEqn);
145 }
146 
147 
149 (
150  const fvVectorMatrix& UEqn,
151  volTensorField& AU,
152  bool correctAUprocBC
153 )
154 {
155  this->correct(UEqn, AU);
156 
157  if (correctAUprocBC)
158  {
159  // Correct the boundary conditions of the tensorial diagonal to
160  // ensure processor boundaries are correctly handled when AU^-1 is
161  // interpolated for the pressure equation.
163  }
164 }
165 
166 
168 {
169  calcTransformModelData();
170 
171  return true;
172 }
173 
174 
176 {
177  calcTransformModelData();
178 }
179 
180 
182 {
183  calcTransformModelData();
184 }
185 
186 
188 (
189  const polyDistributionMap& map
190 )
191 {
192  calcTransformModelData();
193 }
194 
195 
197 {
198  dict.lookup("cellZone") >> zoneName_;
199 
200  return true;
201 }
202 
203 
204 // ************************************************************************* //
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:105
Base class for other coordinate system specifications.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const Type & value() const
Return const reference to value.
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
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Top level model for porosity models.
Definition: porosityModel.H:57
virtual void addResistance(fvVectorMatrix &UEqn)
Add resistance.
virtual bool movePoints()
Update for mesh motion.
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:75
virtual ~porosityModel()
Destructor.
word zoneName_
Name of cellZone.
Definition: porosityModel.H:78
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:40
virtual bool read()
Inherit read from regIOobject.
label fieldIndex(const label index) const
Return label index.
Definition: porosityModel.C:66
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
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
U
Definition: pEqn.H:72
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
const dimensionedScalar mu
Atomic mass unit.
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
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dictionary dict