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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
34namespace 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:
80 (
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
112 bool foundZone = !cellZoneIDs_.empty();
113 reduce(foundZone, orOp<bool>());
114
115 if (!foundZone && Pstream::master())
116 {
118 << "Cannot find porous cellZone " << zoneName_ << endl
119 << "Valid zones : "
121 << "Valid groups: "
123 << exit(FatalError);
124 }
125
126 Info<< incrIndent << indent << csys() << decrIndent << endl;
127
128 const pointField& points = mesh_.points();
129 const cellList& cells = mesh_.cells();
130 const faceList& faces = mesh_.faces();
131
132 for (const label zonei : cellZoneIDs_)
133 {
134 const cellZone& cZone = mesh_.cellZones()[zonei];
135
136 boundBox bb;
137
138 for (const label celli : cZone)
139 {
140 const cell& c = cells[celli];
141 const pointField cellPoints(c.points(faces, points));
142
143 for (const point& pt : cellPoints)
144 {
145 bb.add(csys().localPosition(pt));
146 }
147 }
148
149 bb.reduce();
150
151 Info<< " local bounds: " << bb.span() << nl << endl;
152 }
153}
154
155
156// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157
159{
160 if (!mesh_.upToDatePoints(*this))
161 {
162 calcTransformModelData();
163
164 // set model up-to-date wrt points
165 mesh_.setUpToDatePoints(*this);
166 }
167}
168
169
170Foam::tmp<Foam::vectorField> Foam::porosityModel::porosityModel::force
171(
172 const volVectorField& U,
173 const volScalarField& rho,
174 const volScalarField& mu
175)
176{
177 transformModelData();
178
179 tmp<vectorField> tforce(new vectorField(U.size(), Zero));
180
181 if (!cellZoneIDs_.empty())
182 {
183 this->calcForce(U, rho, mu, tforce.ref());
184 }
185
186 return tforce;
187}
188
189
191{
192 if (cellZoneIDs_.empty())
193 {
194 return;
195 }
196
197 transformModelData();
198 this->correct(UEqn);
199}
200
201
203(
205 const volScalarField& rho,
206 const volScalarField& mu
207)
208{
209 if (cellZoneIDs_.empty())
210 {
211 return;
212 }
213
214 transformModelData();
215 this->correct(UEqn, rho, mu);
216}
217
218
220(
221 const fvVectorMatrix& UEqn,
222 volTensorField& AU,
223 bool correctAUprocBC
224)
225{
226 if (cellZoneIDs_.empty())
227 {
228 return;
229 }
230
231 transformModelData();
232 this->correct(UEqn, AU);
233
234 if (correctAUprocBC)
235 {
236 // Correct the boundary conditions of the tensorial diagonal to ensure
237 // processor boundaries are correctly handled when AU^-1 is interpolated
238 // for the pressure equation.
240 }
241}
242
243
245{
246 return true;
247}
248
249
251{
252 dict.readIfPresent("active", active_);
253
254 coeffs_ = dict.optionalSubDict(type() + "Coeffs");
255
256 dict.readEntry("cellZone", zoneName_);
257 cellZoneIDs_ = mesh_.cellZones().indices(zoneName_);
258
259 return true;
260}
261
262
263// ************************************************************************* //
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
wordList groupNames() const
A list of the zone group names (if any)
Definition: ZoneMesh.C:311
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:377
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
void reduce()
Parallel reduction of min/max values.
Definition: boundBox.C:184
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:191
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
A subset of mesh cells.
Definition: cellZone.H:65
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Base class for coordinate system specification, the default coordinate system type is cartesian .
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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:577
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
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:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
Top level model for porosity models.
Definition: porosityModel.H:60
virtual void addResistance(fvVectorMatrix &UEqn)
Add resistance.
bool active_
Porosity active flag.
Definition: porosityModel.H:78
const fvMesh & mesh_
Reference to the mesh database.
Definition: porosityModel.H:69
virtual bool writeData(Ostream &os) const
Write.
virtual void transformModelData()
Transform the model data wrt mesh changes.
const dictionary dict_
Dictionary used for model construction.
Definition: porosityModel.H:72
wordRe zoneName_
Name(s) of cell-zone.
Definition: porosityModel.H:81
const coordinateSystem & csys() const
Local coordinate system.
const dictionary & dict() const
Return dictionary used for model construction.
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
Definition: porosityModel.C:43
labelList cellZoneIDs_
Cell zone IDs.
Definition: porosityModel.H:84
virtual bool read()
Inherit read from regIOobject.
const cellList & cells() const
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
thermo correct()
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & mu
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
Field< vector > vectorField
Specialisation of Field<T> for vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict