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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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"
30 #include "mathematicalConstants.H"
31 #include "surfaceFilmRegionModel.H"
32 #include "liquidFilmBase.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class CloudType>
40 :
42  g_(owner.g()),
43  ejectedParcelType_(0),
44  injectionOffset_(1.1),
45  minDiameter_(0),
46  massParcelPatch_(),
47  diameterParcelPatch_(),
48  UFilmPatch_(),
49  rhoFilmPatch_(),
50  deltaFilmPatch_(),
51  nParcelsTransferred_(0),
52  nParcelsInjected_(0),
53  totalMassTransferred_(0)
54 {}
55 
56 
57 template<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_(),
80  diameterParcelPatch_(),
81  UFilmPatch_(),
82  rhoFilmPatch_(),
83  deltaFilmPatch_(owner.mesh().boundary().size()),
84  nParcelsTransferred_(0),
85  nParcelsInjected_(0),
86  totalMassTransferred_()
87 {}
88 
89 
90 template<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 {}
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class CloudType>
115 template<class CloudTrackType>
117 (
118  const label primaryPatchi,
119  const labelUList& injectorCells,
120  CloudTrackType& cloud
121 )
122 {
123  const fvMesh& mesh = this->owner().mesh();
124 
125  const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
126  const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
127 
128  // Cached sizes are identical to the patch size
129  forAll(injectorCells, facei)
130  {
131  const label celli = injectorCells[facei];
132 
133  if (diameterParcelPatch_[facei] > 0)
134  {
135  const scalar offset =
136  (
137  injectionOffset_
138  * max
139  (
140  diameterParcelPatch_[facei],
141  deltaFilmPatch_[primaryPatchi][facei]
142  )
143  );
144 
145  const point pos = Cf[facei] - offset*normalised(Sf[facei]);
146 
147  // Create a new parcel
148  parcelType* pPtr =
149  new parcelType(this->owner().pMesh(), pos, celli);
150 
151  // Check/set new parcel thermo properties
152  cloud.setParcelThermoProperties(*pPtr, 0.0);
153 
154  setParcelProperties(*pPtr, facei);
155 
156  if (pPtr->nParticle() > 0.001)
157  {
158  // Check new parcel properties
159  cloud.checkParcelProperties(*pPtr, 0.0, false);
160 
161  // Add the new parcel to the cloud
162  cloud.addParticle(pPtr);
163 
164  ++nParcelsInjected_;
165  }
166  else
167  {
168  // TODO: cache mass and re-distribute?
169  delete pPtr;
170  }
171  }
172  }
173 }
174 
175 
176 template<class CloudType>
177 template<class CloudTrackType>
179 (
180  const UList<labelPair>& patchFaces,
181  CloudTrackType& cloud
182 )
183 {
184  const fvMesh& mesh = this->owner().mesh();
186 
187  const auto& Cf = mesh.C().boundaryField();
188  const auto& Sf = mesh.Sf().boundaryField();
189 
190  forAll(patchFaces, filmFacei)
191  {
192  const labelPair& patchAndFace = patchFaces[filmFacei];
193  const label patchi = patchAndFace.first();
194  const label facei = patchAndFace.second();
195 
196  if (patchi < 0) continue; // extra safety
197 
198  const label celli = pbm[patchi].faceCells()[facei];
199 
200  if (diameterParcelPatch_[filmFacei] > 0)
201  {
202  const scalar offset =
203  injectionOffset_ * max
204  (
205  diameterParcelPatch_[filmFacei],
206  deltaFilmPatch_[patchi][facei]
207  );
208 
209  const point pos =
210  (
211  Cf[patchAndFace]
212  - offset * normalised(Sf[patchAndFace])
213  );
214 
215  // Create a new parcel
216  parcelType* pPtr =
217  new parcelType(this->owner().pMesh(), pos, celli);
218 
219  // Check/set new parcel thermo properties
220  cloud.setParcelThermoProperties(*pPtr, 0.0);
221 
222  setParcelProperties(*pPtr, filmFacei);
223 
224  if (pPtr->nParticle() > 0.001)
225  {
226  // Check new parcel properties
227  cloud.checkParcelProperties(*pPtr, 0.0, false);
228 
229  // Add the new parcel to the cloud
230  cloud.addParticle(pPtr);
231 
232  ++nParcelsInjected_;
233  }
234  else
235  {
236  // TODO: cache mass and re-distribute?
237  delete pPtr;
238  }
239  }
240  }
241 }
242 
243 
244 template<class CloudType>
245 template<class TrackCloudType>
247 {
248  if (!this->active())
249  {
250  return;
251  }
252 
253  const fvMesh& mesh = this->owner().mesh();
255 
256  // Check the singleLayer type of films
257  {
258  const auto* filmPtr =
259  mesh.time().objectRegistry::template findObject<regionFilm>
260  (
261  "surfaceFilmProperties"
262  );
263 
264  if (filmPtr && filmPtr->active())
265  {
266  const auto& film = *filmPtr;
267  const labelList& filmPatches = film.intCoupledPatchIDs();
268  const labelList& primaryPatches = film.primaryPatchIDs();
269 
270  forAll(filmPatches, i)
271  {
272  const label filmPatchi = filmPatches[i];
273  const label primaryPatchi = primaryPatches[i];
274 
275  cacheFilmFields(filmPatchi, primaryPatchi, film);
276 
277  injectParticles
278  (
279  primaryPatchi,
280  pbm[primaryPatchi].faceCells(), // injector cells
281  cloud
282  );
283  }
284  }
285  }
286 
287  // Check finite area films
288  for
289  (
290  const areaFilm& regionFa
291  : mesh.time().objectRegistry::template csorted<areaFilm>()
292  )
293  {
294  if (regionFa.active())
295  {
296  auto& film = const_cast<areaFilm&>(regionFa);
297 
298  const List<labelPair>& patchFaces =
299  film.regionMesh().whichPatchFaces();
300 
301 
302  cacheFilmFields(film);
303 
304  injectParticles(patchFaces, cloud);
305 
306  forAll(patchFaces, filmFacei)
307  {
308  const label patchi = patchFaces[filmFacei].first();
309  const label facei = patchFaces[filmFacei].second();
310 
311  if (diameterParcelPatch_[filmFacei] > 0)
312  {
313  film.addSources
314  (
315  patchi,
316  facei,
317  -massParcelPatch_[filmFacei],// mass
318  Zero, // tangential momentum
319  Zero, // impingement
320  Zero // energy
321  );
322  }
323  }
324  }
325  }
326 }
327 
328 
329 template<class CloudType>
331 (
332  const label filmPatchi,
333  const label primaryPatchi,
335 )
336 {
337  massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
338  filmModel.toPrimary(filmPatchi, massParcelPatch_);
339 
340  diameterParcelPatch_ =
341  filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
342  filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
343 
344  UFilmPatch_ = filmModel.Us().boundaryField()[filmPatchi];
345  filmModel.toPrimary(filmPatchi, UFilmPatch_);
346 
347  rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
348  filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
349 
350  deltaFilmPatch_[primaryPatchi] =
351  filmModel.delta().boundaryField()[filmPatchi];
352  filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
353 }
354 
355 
356 template<class CloudType>
358 (
360 )
361 {
362  const polyBoundaryMesh& pbm = this->owner().mesh().boundaryMesh();
363  const volSurfaceMapping& map = film.region().vsm();
364 
365  // The polyPatch/local-face for each of the faceLabels
366  const List<labelPair>& patchFaces =
367  film.regionMesh().whichPatchFaces();
368 
369  const label nFaces = film.Uf().size(); // or regionMesh().nFaces()
370 
371 
372  // Flat fields
373 
374  massParcelPatch_.resize(nFaces, Zero);
375  map.mapToSurface(film.cloudMassTrans(), massParcelPatch_);
376 
377  diameterParcelPatch_.resize(nFaces, Zero);
378  map.mapToSurface(film.cloudDiameterTrans(), diameterParcelPatch_);
379 
380  // Direct copy (one-to-one mapping)
381  UFilmPatch_ = film.Uf().primitiveField();
382 
383  // Direct copy (one-to-one mapping)
384  rhoFilmPatch_ = film.rho().primitiveField();
385 
386 
387  // Per-patch fields
388 
389  // Same as film.region().primaryPatchIDs()
390  for (const label patchi : film.regionMesh().whichPolyPatches())
391  {
392  deltaFilmPatch_[patchi].resize(pbm[patchi].size(), Zero);
393  }
394 
395  forAll(patchFaces, i)
396  {
397  const label patchi = patchFaces[i].first();
398  const label facei = patchFaces[i].second();
399 
400  if (patchi >= 0)
401  {
402  deltaFilmPatch_[patchi][facei] = film.h()[i];
403  }
404  }
405 }
406 
407 
408 template<class CloudType>
410 (
411  parcelType& p,
412  const label filmFacei
413 ) const
414 {
415  // Set parcel properties
416  scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
417  p.d() = diameterParcelPatch_[filmFacei];
418  p.U() = UFilmPatch_[filmFacei];
419  p.rho() = rhoFilmPatch_[filmFacei];
420 
421  p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
422 
423  if (minDiameter_ != -1)
424  {
425  if (p.d() < minDiameter_)
426  {
427  p.nParticle() = 0;
428  }
429  }
430 
431  if (ejectedParcelType_ >= 0)
432  {
433  p.typeId() = ejectedParcelType_;
434  }
435 }
436 
437 
438 template<class CloudType>
440 {
442 
443  label nTrans0 =
444  this->template getModelProperty<label>("nParcelsTransferred");
445 
446  label nInject0 =
447  this->template getModelProperty<label>("nParcelsInjected");
448 
449  scalar massTransferred0 =
450  this->template getModelProperty<scalar>("massTransferred");
451 
452  label nTransTotal =
453  nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
454 
455  label nInjectTotal =
456  nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
457 
458  scalar massTransferredTotal =
459  massTransferred0 + returnReduce(totalMassTransferred_, sumOp<scalar>());
460 
461 
462  Log_<< " Surface film:" << nl
463  << " - parcels absorbed = " << nTransTotal << nl
464  << " - mass absorbed = " << massTransferredTotal << nl
465  << " - parcels ejected = " << nInjectTotal << endl;
466 
467  if (this->writeTime())
468  {
469  this->setModelProperty("nParcelsTransferred", nTransTotal);
470  this->setModelProperty("nParcelsInjected", nInjectTotal);
471  this->setModelProperty("massTransferred", massTransferredTotal);
472 
473  nParcelsTransferred_ = 0;
474  nParcelsInjected_ = 0;
475  totalMassTransferred_ = 0;
476  }
477 }
478 
479 
480 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
481 
482 #include "SurfaceFilmModelNew.C"
483 
484 // ************************************************************************* //
Different types of constants.
faceListList boundary
const polyBoundaryMesh & pbm
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
const T & first() const noexcept
Access the first element.
Definition: Pair.H:137
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const surfaceVectorField & Sf() const
Return cell face area vectors.
void injectParticles(const label primaryPatchi, const labelUList &injectorCells, TrackCloudType &cloud)
Inject particles in cloud.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
const areaScalarField & h() const
Access const reference h.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition: faMeshI.H:135
Base class for cloud sub-models.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
virtual const volScalarField & cloudMassTrans() const =0
Return mass transfer source - Eulerian phase only.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Volume to surface and surface to volume mapping.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dimensionedScalar pos(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: POSIX.C:799
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
dynamicFvMesh & mesh
#define Log_
Report write to Foam::Info if the class log switch is true.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
const areaVectorField & Uf() const
Access const reference Uf.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const faMesh & regionMesh() const
Return the region mesh database.
constexpr scalar pi(M_PI)
void inject(TrackCloudType &cloud)
Inject parcels into the cloud.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual void info()
Write surface film info.
Templated wall surface film model class.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
const labelList & primaryPatchIDs() const
List of patch IDs on the primary region coupled to this region.
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionFilm &)
Cache the film fields in preparation for injection.
Foam::KinematicCloud< Cloud< basicKinematicCollidingParcel > > ::parcelType parcelType
Convenience typedef to the cloud&#39;s parcel type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const areaScalarField & rho() const =0
Access const reference rho.
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
Definition: faMeshI.H:145
const volVectorField & C() const
Return cell centres as volVectorField.
volScalarField & p
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
const regionFaModel & region() const
Access to this region.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film to cloud.
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
Add sources.