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  (
89  "void Foam::fv::interRegionExplicitPorositySource::initialise()"
90  )
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 
113 Foam::fv::interRegionExplicitPorositySource::interRegionExplicitPorositySource
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_(NULL),
123  firstIter_(-1),
124  UName_(coeffs_.lookupOrDefault<word>("UName", "U")),
125  muName_(coeffs_.lookupOrDefault<word>("muName", "thermo:mu"))
126 {
127  if (active_)
128  {
129  fieldNames_.setSize(1, UName_);
130  applied_.setSize(1, false);
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 (
139  fvMatrix<vector>& eqn,
140  const label fieldI
141 )
142 {
143  initialise();
144 
145  const fvMesh& 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("zero", U.dimensions(), vector::zero)
161  );
162 
163  // map local velocity onto neighbour region
164  meshInterp().mapSrcToTgt
165  (
166  U.internalField(),
168  UNbr.internalField()
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(), vector::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 fvMesh& 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("zero", U.dimensions(), vector::zero)
215  );
216 
217  // map local velocity onto neighbour region
218  meshInterp().mapSrcToTgt
219  (
220  U.internalField(),
222  UNbr.internalField()
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,
238  dimensionedScalar("zero", dimDensity, 0.0)
239  );
240 
241  volScalarField muNbr
242  (
243  IOobject
244  (
245  "mu:UNbr",
246  nbrMesh.time().timeName(),
247  nbrMesh,
250  ),
251  nbrMesh,
252  dimensionedScalar("zero", dimViscosity, 0.0)
253  );
254 
255  const volScalarField& mu =
256  mesh_.lookupObject<volScalarField>(muName_);
257 
258  // map local rho onto neighbour region
259  meshInterp().mapSrcToTgt
260  (
261  rho.internalField(),
263  rhoNbr.internalField()
264  );
265 
266  // map local mu onto neighbour region
267  meshInterp().mapSrcToTgt
268  (
269  mu.internalField(),
271  muNbr.internalField()
272  );
273 
274  porosityPtr_->addResistance(nbrEqn, rhoNbr, muNbr);
275 
276  // convert source from neighbour to local region
277  fvMatrix<vector> porosityEqn(U, eqn.dimensions());
278  scalarField& Udiag = porosityEqn.diag();
279  vectorField& Usource = porosityEqn.source();
280 
281  Udiag.setSize(eqn.diag().size(), 0.0);
282  Usource.setSize(eqn.source().size(), vector::zero);
283 
284  meshInterp().mapTgtToSrc(nbrEqn.diag(), plusEqOp<scalar>(), Udiag);
285  meshInterp().mapTgtToSrc(nbrEqn.source(), plusEqOp<vector>(), Usource);
286 
287  eqn -= porosityEqn;
288 }
289 
290 
292 {
293  if (interRegionOption::read(dict))
294  {
295  coeffs_.readIfPresent("UName", UName_);
296  coeffs_.readIfPresent("muName", muName_);
297 
298  // reset the porosity model?
299 
300  return true;
301  }
302  else
303  {
304  return false;
305  }
306 }
307 
308 
309 // ************************************************************************* //
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
defineTypeNameAndDebug(cellSetOption, 0)
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet & dimensions() const
Definition: fvMatrix.H:286
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
bool set(const label) const
Is element set.
Definition: PtrListI.H:107
virtual bool read(const dictionary &dict)
Read dictionary.
A subset of mesh cells.
Definition: cellZone.H:61
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Selector.
A class for handling words, derived from string.
Definition: word.H:59
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 word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
InternalField & internalField()
Return internal field.
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:68
const dimensionSet dimDensity
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Base class for inter-region exchange.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Vector.
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:400
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
void setSize(const label)
Reset size of List.
Definition: List.C:318
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
const dimensionSet dimViscosity
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.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalarField & diag()
Definition: lduMatrix.C:183
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
labelList fv(nPoints)
static const Vector zero
Definition: Vector.H:80
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:354
Field< Type > & source()
Definition: fvMatrix.H:291
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1079
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read(const dictionary &dict)
Read dictionary.
U
Definition: pEqn.H:82
A special matrix type and solver, designed for finite volume solutions of scalar equations.