interRegionExplicitPorositySource.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-2021 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 {
38  defineTypeNameAndDebug(interRegionExplicitPorositySource, 0);
40  (
41  fvModel,
42  interRegionExplicitPorositySource,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
50 
51 void Foam::fv::interRegionExplicitPorositySource::readCoeffs()
52 {
53  UName_ = coeffs().lookupOrDefault<word>("U", "U");
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const word& name,
62  const word& modelType,
63  const dictionary& dict,
64  const fvMesh& mesh
65 )
66 :
67  interRegionModel(name, modelType, dict, mesh),
68  UName_(word::null),
69  filter_
70  (
72  (
73  "filter",
74  mesh,
76  )
77  ),
78  porosityPtr_(nullptr)
79 {
80  readCoeffs();
81 
82  const fvMesh& nbrMesh = mesh.time().lookupObject<fvMesh>(nbrRegionName());
83 
84  meshInterp().mapTgtToSrc
85  (
86  scalarField(nbrMesh.nCells(), 1),
88  filter_
89  );
90 
91  const word zoneName(name + ":porous");
92 
93  const meshCellZones& cellZones = mesh.cellZones();
94  label zoneID = cellZones.findZoneID(zoneName);
95 
96  if (zoneID == -1)
97  {
98  meshCellZones& cz = const_cast<meshCellZones&>(cellZones);
99 
100  zoneID = cz.size();
101 
102  cz.setSize(zoneID + 1);
103 
104  // Scan the porous region filter for all cells containing porosity
105  labelList porousCells(mesh.nCells());
106 
107  label i = 0;
108  forAll(filter_, celli)
109  {
110  if (filter_[celli] > small)
111  {
112  porousCells[i++] = celli;
113  }
114  }
115  porousCells.setSize(i);
116 
117  cz.set
118  (
119  zoneID,
120  new cellZone
121  (
122  zoneName,
123  porousCells,
124  zoneID,
125  cellZones
126  )
127  );
128 
129  cz.clearAddressing();
130  }
131  else
132  {
134  << "Unable to create porous cellZone " << zoneName
135  << ": zone already exists"
136  << abort(FatalError);
137  }
138 
139  porosityPtr_ = porosityModel::New
140  (
141  name,
142  mesh,
143  coeffs(),
144  zoneName
145  );
146 }
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
153 {
154  return wordList(1, UName_);
155 }
156 
157 
159 (
160  const volScalarField& rho,
161  fvMatrix<vector>& eqn,
162  const word& fieldName
163 ) const
164 {
165  fvMatrix<vector> porosityEqn(eqn.psi(), eqn.dimensions());
166  porosityPtr_->addResistance(porosityEqn);
167  eqn -= filter_*porosityEqn;
168 }
169 
170 
172 {
173  if (interRegionModel::read(dict))
174  {
175  readCoeffs();
176  return true;
177  }
178  else
179  {
180  return false;
181  }
182 }
183 
184 
185 // ************************************************************************* //
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
virtual bool read(const dictionary &dict)
Read dictionary.
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:482
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
label nCells() const
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
Base class for inter-region exchange.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const word null
An empty word.
Definition: word.H:77
interRegionExplicitPorositySource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
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:72
virtual bool read(const dictionary &dict)
Read dictionary.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
A subset of mesh cells.
Definition: cellZone.H:61
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite volume solutions of scalar equations.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Selector.
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
Namespace for OpenFOAM.
virtual void addSup(const volScalarField &rho, fvMatrix< vector > &eqn, const word &fieldName) const
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287