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 -------------------------------------------------------------------------------
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 "surfaceFilmRegionModel.H"
31 #include "mathematicalConstants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class CloudType>
39 :
41  g_(owner.g()),
42  ejectedParcelType_(0),
43  massParcelPatch_(0),
44  diameterParcelPatch_(0),
45  UFilmPatch_(0),
46  rhoFilmPatch_(0),
47  deltaFilmPatch_(0),
48  nParcelsTransferred_(0),
49  nParcelsInjected_(0)
50 {}
51 
52 
53 template<class CloudType>
55 (
56  const dictionary& dict,
57  CloudType& owner,
58  const word& type
59 )
60 :
61  CloudSubModelBase<CloudType>(owner, dict, typeName, type),
62  g_(owner.g()),
63  ejectedParcelType_
64  (
65  this->coeffDict().getOrDefault("ejectedParcelType", -1)
66  ),
67  massParcelPatch_(0),
68  diameterParcelPatch_(0),
69  UFilmPatch_(0),
70  rhoFilmPatch_(0),
71  deltaFilmPatch_(owner.mesh().boundary().size()),
72  nParcelsTransferred_(0),
73  nParcelsInjected_(0)
74 {}
75 
76 
77 template<class CloudType>
79 (
81 )
82 :
84  g_(sfm.g_),
85  ejectedParcelType_(sfm.ejectedParcelType_),
86  massParcelPatch_(sfm.massParcelPatch_),
87  diameterParcelPatch_(sfm.diameterParcelPatch_),
88  UFilmPatch_(sfm.UFilmPatch_),
89  rhoFilmPatch_(sfm.rhoFilmPatch_),
90  deltaFilmPatch_(sfm.deltaFilmPatch_),
91  nParcelsTransferred_(sfm.nParcelsTransferred_),
92  nParcelsInjected_(sfm.nParcelsInjected_)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
98 template<class CloudType>
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class CloudType>
106 template<class TrackCloudType>
108 {
109  if (!this->active())
110  {
111  return;
112  }
113 
114  // Retrieve the film model from the owner database
116  this->owner().mesh().time().objectRegistry::template lookupObject
118  (
119  "surfaceFilmProperties"
120  );
121 
122  if (!filmModel.active())
123  {
124  return;
125  }
126 
127  const labelList& filmPatches = filmModel.intCoupledPatchIDs();
128  const labelList& primaryPatches = filmModel.primaryPatchIDs();
129 
130  const fvMesh& mesh = this->owner().mesh();
131  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
132 
133  forAll(filmPatches, i)
134  {
135  const label filmPatchi = filmPatches[i];
136  const label primaryPatchi = primaryPatches[i];
137 
138  const labelList& injectorCellsPatch = pbm[primaryPatchi].faceCells();
139 
140  cacheFilmFields(filmPatchi, primaryPatchi, filmModel);
141 
142  const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
143  const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
144  const scalarField& magSf = mesh.magSf().boundaryField()[primaryPatchi];
145 
146  forAll(injectorCellsPatch, j)
147  {
148  if (diameterParcelPatch_[j] > 0)
149  {
150  const label celli = injectorCellsPatch[j];
151 
152  const scalar offset =
153  max
154  (
155  diameterParcelPatch_[j],
156  deltaFilmPatch_[primaryPatchi][j]
157  );
158  const point pos = Cf[j] - 1.1*offset*Sf[j]/magSf[j];
159 
160  // Create a new parcel
161  parcelType* pPtr =
162  new parcelType(this->owner().pMesh(), pos, celli);
163 
164  // Check/set new parcel thermo properties
165  cloud.setParcelThermoProperties(*pPtr, 0.0);
166 
167  setParcelProperties(*pPtr, j);
168 
169  if (pPtr->nParticle() > 0.001)
170  {
171  // Check new parcel properties
172  // cloud.checkParcelProperties(*pPtr, 0.0, true);
173  cloud.checkParcelProperties(*pPtr, 0.0, false);
174 
175  // Add the new parcel to the cloud
176  cloud.addParticle(pPtr);
177 
178  nParcelsInjected_++;
179  }
180  else
181  {
182  // TODO: cache mass and re-distribute?
183  delete pPtr;
184  }
185  }
186  }
187  }
188 }
189 
190 
191 template<class CloudType>
193 (
194  const label filmPatchi,
195  const label primaryPatchi,
197 )
198 {
199  massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
200  filmModel.toPrimary(filmPatchi, massParcelPatch_);
201 
202  diameterParcelPatch_ =
203  filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
204  filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
205 
206  UFilmPatch_ = filmModel.Us().boundaryField()[filmPatchi];
207  filmModel.toPrimary(filmPatchi, UFilmPatch_);
208 
209  rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
210  filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
211 
212  deltaFilmPatch_[primaryPatchi] =
213  filmModel.delta().boundaryField()[filmPatchi];
214  filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
215 }
216 
217 
218 template<class CloudType>
220 (
221  parcelType& p,
222  const label filmFacei
223 ) const
224 {
225  // Set parcel properties
226  scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
227  p.d() = diameterParcelPatch_[filmFacei];
228  p.U() = UFilmPatch_[filmFacei];
229  p.rho() = rhoFilmPatch_[filmFacei];
230 
231  p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
232 
233  if (ejectedParcelType_ >= 0)
234  {
235  p.typeId() = ejectedParcelType_;
236  }
237 }
238 
239 
240 template<class CloudType>
242 {
243  label nTrans0 =
244  this->template getModelProperty<label>("nParcelsTransferred");
245 
246  label nInject0 =
247  this->template getModelProperty<label>("nParcelsInjected");
248 
249  label nTransTotal =
250  nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
251 
252  label nInjectTotal =
253  nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
254 
255  os << " Surface film:" << nl
256  << " - parcels absorbed = " << nTransTotal << nl
257  << " - parcels ejected = " << nInjectTotal << endl;
258 
259  if (this->writeTime())
260  {
261  this->setModelProperty("nParcelsTransferred", nTransTotal);
262  this->setModelProperty("nParcelsInjected", nInjectTotal);
263  nParcelsTransferred_ = 0;
264  nParcelsInjected_ = 0;
265  }
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 #include "SurfaceFilmModelNew.C"
272 
273 // ************************************************************************* //
Foam::SurfaceFilmModel::cacheFilmFields
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmRegionModel &)
Cache the film fields in preparation for injection.
Definition: SurfaceFilmModel.C:193
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::regionModels::regionModel::active
Switch active() const
Return the active flag.
Definition: regionModelI.H:46
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::SurfaceFilmModel::g_
const dimensionedVector & g_
Gravitational acceleration constant.
Definition: SurfaceFilmModel.H:79
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::SurfaceFilmModel::deltaFilmPatch_
scalarListList deltaFilmPatch_
Film height of all film patches / patch face.
Definition: SurfaceFilmModel.H:102
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::delta
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Foam::SurfaceFilmModel::massParcelPatch_
scalarList massParcelPatch_
Parcel mass / patch face.
Definition: SurfaceFilmModel.H:90
Foam::CloudSubModelBase
Base class for cloud sub-models.
Definition: CloudSubModelBase.H:51
Foam::SurfaceFilmModel::nParcelsTransferred_
label nParcelsTransferred_
Number of parcels transferred to the film model.
Definition: SurfaceFilmModel.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFilmRegionModel.H
Foam::SurfaceFilmModel::~SurfaceFilmModel
virtual ~SurfaceFilmModel()
Destructor.
Definition: SurfaceFilmModel.C:99
Foam::SurfaceFilmModel::nParcelsInjected_
label nParcelsInjected_
Number of parcels injected from the film model.
Definition: SurfaceFilmModel.H:111
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::cloudDiameterTrans
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::Us
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
Foam::regionModels::regionModel::toPrimary
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
Definition: regionModelTemplates.C:147
Foam::SurfaceFilmModel< Foam::KinematicCloud< Cloud< basicKinematicCollidingParcel > > >::parcelType
Foam::KinematicCloud< Cloud< basicKinematicCollidingParcel > > ::parcelType parcelType
Convenience typedef to the cloud's parcel type.
Definition: SurfaceFilmModel.H:76
Foam::Field< vector >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::cloudMassTrans
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::SurfaceFilmModel
Templated wall surface film model class.
Definition: KinematicCloud.H:92
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::SurfaceFilmModel::inject
void inject(TrackCloudType &cloud)
Inject parcels into the cloud.
Definition: SurfaceFilmModel.C:107
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::SurfaceFilmModel::diameterParcelPatch_
scalarList diameterParcelPatch_
Parcel diameter / patch face.
Definition: SurfaceFilmModel.H:93
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::rho
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam::SurfaceFilmModel::setParcelProperties
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
Definition: SurfaceFilmModel.C:220
Foam::maxEqOp
Definition: ops.H:80
Foam::regionModels::regionModel::primaryPatchIDs
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:168
Foam::SurfaceFilmModel::rhoFilmPatch_
scalarList rhoFilmPatch_
Film density / patch face.
Definition: SurfaceFilmModel.H:99
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:679
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::List< label >
SurfaceFilmModelNew.C
Foam::SurfaceFilmModel::info
virtual void info(Ostream &os)
Write surface film info to stream.
Definition: SurfaceFilmModel.C:241
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::SurfaceFilmModel::ejectedParcelType_
label ejectedParcelType_
Ejected parcel type label - id assigned to identify parcel for.
Definition: SurfaceFilmModel.H:84
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
SurfaceFilmModel.H
Foam::SurfaceFilmModel::UFilmPatch_
List< vector > UFilmPatch_
Film velocity / patch face.
Definition: SurfaceFilmModel.H:96
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::SurfaceFilmModel::SurfaceFilmModel
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
Definition: SurfaceFilmModel.C:38
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:175
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177