interRegionPorosityForce.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 "fvMesh.H"
28 #include "fvMatrices.H"
29 #include "porosityModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
41  (
42  fvModel,
44  dictionary,
45  interRegionExplicitPorositySource,
46  "interRegionExplicitPorositySource"
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
53 
54 void Foam::fv::interRegionPorosityForce::readCoeffs()
55 {
56  UName_ = coeffs().lookupOrDefault<word>("U", "U");
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const word& modelType,
66  const fvMesh& mesh,
67  const dictionary& dict
68 )
69 :
70  interRegionModel(name, modelType, mesh, dict),
71  UName_(word::null),
72  filter_
73  (
74  volScalarField::Internal::New
75  (
76  "filter",
77  mesh,
79  )
80  ),
81  porosityPtr_(nullptr)
82 {
83  readCoeffs();
84 
86 
87  interpolate(scalarField(nbrMesh.nCells(), 1), filter_);
88 
89  const word zoneName(name + ":porous");
90 
91  if (!mesh.cellZones().found(zoneName))
92  {
93  // Scan the porous region filter for all cells containing porosity
94  labelList porousCells(mesh.nCells());
95 
96  label i = 0;
97  forAll(filter_, celli)
98  {
99  if (filter_[celli] > small)
100  {
101  porousCells[i++] = celli;
102  }
103  }
104  porousCells.setSize(i);
105 
107  (
108  new cellZone
109  (
110  zoneName,
111  porousCells,
112  mesh.cellZones()
113  )
114  );
115  }
116  else
117  {
119  << "Unable to create porous cellZone " << zoneName
120  << ": zone already exists"
121  << abort(FatalError);
122  }
123 
124  porosityPtr_ = porosityModel::New
125  (
126  name,
127  mesh,
128  coeffs(),
129  zoneName
130  );
131 }
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 {
138  return wordList(1, UName_);
139 }
140 
141 
143 (
144  const volVectorField& U,
145  fvMatrix<vector>& eqn
146 ) const
147 {
148  fvMatrix<vector> porosityEqn(eqn.psi(), eqn.dimensions());
149  porosityPtr_->addResistance(porosityEqn);
150  eqn -= filter_*porosityEqn;
151 }
152 
153 
155 (
156  const volScalarField& rho,
157  const volVectorField& U,
158  fvMatrix<vector>& eqn
159 ) const
160 {
161  fvMatrix<vector> porosityEqn(eqn.psi(), eqn.dimensions());
162  porosityPtr_->addResistance(porosityEqn);
163  eqn -= filter_*porosityEqn;
164 }
165 
166 
168 (
169  const volScalarField& alpha,
170  const volScalarField& rho,
171  const volVectorField& U,
172  fvMatrix<vector>& eqn
173 ) const
174 {
175  fvMatrix<vector> porosityEqn(eqn.psi(), eqn.dimensions());
176  porosityPtr_->addResistance(porosityEqn);
177  eqn -= alpha*filter_*porosityEqn;
178 }
179 
180 
182 {
184  return true;
185 }
186 
187 
189 {
191 }
192 
193 
195 {
197 }
198 
199 
201 {
203 }
204 
205 
207 {
209  {
210  readCoeffs();
211  return true;
212  }
213  else
214  {
215  return false;
216  }
217 }
218 
219 
220 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
void setSize(const label)
Reset size of List.
Definition: List.C:281
void append(ZoneType *) const
Append or update a zone.
Definition: ZoneList.C:176
bool found(const label objectIndex) const
Return true if objectIndex is in any zone.
Definition: ZoneList.C:122
A subset of mesh cells.
Definition: cellZone.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
const dimensionSet & dimensions() const
Definition: fvMatrix.H:302
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
Base class for inter-region exchange.
virtual bool read(const dictionary &dict)
Read dictionary.
const word & nbrRegionName() const
Return const access to the neighbour region name.
const fvMesh & nbrMesh() const
Return const access to the neighbour mesh.
tmp< Field< Type > > interpolate(const Field< Type > &field) const
Interpolate field.
This model applies the force exerted on the fluid by a porous media, the extent of which is defined b...
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Add source to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
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 bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
interRegionPorosityForce(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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:452
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Selector.
label nCells() const
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#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.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
labelList fv(nPoints)
dictionary dict