Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
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 (
100 )
101{
102 if (zoneName_.empty())
103 {
104 dict.readIfPresent("active", active_);
105 dict_.readEntry("cellZone", zoneName_);
106 }
107
108 cellZoneIDs_ = mesh_.cellZones().indices(zoneName_);
109
110 Info<< " creating porous zone: " << zoneName_ << endl;
111
113 {
114 FatalErrorInFunction
115 << "Cannot find porous cellZone " << zoneName_ << endl
116 << "Valid zones : "
117 << flatOutput(mesh_.cellZones().names()) << nl
118 << "Valid groups: "
119 << flatOutput(mesh_.cellZones().groupNames()) << nl
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
167Foam::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 auto tforce = tmp<vectorField>::New(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;
193
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;
210
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.
236 AU.correctBoundaryConditions();
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// ************************************************************************* //
if(maxValue - minValue< SMALL)
if(patchID !=-1)
void correctBoundaryConditions()
Correct boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static constexpr direction nComponents
Number of components in this vector space.
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
void reduce()
Inplace parallel reduction of min/max values, using UPstream::worldComm.
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
vector span() const
The bounding box span (from minimum to maximum).
Definition boundBoxI.H:192
A subset of mesh cells.
Definition cellZone.H:61
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
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:133
const Type & value() const noexcept
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Top level model for porosity models.
virtual void addResistance(fvVectorMatrix &UEqn)
Add resistance.
virtual void calcTransformModelData()=0
Transform the model data wrt mesh changes.
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const =0
Calculate the porosity force.
bool active_
Porosity active flag.
const fvMesh & mesh_
Reference to the mesh database.
virtual bool writeData(Ostream &os) const
Write.
word name_
Porosity name.
autoPtr< coordinateSystem > csysPtr_
Local coordinate system.
dictionary coeffs_
Model coefficients dictionary.
porosityModel(const porosityModel &)=delete
No copy construct.
virtual void transformModelData()
Transform the model data wrt mesh changes.
const dictionary dict_
Dictionary used for model construction.
wordRe zoneName_
Name(s) of cell-zone.
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.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const wordRe &cellZoneName=wordRe::null)
Selector.
const word & name() const
Return const access to the porosity model name.
labelList cellZoneIDs_
Cell zone IDs.
virtual bool read()
Inherit read from regIOobject.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
thermo correct()
fvVectorMatrix & UEqn
Definition UEqn.H:13
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
GeometricField< tensor, fvPatchField, volMesh > volTensorField
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
fvMatrix< vector > fvVectorMatrix
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vector point
Point is a vector.
Definition point.H:37
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:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict