SurfaceFilmModel.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 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
29#include "SurfaceFilmModel.H"
32#include "liquidFilmBase.H"
33
34using namespace Foam::constant;
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
38template<class CloudType>
40:
42 g_(owner.g()),
43 ejectedParcelType_(0),
44 injectionOffset_(1.1),
45 minDiameter_(0),
46 massParcelPatch_(0),
47 diameterParcelPatch_(0),
48 UFilmPatch_(0),
49 rhoFilmPatch_(0),
50 deltaFilmPatch_(0),
51 nParcelsTransferred_(0),
52 nParcelsInjected_(0),
53 totalMassTransferred_(0)
54{}
55
56
57template<class CloudType>
59(
60 const dictionary& dict,
61 CloudType& owner,
62 const word& type
63)
64:
65 CloudSubModelBase<CloudType>(owner, dict, typeName, type),
66 g_(owner.g()),
67 ejectedParcelType_
68 (
69 this->coeffDict().template getOrDefault<label>("ejectedParcelType", -1)
70 ),
71 injectionOffset_
72 (
73 this->coeffDict().template getOrDefault<scalar>("injectionOffset", 1.1)
74 ),
75 minDiameter_
76 (
77 this->coeffDict().template getOrDefault<scalar>("minDiameter", -1)
78 ),
79 massParcelPatch_(0),
80 diameterParcelPatch_(0),
81 UFilmPatch_(0),
82 rhoFilmPatch_(0),
83 deltaFilmPatch_(owner.mesh().boundary().size()),
84 nParcelsTransferred_(0),
85 nParcelsInjected_(0),
86 totalMassTransferred_()
87{}
88
89
90template<class CloudType>
92(
94)
95:
97 g_(sfm.g_),
98 ejectedParcelType_(sfm.ejectedParcelType_),
99 injectionOffset_(sfm.injectionOffset_),
100 minDiameter_(sfm.minDiameter_),
101 massParcelPatch_(sfm.massParcelPatch_),
102 diameterParcelPatch_(sfm.diameterParcelPatch_),
103 UFilmPatch_(sfm.UFilmPatch_),
104 rhoFilmPatch_(sfm.rhoFilmPatch_),
105 deltaFilmPatch_(sfm.deltaFilmPatch_),
106 nParcelsTransferred_(sfm.nParcelsTransferred_),
107 nParcelsInjected_(sfm.nParcelsInjected_),
108 totalMassTransferred_(sfm.totalMassTransferred_)
109{}
110
111
112// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113
114template<class CloudType>
116{}
117
118
119// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120
121template<class CloudType>
122template<class CloudTrackType>
124(
125 const label primaryPatchi,
126 const labelList& injectorCellsPatch,
127 CloudTrackType& cloud
128)
129{
130 const fvMesh& mesh = this->owner().mesh();
131 const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
132 const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
133 const scalarField& magSf =
134 mesh.magSf().boundaryField()[primaryPatchi];
135
136 forAll(injectorCellsPatch, j)
137 {
138 if (diameterParcelPatch_[j] > 0)
139 {
140 const label celli = injectorCellsPatch[j];
141
142 const scalar offset =
144 (
145 diameterParcelPatch_[j],
146 deltaFilmPatch_[primaryPatchi][j]
147 );
148 const point pos = Cf[j] - injectionOffset_*offset*Sf[j]/magSf[j];
149
150 // Create a new parcel
152 new parcelType(this->owner().pMesh(), pos, celli);
153
154 // Check/set new parcel thermo properties
155 cloud.setParcelThermoProperties(*pPtr, 0.0);
156
157 setParcelProperties(*pPtr, j);
158
159 if (pPtr->nParticle() > 0.001)
160 {
161 // Check new parcel properties
162 cloud.checkParcelProperties(*pPtr, 0.0, false);
163
164 // Add the new parcel to the cloud
165 cloud.addParticle(pPtr);
166
167 nParcelsInjected_++;
168 }
169 else
170 {
171 // TODO: cache mass and re-distribute?
172 delete pPtr;
173 }
174 }
175 }
176}
177
178template<class CloudType>
179template<class TrackCloudType>
181{
182 if (!this->active())
183 {
184 return;
185 }
186
187 const fvMesh& mesh = this->owner().mesh();
188 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
189
190 // Retrieve the film model from the owner database
192 mesh.time().objectRegistry::template findObject
194 (
195 "surfaceFilmProperties"
196 );
197
198 // Check the singleLayer type of films
199 if (filmModel && filmModel->active())
200 {
201
202 const labelList& filmPatches = filmModel->intCoupledPatchIDs();
203 const labelList& primaryPatches = filmModel->primaryPatchIDs();
204
205 forAll(filmPatches, i)
206 {
207 const label filmPatchi = filmPatches[i];
208 const label primaryPatchi = primaryPatches[i];
209
210 const labelList& injectorCellsPatch = pbm[primaryPatchi].faceCells();
211
212 cacheFilmFields(filmPatchi, primaryPatchi, *filmModel);
213
214 injectParticles(primaryPatchi, injectorCellsPatch, cloud);
215 }
216 }
217
218 // Check finite area films
219 wordList names =
220 mesh.time().objectRegistry::template
221 sortedNames<regionModels::regionFaModel>();
222
223 forAll (names, i)
224 {
225 const regionModels::regionFaModel* regionFa =
226 mesh.time().objectRegistry::template cfindObject
227 <
229 >(names[i]);
230
231 // Check that it is a type areaFilm
232 if (regionFa && isA<areaFilm>(*regionFa) && regionFa->active())
233 {
234 areaFilm& film =
235 const_cast<areaFilm&>(refCast<const areaFilm>(*regionFa));
236
237 const label patchId = regionFa->patchID();
238
239 const labelList& injectorCellsPatch = pbm[patchId].faceCells();
240
241 cacheFilmFields(patchId, film);
242
243 injectParticles(patchId, injectorCellsPatch, cloud);
244
245 forAll(injectorCellsPatch, facei)
246 {
247 if (diameterParcelPatch_[facei] > 0)
248 {
249 film.addSources
250 (
251 patchId,
252 facei,
253 -massParcelPatch_[facei], // mass
254 Zero, // tangential momentum
255 Zero, // impingement
256 Zero // energy
257 );
258 }
259 }
260 }
261 }
262}
263
264
265template<class CloudType>
267(
268 const label filmPatchi,
269 const label primaryPatchi,
271)
272{
273 massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
274 filmModel.toPrimary(filmPatchi, massParcelPatch_);
275
276 diameterParcelPatch_ =
277 filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
278 filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
279
280 UFilmPatch_ = filmModel.Us().boundaryField()[filmPatchi];
281 filmModel.toPrimary(filmPatchi, UFilmPatch_);
282
283 rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
284 filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
285
286 deltaFilmPatch_[primaryPatchi] =
287 filmModel.delta().boundaryField()[filmPatchi];
288 filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
289}
290
291
292template<class CloudType>
294(
295 const label filmPatchi,
297)
298{
299 const volSurfaceMapping& map = filmModel.region().vsm();
300
301 massParcelPatch_.setSize(filmModel.Uf().size(), Zero);
302
303 const scalarField& massParcelPatch =
304 filmModel.cloudMassTrans().boundaryField()[filmPatchi];
305
306 map.mapToField(massParcelPatch, massParcelPatch_);
307
308 const scalarField& diameterParcelPatch =
309 filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
310
311 diameterParcelPatch_.setSize(filmModel.Uf().size(), Zero);
312 map.mapToField(diameterParcelPatch, diameterParcelPatch_);
313
314 UFilmPatch_.setSize(filmModel.Uf().size(), vector::zero);
315 map.mapToField(filmModel.Uf(), UFilmPatch_);
316
317 rhoFilmPatch_.setSize(UFilmPatch_.size(), Zero);
318 map.mapToField(filmModel.rho(), rhoFilmPatch_);
319
320 deltaFilmPatch_[filmPatchi].setSize(UFilmPatch_.size(), Zero);
321 map.mapToField(filmModel.h(), deltaFilmPatch_[filmPatchi]);
322}
323
324
325template<class CloudType>
327(
328 parcelType& p,
329 const label filmFacei
330) const
331{
332 // Set parcel properties
333 scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
334 p.d() = diameterParcelPatch_[filmFacei];
335 p.U() = UFilmPatch_[filmFacei];
336 p.rho() = rhoFilmPatch_[filmFacei];
337
338 p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
339
340 if (minDiameter_ != -1)
341 {
342 if (p.d() < minDiameter_)
343 {
344 p.nParticle() = 0;
345 }
346 }
347
348 if (ejectedParcelType_ >= 0)
349 {
350 p.typeId() = ejectedParcelType_;
351 }
352}
353
354
355template<class CloudType>
357{
358 label nTrans0 =
359 this->template getModelProperty<label>("nParcelsTransferred");
360
361 label nInject0 =
362 this->template getModelProperty<label>("nParcelsInjected");
363
364 scalar massTransferred0 =
365 this->template getModelProperty<scalar>("massTransferred");
366
367 label nTransTotal =
368 nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
369
370 label nInjectTotal =
371 nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
372
373 scalar massTransferredTotal =
374 massTransferred0 + returnReduce(totalMassTransferred_, sumOp<scalar>());
375
376
377 os << " Surface film:" << nl
378 << " - parcels absorbed = " << nTransTotal << nl
379 << " - mass absorbed = " << massTransferredTotal << nl
380 << " - parcels ejected = " << nInjectTotal << endl;
381
382 if (this->writeTime())
383 {
384 this->setModelProperty("nParcelsTransferred", nTransTotal);
385 this->setModelProperty("nParcelsInjected", nInjectTotal);
386 this->setModelProperty("massTransferred", massTransferredTotal);
387
388 nParcelsTransferred_ = 0;
389 nParcelsInjected_ = 0;
390 totalMassTransferred_ = 0;
391 }
392}
393
394
395// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396
397#include "SurfaceFilmModelNew.C"
398
399// ************************************************************************* //
const uniformDimensionedVectorField & g
Base class for cloud sub-models.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
Templated wall surface film model class.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmRegionModel &)
Cache the film fields in preparation for injection.
void inject(TrackCloudType &cloud)
Inject parcels into the cloud.
void injectParticles(const label primaryPatchi, const labelList &injectorCellsPatch, TrackCloudType &cloud)
Inject particles in cloud.
virtual ~SurfaceFilmModel()
Destructor.
CloudType::parcelType parcelType
Convenience typedef to the cloud's parcel type.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const regionFaModel & region() const
Access to this region.
const areaVectorField & Uf() const
Access const reference Uf.
const areaScalarField & h() const
Access const reference h.
virtual const areaScalarField & rho() const =0
Access const reference rho.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film to cloud.
virtual const volScalarField & cloudMassTrans() const =0
Return mass transfer source - Eulerian phase only.
Base class for area region models.
label patchID() const
Return patch ID.
const Switch & active() const
Return the active flag.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:168
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:175
Switch active() const
Return the active flag.
Definition: regionModelI.H:46
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
Volume to surface and surface to volume mapping.
void mapToField(const GeometricField< Type, faPatchField, areaMesh > &af, Field< Type > &f) const
Map surface field to field.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
faceListList boundary
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
label patchId(-1)
constexpr scalar pi(M_PI)
Different types of constants.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pow3(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59