porosityModel.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-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 
29 #include "porosityModel.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
38 }
39 
40 
41 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
42 
44 {
45  scalar maxCmpt = cmptMax(resist.value());
46 
47  if (maxCmpt < 0)
48  {
50  << "Cannot have all resistances set to negative, resistance = "
51  << resist
52  << exit(FatalError);
53  }
54  else
55  {
56  vector& val = resist.value();
57  for (label cmpt = 0; cmpt < vector::nComponents; cmpt++)
58  {
59  if (val[cmpt] < 0)
60  {
61  val[cmpt] *= -maxCmpt;
62  }
63  }
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const word& name,
73  const word& modelType,
74  const fvMesh& mesh,
75  const dictionary& dict,
76  const wordRe& cellZoneName
77 )
78 :
79  regIOobject
80  (
81  IOobject
82  (
83  name,
84  mesh.time().timeName(),
85  mesh,
86  IOobject::NO_READ,
87  IOobject::NO_WRITE
88  )
89  ),
90  name_(name),
91  mesh_(mesh),
92  dict_(dict),
93  coeffs_(dict.optionalSubDict(modelType + "Coeffs")),
94  active_(true),
95  zoneName_(cellZoneName),
96  cellZoneIDs_(),
97  csysPtr_
98  (
99  coordinateSystem::New(mesh, coeffs_, coordinateSystem::typeName)
100  )
101 {
102  if (zoneName_.empty())
103  {
104  dict.readIfPresent("active", active_);
105  dict_.readEntry("cellZone", zoneName_);
106  }
107 
109 
110  Info<< " creating porous zone: " << zoneName_ << endl;
111 
113  {
115  << "Cannot find porous cellZone " << zoneName_ << endl
116  << "Valid zones : "
117  << flatOutput(mesh_.cellZones().names()) << nl
118  << "Valid groups: "
120  << exit(FatalError);
121  }
122 
123  Info<< incrIndent << indent << csys() << decrIndent << endl;
124 
125  const pointField& points = mesh_.points();
126  const cellList& cells = mesh_.cells();
127  const faceList& faces = mesh_.faces();
128 
129  for (const label zonei : cellZoneIDs_)
130  {
131  const cellZone& cZone = mesh_.cellZones()[zonei];
132 
133  boundBox bb;
134 
135  for (const label celli : cZone)
136  {
137  const cell& c = cells[celli];
138  const pointField cellPoints(c.points(faces, points));
139 
140  for (const point& pt : cellPoints)
141  {
142  bb.add(csys().localPosition(pt));
143  }
144  }
145 
146  bb.reduce();
147 
148  Info<< " local bounds: " << bb.span() << nl << endl;
149  }
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
156 {
157  if (!mesh_.upToDatePoints(*this))
158  {
159  calcTransformModelData();
160 
161  // set model up-to-date wrt points
162  mesh_.setUpToDatePoints(*this);
163  }
164 }
165 
166 
167 Foam::tmp<Foam::vectorField> Foam::porosityModel::porosityModel::force
168 (
169  const volVectorField& U,
170  const volScalarField& rho,
171  const volScalarField& mu
172 )
173 {
174  transformModelData();
175 
176  tmp<vectorField> tforce(new vectorField(U.size(), Zero));
177 
178  if (!cellZoneIDs_.empty())
179  {
180  this->calcForce(U, rho, mu, tforce.ref());
181  }
182 
183  return tforce;
184 }
185 
186 
188 {
189  if (cellZoneIDs_.empty())
190  {
191  return;
192  }
194  transformModelData();
195  this->correct(UEqn);
196 }
197 
198 
200 (
202  const volScalarField& rho,
203  const volScalarField& mu
204 )
205 {
206  if (cellZoneIDs_.empty())
207  {
208  return;
209  }
211  transformModelData();
212  this->correct(UEqn, rho, mu);
213 }
214 
215 
217 (
218  const fvVectorMatrix& UEqn,
219  volTensorField& AU,
220  bool correctAUprocBC
221 )
222 {
223  if (cellZoneIDs_.empty())
224  {
225  return;
226  }
227 
228  transformModelData();
229  this->correct(UEqn, AU);
230 
231  if (correctAUprocBC)
232  {
233  // Correct the boundary conditions of the tensorial diagonal to ensure
234  // processor boundaries are correctly handled when AU^-1 is interpolated
235  // for the pressure equation.
237  }
238 }
239 
242 {
243  return true;
244 }
245 
246 
248 {
249  dict.readIfPresent("active", active_);
250 
251  coeffs_ = dict.optionalSubDict(type() + "Coeffs");
252 
253  dict.readEntry("cellZone", zoneName_);
254  cellZoneIDs_ = mesh_.cellZones().indices(zoneName_);
255 
256  return true;
257 }
258 
259 
260 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
virtual bool writeData(Ostream &os) const
Write.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const dictionary & dict() const
Return dictionary used for model construction.
virtual bool read()
Read object.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const cellList & cells() const
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
bool active_
Porosity active flag.
Definition: porosityModel.H:84
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
porosityModel(const porosityModel &)=delete
No copy construct.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
word timeName
Definition: getTimeIndex.H:3
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:560
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
wordRe zoneName_
Name(s) of cell-zone.
Definition: porosityModel.H:89
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:36
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
const pointField & points
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:435
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
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:69
virtual void transformModelData()
Transform the model data wrt mesh changes.
const coordinateSystem & csys() const
Local coordinate system.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:94
const dictionary dict_
Dictionary used for model construction.
Definition: porosityModel.H:74
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const dimensionedScalar c
Speed of light in a vacuum.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
wordList groupNames() const
A list of the zone group names (if any)
Definition: ZoneMesh.C:369
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
void correctBoundaryConditions()
Correct boundary field.
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:362
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
virtual void addResistance(fvVectorMatrix &UEqn)
Add resistance.
Top level model for porosity models.
Definition: porosityModel.H:53
Namespace for OpenFOAM.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127