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 | Copyright (C) 2012-2015 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  option,
42  interRegionExplicitPorositySource,
43  dictionary
44  );
45 }
46 }
47 
48 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
49 
51 {
52  if (!firstIter_)
53  {
54  return;
55  }
56 
57  const word zoneName(name_ + ":porous");
58 
59  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
60  const cellZoneMesh& cellZones = nbrMesh.cellZones();
61  label zoneID = cellZones.findZoneID(zoneName);
62 
63  if (zoneID == -1)
64  {
65  cellZoneMesh& cz = const_cast<cellZoneMesh&>(cellZones);
66 
67  zoneID = cz.size();
68 
69  cz.setSize(zoneID + 1);
70 
71  cz.set
72  (
73  zoneID,
74  new cellZone
75  (
76  zoneName,
77  nbrMesh.faceNeighbour(), // neighbour internal cells
78  zoneID,
79  cellZones
80  )
81  );
82 
83  cz.clearAddressing();
84  }
85  else
86  {
88  << "Unable to create porous cellZone " << zoneName
89  << ": zone already exists"
90  << abort(FatalError);
91  }
92 
93  porosityPtr_.reset
94  (
96  (
97  name_,
98  nbrMesh,
99  coeffs_,
100  zoneName
101  ).ptr()
102  ),
103 
104  firstIter_ = false;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const word& name,
113  const word& modelType,
114  const dictionary& dict,
115  const fvMesh& mesh
116 )
117 :
118  interRegionOption(name, modelType, dict, mesh),
119  porosityPtr_(NULL),
120  firstIter_(-1),
121  UName_(coeffs_.lookupOrDefault<word>("UName", "U")),
122  muName_(coeffs_.lookupOrDefault<word>("muName", "thermo:mu"))
123 {
124  if (active_)
125  {
126  fieldNames_.setSize(1, UName_);
127  applied_.setSize(1, false);
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 (
136  fvMatrix<vector>& eqn,
137  const label fieldI
138 )
139 {
140  initialise();
141 
142  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
143 
144  const volVectorField& U = eqn.psi();
145 
146  volVectorField UNbr
147  (
148  IOobject
149  (
150  name_ + ":UNbr",
151  nbrMesh.time().timeName(),
152  nbrMesh,
155  ),
156  nbrMesh,
157  dimensionedVector("zero", U.dimensions(), vector::zero)
158  );
159 
160  // map local velocity onto neighbour region
161  meshInterp().mapSrcToTgt
162  (
163  U.internalField(),
165  UNbr.internalField()
166  );
167 
168  fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
169 
170  porosityPtr_->addResistance(nbrEqn);
171 
172  // convert source from neighbour to local region
173  fvMatrix<vector> porosityEqn(U, eqn.dimensions());
174  scalarField& Udiag = porosityEqn.diag();
175  vectorField& Usource = porosityEqn.source();
176 
177  Udiag.setSize(eqn.diag().size(), 0.0);
178  Usource.setSize(eqn.source().size(), vector::zero);
179 
180  meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
181  meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
182 
183  eqn -= porosityEqn;
184 }
185 
186 
188 (
189  const volScalarField& rho,
190  fvMatrix<vector>& eqn,
191  const label fieldI
192 )
193 {
194  initialise();
195 
196  const fvMesh& nbrMesh = mesh_.time().lookupObject<fvMesh>(nbrRegionName_);
197 
198  const volVectorField& U = eqn.psi();
199 
200  volVectorField UNbr
201  (
202  IOobject
203  (
204  name_ + ":UNbr",
205  nbrMesh.time().timeName(),
206  nbrMesh,
209  ),
210  nbrMesh,
211  dimensionedVector("zero", U.dimensions(), vector::zero)
212  );
213 
214  // map local velocity onto neighbour region
215  meshInterp().mapSrcToTgt
216  (
217  U.internalField(),
219  UNbr.internalField()
220  );
221 
222  fvMatrix<vector> nbrEqn(UNbr, eqn.dimensions());
223 
224  volScalarField rhoNbr
225  (
226  IOobject
227  (
228  "rho:UNbr",
229  nbrMesh.time().timeName(),
230  nbrMesh,
233  ),
234  nbrMesh,
235  dimensionedScalar("zero", dimDensity, 0.0)
236  );
237 
238  volScalarField muNbr
239  (
240  IOobject
241  (
242  "mu:UNbr",
243  nbrMesh.time().timeName(),
244  nbrMesh,
247  ),
248  nbrMesh,
249  dimensionedScalar("zero", dimViscosity, 0.0)
250  );
251 
252  const volScalarField& mu =
253  mesh_.lookupObject<volScalarField>(muName_);
254 
255  // map local rho onto neighbour region
256  meshInterp().mapSrcToTgt
257  (
258  rho.internalField(),
260  rhoNbr.internalField()
261  );
262 
263  // map local mu onto neighbour region
264  meshInterp().mapSrcToTgt
265  (
266  mu.internalField(),
268  muNbr.internalField()
269  );
270 
271  porosityPtr_->addResistance(nbrEqn, rhoNbr, muNbr);
272 
273  // convert source from neighbour to local region
274  fvMatrix<vector> porosityEqn(U, eqn.dimensions());
275  scalarField& Udiag = porosityEqn.diag();
276  vectorField& Usource = porosityEqn.source();
277 
278  Udiag.setSize(eqn.diag().size(), 0.0);
279  Usource.setSize(eqn.source().size(), vector::zero);
280 
281  meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
282  meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
283 
284  eqn -= porosityEqn;
285 }
286 
287 
289 {
291  {
292  coeffs_.readIfPresent("UName", UName_);
293  coeffs_.readIfPresent("muName", muName_);
294 
295  // reset the porosity model?
296 
297  return true;
298  }
299  else
300  {
301  return false;
302  }
303 }
304 
305 
306 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fv::interRegionOption::nbrRegionName_
word nbrRegionName_
Name of the neighbour region to map.
Definition: interRegionOption.H:63
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::dimDensity
const dimensionSet dimDensity
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:286
Foam::fv::interRegionExplicitPorositySource::firstIter_
bool firstIter_
First iteration.
Definition: interRegionExplicitPorositySource.H:88
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:72
porosityModel.H
Foam::fv::interRegionExplicitPorositySource::addSup
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Vector.
Definition: interRegionExplicitPorositySource.C:135
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:78
Foam::fv::interRegionExplicitPorositySource::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: interRegionExplicitPorositySource.C:288
U
U
Definition: pEqn.H:46
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:61
interRegionExplicitPorositySource.H
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fv::interRegionOption
Base class for inter-region exchange.
Definition: interRegionOption.H:51
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fv::option::coeffs_
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:84
Foam::ZoneMesh< cellZone, polyMesh >
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::fv::interRegionExplicitPorositySource::porosityPtr_
autoPtr< porosityModel > porosityPtr_
Run-time selectable porosity model.
Definition: interRegionExplicitPorositySource.H:85
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::plusEqOp
Definition: ops.H:71
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
rho
rho
Definition: pEqn.H:3
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
Foam::fv::interRegionExplicitPorositySource::interRegionExplicitPorositySource
interRegionExplicitPorositySource(const interRegionExplicitPorositySource &)
Disallow default bitwise copy construct.
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::porosityModel::New
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Selector.
Definition: porosityModelNew.C:31
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:394
Foam::fv::interRegionExplicitPorositySource::initialise
void initialise()
Initialise.
Definition: interRegionExplicitPorositySource.C:50
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fv::interRegionOption::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: interRegionOptionIO.C:30
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023