SurfaceFilmModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 
26 #include "SurfaceFilmModel.H"
27 #include "surfaceFilmRegionModel.H"
28 #include "mathematicalConstants.H"
29 
30 using namespace Foam::constant;
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 :
38  g_(owner.g()),
39  ejectedParcelType_(0),
40  massParcelPatch_(0),
41  diameterParcelPatch_(0),
42  UFilmPatch_(0),
43  rhoFilmPatch_(0),
44  deltaFilmPatch_(0),
45  nParcelsTransferred_(0),
46  nParcelsInjected_(0)
47 {}
48 
49 
50 template<class CloudType>
52 (
53  const dictionary& dict,
54  CloudType& owner,
55  const word& type
56 )
57 :
58  CloudSubModelBase<CloudType>(owner, dict, typeName, type),
59  g_(owner.g()),
60  ejectedParcelType_
61  (
62  this->coeffDict().lookupOrDefault("ejectedParcelType", -1)
63  ),
64  massParcelPatch_(0),
65  diameterParcelPatch_(0),
66  UFilmPatch_(0),
67  rhoFilmPatch_(0),
68  deltaFilmPatch_(owner.mesh().boundary().size()),
69  nParcelsTransferred_(0),
70  nParcelsInjected_(0)
71 {}
72 
73 
74 template<class CloudType>
76 (
78 )
79 :
81  g_(sfm.g_),
82  ejectedParcelType_(sfm.ejectedParcelType_),
83  massParcelPatch_(sfm.massParcelPatch_),
84  diameterParcelPatch_(sfm.diameterParcelPatch_),
85  UFilmPatch_(sfm.UFilmPatch_),
86  rhoFilmPatch_(sfm.rhoFilmPatch_),
87  deltaFilmPatch_(sfm.deltaFilmPatch_),
88  nParcelsTransferred_(sfm.nParcelsTransferred_),
89  nParcelsInjected_(sfm.nParcelsInjected_)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
95 template<class CloudType>
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 template<class CloudType>
103 template<class TrackCloudType>
105 {
106  if
107  (
108  !this->owner().mesh().time().objectRegistry::template foundObject
110  (
111  "surfaceFilmProperties"
112  )
113  )
114  {
115  return;
116  }
117 
118  // Retrieve the film model from the owner database
120  this->owner().mesh().time().objectRegistry::template lookupObject
122  (
123  "surfaceFilmProperties"
124  );
125 
126  const labelList& filmPatches = filmModel.intCoupledPatchIDs();
127  const labelList& primaryPatches = filmModel.primaryPatchIDs();
128 
129  const fvMesh& mesh = this->owner().mesh();
130  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
131 
132  forAll(filmPatches, i)
133  {
134  const label filmPatchi = filmPatches[i];
135  const label primaryPatchi = primaryPatches[i];
136 
137  const labelList& injectorCellsPatch = pbm[primaryPatchi].faceCells();
138 
139  cacheFilmFields(filmPatchi, primaryPatchi, filmModel);
140 
141  const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
142  const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
143  const scalarField& magSf = mesh.magSf().boundaryField()[primaryPatchi];
144 
145  forAll(injectorCellsPatch, j)
146  {
147  if (diameterParcelPatch_[j] > 0)
148  {
149  const label celli = injectorCellsPatch[j];
150 
151  const scalar offset =
152  max
153  (
154  diameterParcelPatch_[j],
155  deltaFilmPatch_[primaryPatchi][j]
156  );
157  const point pos = Cf[j] - 1.1*offset*Sf[j]/magSf[j];
158 
159  // Create a new parcel
160  parcelType* pPtr =
161  new parcelType(this->owner().pMesh(), pos, celli);
162 
163  // Check/set new parcel thermo properties
164  cloud.setParcelThermoProperties(*pPtr, 0.0);
165 
166  setParcelProperties(*pPtr, j);
167 
168  if (pPtr->nParticle() > 0.001)
169  {
170  // Check new parcel properties
171  // cloud.checkParcelProperties(*pPtr, 0.0, true);
172  cloud.checkParcelProperties(*pPtr, 0.0, false);
173 
174  // Add the new parcel to the cloud
175  cloud.addParticle(pPtr);
176 
177  nParcelsInjected_++;
178  }
179  else
180  {
181  // TODO: cache mass and re-distribute?
182  delete pPtr;
183  }
184  }
185  }
186  }
187 }
188 
189 
190 template<class CloudType>
192 (
193  const label filmPatchi,
194  const label primaryPatchi,
196 )
197 {
198  massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
199  filmModel.toPrimary(filmPatchi, massParcelPatch_);
200 
201  diameterParcelPatch_ =
202  filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
203  filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
204 
205  UFilmPatch_ = filmModel.U().boundaryField()[filmPatchi];
206  filmModel.toPrimary(filmPatchi, UFilmPatch_);
207 
208  rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
209  filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
210 
211  deltaFilmPatch_[primaryPatchi] =
212  filmModel.delta().boundaryField()[filmPatchi];
213  filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
214 }
215 
216 
217 template<class CloudType>
219 (
220  parcelType& p,
221  const label filmFacei
222 ) const
223 {
224  // Set parcel properties
225  scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
226  p.d() = diameterParcelPatch_[filmFacei];
227  p.U() = UFilmPatch_[filmFacei];
228  p.rho() = rhoFilmPatch_[filmFacei];
229 
230  p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
231 
232  if (ejectedParcelType_ >= 0)
233  {
234  p.typeId() = ejectedParcelType_;
235  }
236 }
237 
238 
239 template<class CloudType>
241 {
242  label nTrans0 =
243  this->template getModelProperty<label>("nParcelsTransferred");
244 
245  label nInject0 =
246  this->template getModelProperty<label>("nParcelsInjected");
247 
248  label nTransTotal =
249  nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
250 
251  label nInjectTotal =
252  nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
253 
254  os << " Parcels absorbed into film = " << nTransTotal << nl
255  << " New film detached parcels = " << nInjectTotal << endl;
256 
257  if (this->writeTime())
258  {
259  this->setModelProperty("nParcelsTransferred", nTransTotal);
260  this->setModelProperty("nParcelsInjected", nInjectTotal);
261  nParcelsTransferred_ = 0;
262  nParcelsInjected_ = 0;
263  }
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #include "SurfaceFilmModelNew.C"
270 
271 // ************************************************************************* //
Collection of constants.
virtual void info(Ostream &os)
Write surface film info to stream.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label nParcelsInjected_
Number of parcels injected from the film model.
Base class for cloud sub-models.
fvMesh & mesh
virtual const volScalarField & rho() const =0
Return the film density [kg/m^3].
dimensionedScalar pos(const dimensionedScalar &ds)
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
scalarListList deltaFilmPatch_
Film height of all film patches / patch face.
A class for handling words, derived from string.
Definition: word.H:59
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
virtual const volScalarField & delta() const =0
Return the film thickness [m].
void inject(TrackCloudType &cloud)
Inject parcels into the cloud.
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
Foam::polyBoundaryMesh.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:166
static const char nl
Definition: Ostream.H:260
List< vector > UFilmPatch_
Film velocity / patch face.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
virtual ~SurfaceFilmModel()
Destructor.
dimensionedScalar pow3(const dimensionedScalar &ds)
Foam::MomentumCloud< CloudType > ::parcelType parcelType
Convenience typedef to the cloud&#39;s parcel type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalarList diameterParcelPatch_
Parcel diameter / patch face.
label nParcelsTransferred_
Number of parcels transferred to the film model.
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmRegionModel &)
Cache the film fields in preparation for injection.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
scalarList rhoFilmPatch_
Film density / patch face.
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
const volVectorField & C() const
Return cell centres.
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:173
const dimensionedVector & g
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
label ejectedParcelType_
Ejected parcel type label - id assigned to identify parcel for.
const dimensionedVector & g_
Gravitational acceleration constant.
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.
scalarList massParcelPatch_
Parcel mass / patch face.