interRegionExplicitPorositySource.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-2016 OpenFOAM Foundation
9  Copyright (C) 2020 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 
30 #include "fvMesh.H"
31 #include "fvMatrices.H"
32 #include "porosityModel.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
41  defineTypeNameAndDebug(interRegionExplicitPorositySource, 0);
43  (
44  option,
45  interRegionExplicitPorositySource,
47  );
48 }
49 }
50 
51 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
52 
54 {
55  if (!firstIter_)
56  {
57  return;
58  }
59 
60  const word zoneName(name_ + ":porous");
61 
62  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
63  const cellZoneMesh& cellZones = nbrMesh.cellZones();
64  label zoneID = cellZones.findZoneID(zoneName);
65 
66  if (zoneID == -1)
67  {
68  cellZoneMesh& cz = const_cast<cellZoneMesh&>(cellZones);
69 
70  zoneID = cz.size();
71 
72  cz.setSize(zoneID + 1);
73 
74  cz.set
75  (
76  zoneID,
77  new cellZone
78  (
79  zoneName,
80  nbrMesh.faceNeighbour(), // Neighbour internal cells
81  zoneID,
82  cellZones
83  )
84  );
85 
86  cz.clearAddressing();
87  }
88  else
89  {
91  << "Unable to create porous cellZone " << zoneName
92  << ": zone already exists"
93  << abort(FatalError);
94  }
95 
96  porosityPtr_.reset
97  (
99  (
100  name_,
101  nbrMesh,
102  coeffs_,
103  zoneName
104  ).ptr()
105  ),
106 
107  firstIter_ = false;
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
114 (
115  const word& name,
116  const word& modelType,
117  const dictionary& dict,
118  const fvMesh& mesh
119 )
120 :
121  interRegionOption(name, modelType, dict, mesh),
122  porosityPtr_(nullptr),
123  firstIter_(true),
124  UName_(coeffs_.getOrDefault<word>("U", "U")),
125  muName_(coeffs_.getOrDefault<word>("mu", "thermo:mu"))
126 {
127  if (active_)
128  {
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 (
139  fvMatrix<vector>& eqn,
140  const label fieldi
141 )
142 {
143  initialise();
144 
145  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
146 
147  const volVectorField& U = eqn.psi();
148 
149  volVectorField UNbr
150  (
151  IOobject
152  (
153  name_ + ":UNbr",
154  nbrMesh.time().timeName(),
155  nbrMesh,
158  ),
159  nbrMesh,
160  dimensionedVector(U.dimensions(), Zero)
161  );
162 
163  // Map local velocity onto neighbour region
164  meshInterp().mapSrcToTgt
165  (
166  U.primitiveField(),
168  UNbr.primitiveFieldRef()
169  );
170 
171  fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
172 
173  porosityPtr_->addResistance(nbrEqn);
174 
175  // Convert source from neighbour to local region
176  fvMatrix<vector> porosityEqn(U, eqn.dimensions());
177  scalarField& Udiag = porosityEqn.diag();
178  vectorField& Usource = porosityEqn.source();
179 
180  Udiag.setSize(eqn.diag().size(), 0.0);
181  Usource.setSize(eqn.source().size(), Zero);
182 
183  meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
184  meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
185 
186  eqn -= porosityEqn;
187 }
188 
189 
191 (
192  const volScalarField& rho,
193  fvMatrix<vector>& eqn,
194  const label fieldi
195 )
196 {
197  initialise();
198 
199  const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
200 
201  const volVectorField& U = eqn.psi();
202 
203  volVectorField UNbr
204  (
205  IOobject
206  (
207  name_ + ":UNbr",
208  nbrMesh.time().timeName(),
209  nbrMesh,
212  ),
213  nbrMesh,
214  dimensionedVector(U.dimensions(), Zero)
215  );
216 
217  // Map local velocity onto neighbour region
218  meshInterp().mapSrcToTgt
219  (
220  U.primitiveField(),
222  UNbr.primitiveFieldRef()
223  );
224 
225  fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
226 
227  volScalarField rhoNbr
228  (
229  IOobject
230  (
231  "rho:UNbr",
232  nbrMesh.time().timeName(),
233  nbrMesh,
236  ),
237  nbrMesh,
239  );
240 
241  volScalarField muNbr
242  (
243  IOobject
244  (
245  "mu:UNbr",
246  nbrMesh.time().timeName(),
247  nbrMesh,
250  ),
251  nbrMesh,
253  );
254 
255  const auto& mu = mesh_.lookupObject<volScalarField>(muName_);
256 
257  // Map local rho onto neighbour region
258  meshInterp().mapSrcToTgt
259  (
260  rho.primitiveField(),
262  rhoNbr.primitiveFieldRef()
263  );
264 
265  // Map local mu onto neighbour region
266  meshInterp().mapSrcToTgt
267  (
268  mu.primitiveField(),
270  muNbr.primitiveFieldRef()
271  );
272 
273  porosityPtr_->addResistance(nbrEqn, rhoNbr, muNbr);
274 
275  // Convert source from neighbour to local region
276  fvMatrix<vector> porosityEqn(U, eqn.dimensions());
277  scalarField& Udiag = porosityEqn.diag();
278  vectorField& Usource = porosityEqn.source();
279 
280  Udiag.setSize(eqn.diag().size(), 0.0);
281  Usource.setSize(eqn.source().size(), Zero);
282 
283  meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
284  meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
285 
286  eqn -= porosityEqn;
287 }
288 
289 
291 {
293  {
294  coeffs_.readIfPresent("U", UName_);
295  coeffs_.readIfPresent("mu", muName_);
296 
297  // Reset the porosity model?
298 
299  return true;
300  }
301 
302  return false;
303 }
304 
305 
306 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Intermediate class for handling inter-region exchanges.
virtual bool read(const dictionary &dict)
Read dictionary.
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
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
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Vector.
const dimensionSet dimViscosity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const wordRe &cellZoneName=wordRe::null)
Selector.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &dict)
Read dictionary.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
const word name_
Source name.
Definition: fvOption.H:132
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
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:64
word nbrRegionName_
Name of the neighbour region to map.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
const dimensionSet dimDensity
autoPtr< porosityModel > porosityPtr_
Run-time selectable porosity model.
const dimensionedScalar mu
Atomic mass unit.
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:528
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
Nothing to be read.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
scalarField & diag()
Definition: lduMatrix.C:197
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
bool active_
Source active flag.
Definition: fvOption.H:167
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127