Loading...
Searching...
No Matches
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-2023 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 "fvMesh.H"
31#include "fvMatrices.H"
32#include "porosityModel.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace fv
40{
43 (
44 option,
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(nullptr);
97
100 (
101 name_,
102 nbrMesh,
103 coeffs_,
104 zoneName
105 );
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 auto tUNbr = volVectorField::New
150 (
151 IOobject::scopedName(name_, "UNbr"),
153 nbrMesh,
154 dimensionedVector(U.dimensions(), Zero)
155 );
156 auto& UNbr = tUNbr.ref();
157
158 // Map local velocity onto neighbour region
159 meshInterp().mapSrcToTgt
160 (
161 U.primitiveField(),
163 UNbr.primitiveFieldRef()
164 );
165
166 fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
167
168 porosityPtr_->addResistance(nbrEqn);
169
170 // Convert source from neighbour to local region
171 fvMatrix<vector> porosityEqn(U, eqn.dimensions());
172 scalarField& Udiag = porosityEqn.diag();
173 vectorField& Usource = porosityEqn.source();
174
175 Udiag.setSize(eqn.diag().size(), Zero);
176 Usource.setSize(eqn.source().size(), Zero);
177
178 meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
179 meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
180
181 eqn -= porosityEqn;
182}
183
184
186(
187 const volScalarField& rho,
188 fvMatrix<vector>& eqn,
189 const label fieldi
190)
191{
192 initialise();
193
194 const auto& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
195
196 const volVectorField& U = eqn.psi();
197
198 auto tUNbr = volVectorField::New
199 (
200 IOobject::scopedName(name_, "UNbr"),
202 nbrMesh,
203 dimensionedVector(U.dimensions(), Zero)
204 );
205 auto& UNbr = tUNbr.ref();
206
207 // Map local velocity onto neighbour region
208 meshInterp().mapSrcToTgt
209 (
210 U.primitiveField(),
212 UNbr.primitiveFieldRef()
213 );
214
215 fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
216
217 auto trhoNbr = volScalarField::New
218 (
219 IOobject::scopedName("rho", "UNbr"),
221 nbrMesh,
223 );
224 auto& rhoNbr = trhoNbr.ref();
225
226 auto tmuNbr = volScalarField::New
227 (
228 IOobject::scopedName("mu", "UNbr"),
230 nbrMesh,
232 );
233 auto& muNbr = tmuNbr.ref();
234
235 const auto& mu = mesh_.lookupObject<volScalarField>(muName_);
236
237 // Map local rho onto neighbour region
238 meshInterp().mapSrcToTgt
239 (
240 rho.primitiveField(),
242 rhoNbr.primitiveFieldRef()
243 );
244
245 // Map local mu onto neighbour region
246 meshInterp().mapSrcToTgt
247 (
248 mu.primitiveField(),
250 muNbr.primitiveFieldRef()
251 );
252
253 porosityPtr_->addResistance(nbrEqn, rhoNbr, muNbr);
254
255 // Convert source from neighbour to local region
256 fvMatrix<vector> porosityEqn(U, eqn.dimensions());
257 scalarField& Udiag = porosityEqn.diag();
258 vectorField& Usource = porosityEqn.source();
259
260 Udiag.setSize(eqn.diag().size(), Zero);
261 Usource.setSize(eqn.source().size(), Zero);
262
263 meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
264 meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
265
266 eqn -= porosityEqn;
267}
268
269
271{
273 {
274 coeffs_.readIfPresent("U", UName_);
275 coeffs_.readIfPresent("mu", muName_);
276
277 // Reset the porosity model?
278
279 return true;
280 }
281
282 return false;
283}
284
285
286// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
void setSize(label n)
Alias for resize().
Definition List.H:536
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
void clearAddressing()
Clear addressing.
Definition ZoneMesh.C:944
A subset of mesh cells.
Definition cellZone.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const dimensionSet & dimensions() const noexcept
Definition fvMatrix.H:530
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Field< Type > & source() noexcept
Definition fvMatrix.H:535
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Applies inter-region explicit porosity source.
virtual bool read(const dictionary &dict)
Read dictionary.
word muName_
Name of operand dynamic viscosity field (compressible case only).
autoPtr< porosityModel > porosityPtr_
Run-time selectable porosity model.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Vector.
interRegionExplicitPorositySource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Intermediate class for handling inter-region exchanges.
word nbrRegionName_
Name of the neighbour region to map.
interRegionOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from dictionary.
virtual bool read(const dictionary &dict)
Read dictionary.
const meshToMesh & meshInterp() const
Return const access to the mapToMap pointer.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
bool active_
Source active flag.
Definition fvOption.H:167
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
const scalarField & diag() const
Definition lduMatrix.C:195
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...
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const wordRe &cellZoneName=wordRe::null)
Selector.
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for finite-volume.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict