Loading...
Searching...
No Matches
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-2025 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// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38// The area film models.
39// - 2506 and earlier: registered on Time.
40// - 2512 and later : registered below faMeshesRegistry
41
42template<class CloudType>
45{
46 if (const objectRegistry* parent = faMesh::registry(mesh); parent)
47 {
48 return
49 (
50 parent->thisDb().csorted
51 <
53 >()
54 );
55 }
56 else
57 {
58 return {};
59 }
60}
61
62template<class CloudType>
65{
66 if (const objectRegistry* parent = faMesh::registry(mesh); parent)
67 {
68 return
69 (
70 const_cast<objectRegistry&>(parent->thisDb()).sorted
71 <
73 >()
74 );
75 }
76 else
77 {
78 return {};
79 }
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
85template<class CloudType>
87:
89 g_(owner.g()),
90 ejectedParcelType_(0),
91 injectionOffset_(1.1),
92 minDiameter_(0),
93 massParcelPatch_(),
94 diameterParcelPatch_(),
95 UFilmPatch_(),
96 rhoFilmPatch_(),
97 deltaFilmPatch_(),
101{}
102
103
104template<class CloudType>
106(
107 const dictionary& dict,
108 CloudType& owner,
109 const word& type
110)
111:
113 g_(owner.g()),
114 ejectedParcelType_
115 (
116 this->coeffDict().template getOrDefault<label>("ejectedParcelType", -1)
117 ),
118 injectionOffset_
119 (
120 this->coeffDict().template getOrDefault<scalar>("injectionOffset", 1.1)
121 ),
122 minDiameter_
123 (
124 this->coeffDict().template getOrDefault<scalar>("minDiameter", -1)
125 ),
126 massParcelPatch_(),
127 diameterParcelPatch_(),
128 UFilmPatch_(),
129 rhoFilmPatch_(),
130 deltaFilmPatch_(owner.mesh().boundary().size()),
134{}
135
136
137template<class CloudType>
139(
141)
142:
144 g_(sfm.g_),
145 ejectedParcelType_(sfm.ejectedParcelType_),
146 injectionOffset_(sfm.injectionOffset_),
147 minDiameter_(sfm.minDiameter_),
148 massParcelPatch_(sfm.massParcelPatch_),
149 diameterParcelPatch_(sfm.diameterParcelPatch_),
150 UFilmPatch_(sfm.UFilmPatch_),
151 rhoFilmPatch_(sfm.rhoFilmPatch_),
152 deltaFilmPatch_(sfm.deltaFilmPatch_),
153 nParcelsTransferred_(sfm.nParcelsTransferred_),
154 nParcelsInjected_(sfm.nParcelsInjected_),
155 totalMassTransferred_(sfm.totalMassTransferred_)
157
158
159// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160
161template<class CloudType>
162template<class CloudTrackType>
164(
165 const label primaryPatchi,
166 const labelUList& injectorCells,
167 CloudTrackType& cloud
168)
169{
170 const fvMesh& mesh = this->owner().mesh();
171
172 const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
173 const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
175 // Cached sizes are identical to the patch size
176 forAll(injectorCells, facei)
177 {
178 const label celli = injectorCells[facei];
179
180 if (diameterParcelPatch_[facei] > 0)
181 {
182 const scalar offset =
183 (
185 * max
186 (
188 deltaFilmPatch_[primaryPatchi][facei]
189 )
190 );
191
192 const point pos = Cf[facei] - offset*normalised(Sf[facei]);
193
194 // Create a new parcel
195 parcelType* pPtr =
196 new parcelType(this->owner().pMesh(), pos, celli);
197
198 // Check/set new parcel thermo properties
199 cloud.setParcelThermoProperties(*pPtr, 0.0);
200
201 setParcelProperties(*pPtr, facei);
202
203 if (pPtr->nParticle() > 0.001)
204 {
205 // Check new parcel properties
206 cloud.checkParcelProperties(*pPtr, 0.0, false);
207
208 // Add the new parcel to the cloud
209 cloud.addParticle(pPtr);
212 }
213 else
214 {
215 // TODO: cache mass and re-distribute?
216 delete pPtr;
217 }
219 }
220}
221
222
223template<class CloudType>
224template<class CloudTrackType>
226(
227 const UList<labelPair>& patchFaces,
228 CloudTrackType& cloud
229)
230{
231 const fvMesh& mesh = this->owner().mesh();
232 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
233
234 const auto& Cf = mesh.C().boundaryField();
235 const auto& Sf = mesh.Sf().boundaryField();
236
237 forAll(patchFaces, filmFacei)
238 {
239 const labelPair& patchAndFace = patchFaces[filmFacei];
240 const label patchi = patchAndFace.first();
241 const label facei = patchAndFace.second();
242
243 if (patchi < 0) continue; // extra safety
244
245 const label celli = pbm[patchi].faceCells()[facei];
246
247 if (diameterParcelPatch_[filmFacei] > 0)
248 {
249 const scalar offset =
251 (
252 diameterParcelPatch_[filmFacei],
253 deltaFilmPatch_[patchi][facei]
254 );
255
256 const point pos =
257 (
258 Cf[patchAndFace]
259 - offset * normalised(Sf[patchAndFace])
260 );
261
262 // Create a new parcel
263 parcelType* pPtr =
264 new parcelType(this->owner().pMesh(), pos, celli);
265
266 // Check/set new parcel thermo properties
267 cloud.setParcelThermoProperties(*pPtr, 0.0);
268
269 setParcelProperties(*pPtr, filmFacei);
270
271 if (pPtr->nParticle() > 0.001)
272 {
273 // Check new parcel properties
274 cloud.checkParcelProperties(*pPtr, 0.0, false);
275
276 // Add the new parcel to the cloud
277 cloud.addParticle(pPtr);
278
280 }
281 else
282 {
283 // TODO: cache mass and re-distribute?
284 delete pPtr;
285 }
287 }
288}
289
290
291template<class CloudType>
292template<class TrackCloudType>
294{
295 if (!this->active())
296 {
297 return;
298 }
299
300 const fvMesh& mesh = this->owner().mesh();
301 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
302
303 // Check the singleLayer type of films
304 {
305 const auto* filmPtr =
306 mesh.time().objectRegistry::template findObject<regionFilm>
307 (
308 "surfaceFilmProperties"
309 );
310
311 if (filmPtr && filmPtr->active())
312 {
313 const auto& film = *filmPtr;
314 const labelList& filmPatches = film.intCoupledPatchIDs();
315 const labelList& primaryPatches = film.primaryPatchIDs();
316
317 forAll(filmPatches, i)
318 {
319 const label filmPatchi = filmPatches[i];
320 const label primaryPatchi = primaryPatches[i];
321
322 cacheFilmFields(filmPatchi, primaryPatchi, film);
323
324 injectParticles
325 (
326 primaryPatchi,
327 pbm[primaryPatchi].faceCells(), // injector cells
328 cloud
329 );
330 }
331 }
333
334 // Check finite-area films
335 // - allow non-const access to the area films
337 {
338 if (film.active())
339 {
340 const List<labelPair>& patchFaces =
341 film.regionMesh().whichPatchFaces();
342
343 cacheFilmFields(film);
344
345 injectParticles(patchFaces, cloud);
346
347 forAll(patchFaces, filmFacei)
348 {
349 const label patchi = patchFaces[filmFacei].first();
350 const label facei = patchFaces[filmFacei].second();
351
352 if (diameterParcelPatch_[filmFacei] > 0)
353 {
354 film.addSources
355 (
356 patchi,
357 facei,
358 -massParcelPatch_[filmFacei],// mass
359 Zero, // tangential momentum
360 Zero, // impingement
361 Zero // energy
362 );
363 }
365 }
366 }
368
369
370template<class CloudType>
372(
373 const label filmPatchi,
374 const label primaryPatchi,
376)
377{
378 massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
379 filmModel.toPrimary(filmPatchi, massParcelPatch_);
380
381 diameterParcelPatch_ =
382 filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
383 filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
384
385 UFilmPatch_ = filmModel.Us().boundaryField()[filmPatchi];
386 filmModel.toPrimary(filmPatchi, UFilmPatch_);
387
388 rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
389 filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
390
391 deltaFilmPatch_[primaryPatchi] =
392 filmModel.delta().boundaryField()[filmPatchi];
393 filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
394}
395
396
397template<class CloudType>
399(
401)
402{
403 const polyBoundaryMesh& pbm = this->owner().mesh().boundaryMesh();
404 const volSurfaceMapping& map = film.region().vsm();
405
406 // The polyPatch/local-face for each of the faceLabels
407 const List<labelPair>& patchFaces =
409
410 const label nFaces = film.Uf().size(); // or regionMesh().nFaces()
411
412
413 // Flat fields
414
415 massParcelPatch_.resize(nFaces, Zero);
416 map.mapToSurface(film.cloudMassTrans(), massParcelPatch_);
417
418 diameterParcelPatch_.resize(nFaces, Zero);
419 map.mapToSurface(film.cloudDiameterTrans(), diameterParcelPatch_);
420
421 // Direct copy (one-to-one mapping)
422 UFilmPatch_ = film.Uf().primitiveField();
423
424 // Direct copy (one-to-one mapping)
425 rhoFilmPatch_ = film.rho().primitiveField();
426
427
428 // Per-patch fields
429
430 // Same as film.region().primaryPatchIDs()
431 for (const label patchi : film.regionMesh().whichPolyPatches())
432 {
433 deltaFilmPatch_[patchi].resize(pbm[patchi].size(), Zero);
434 }
435
436 forAll(patchFaces, i)
437 {
438 const label patchi = patchFaces[i].first();
439 const label facei = patchFaces[i].second();
440
441 if (patchi >= 0)
442 {
443 deltaFilmPatch_[patchi][facei] = film.h()[i];
444 }
445 }
446}
447
448
449template<class CloudType>
451(
452 parcelType& p,
453 const label filmFacei
454) const
455{
456 // Set parcel properties
457 scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
458 p.d() = diameterParcelPatch_[filmFacei];
459 p.U() = UFilmPatch_[filmFacei];
460 p.rho() = rhoFilmPatch_[filmFacei];
461
462 p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
463
464 if (minDiameter_ != -1)
465 {
466 if (p.d() < minDiameter_)
467 {
468 p.nParticle() = 0;
469 }
470 }
471
472 if (ejectedParcelType_ >= 0)
474 p.typeId() = ejectedParcelType_;
475 }
476}
477
478
479template<class CloudType>
481{
483
484 label nTrans0 =
485 this->template getModelProperty<label>("nParcelsTransferred");
486
487 label nInject0 =
488 this->template getModelProperty<label>("nParcelsInjected");
489
490 scalar massTransferred0 =
491 this->template getModelProperty<scalar>("massTransferred");
492
493 label nTransTotal =
494 nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
495
496 label nInjectTotal =
497 nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
498
499 scalar massTransferredTotal =
500 massTransferred0 + returnReduce(totalMassTransferred_, sumOp<scalar>());
501
502
503 Log_<< " Surface film:" << nl
504 << " - parcels absorbed = " << nTransTotal << nl
505 << " - mass absorbed = " << massTransferredTotal << nl
506 << " - parcels ejected = " << nInjectTotal << endl;
507
508 if (this->writeTime())
509 {
510 this->setModelProperty("nParcelsTransferred", nTransTotal);
511 this->setModelProperty("nParcelsInjected", nInjectTotal);
512 this->setModelProperty("massTransferred", massTransferredTotal);
513
514 nParcelsTransferred_ = 0;
515 nParcelsInjected_ = 0;
516 totalMassTransferred_ = 0;
517 }
518}
519
520
521// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
522
523#include "SurfaceFilmModelNew.C"
524
525// ************************************************************************* //
const polyBoundaryMesh & pbm
const uniformDimensionedVectorField & g
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
virtual void info()
Write to info.
CloudSubModelBase(CloudType &owner)
Construct null from owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
Templated wall surface film model class.
static UPtrList< const areaFilm > csorted_areaFilms(const polyMesh &)
Registry.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Field< vector > UFilmPatch_
Film velocity / patch face.
scalar minDiameter_
Minimum diameter particle injection.
scalarField rhoFilmPatch_
Film density / patch face.
scalarField massParcelPatch_
Parcel mass / patch face.
void inject(TrackCloudType &cloud)
Inject parcels into the cloud.
label nParcelsInjected_
Number of parcels injected from the film model.
const dimensionedVector & g_
Gravitational acceleration constant.
Field< scalarField > deltaFilmPatch_
Film height of all film patches / patch face.
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionFilm &)
Cache the film fields in preparation for injection.
label ejectedParcelType_
Ejected parcel type label - id assigned to identify parcel for post-processing. If not specified,...
scalarField diameterParcelPatch_
Parcel diameter / patch face.
void injectParticles(const label primaryPatchi, const labelUList &injectorCells, CloudTrackType &cloud)
static UPtrList< areaFilm > sorted_areaFilms(const polyMesh &)
Return a sorted list of area-film objects that are registered on the faMeshesRegistry.
virtual void info()
Write surface film info.
scalar injectionOffset_
Injection offset position.
label nParcelsTransferred_
Number of parcels transferred to the film model.
const dimensionedVector & g() const noexcept
Return gravitational acceleration constant.
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
Foam::KinematicCloud< Foam::DSMCCloud< dsmcParcel > >::parcelType parcelType
void injectParticles(const label primaryPatchi, const labelUList &injectorCells, TrackCloudType &cloud)
Inject particles in cloud.
scalar totalMassTransferred_
Total mass transferred to the film.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels().
Definition faMeshI.H:138
static const objectRegistry * registry(const polyMesh &pMesh)
Find the singleton parent registry (on the polyMesh) that contains all objects related to finite-area...
Definition faMesh.C:145
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition faMeshI.H:128
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Registry of regIOobjects.
UPtrList< Type > sorted()
Return sorted list of objects with a class satisfying isA<Type> or isType<Type> (with Strict).
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const regionFaModel & region() const noexcept
Access to this region.
virtual const areaScalarField & rho() const =0
Access const reference rho.
const areaVectorField & Uf() const noexcept
Access const reference Uf.
const areaScalarField & h() const noexcept
Access const reference h.
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.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const faMesh & regionMesh() const
Return the region mesh database.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
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.
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
virtual bool active() const
Return the model 'active' status - default active = true.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
Volume to surface and surface to volume mapping.
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
faceListList boundary
dynamicFvMesh & mesh
#define Log_
Report write to Foam::Info if the class log switch is true.
constexpr scalar pi(M_PI)
Different types of constants.
dimensionedScalar pos(const dimensionedScalar &ds)
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensionedScalar pow3(const dimensionedScalar &ds)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299