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 | Copyright (C) 2011-2016 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 "surfaceFilmModel.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 :
59  g_(owner.g()),
61  (
62  this->coeffDict().lookupOrDefault("ejectedParcelType", -1)
63  ),
66  UFilmPatch_(0),
67  rhoFilmPatch_(0),
68  deltaFilmPatch_(owner.mesh().boundary().size()),
71 {}
72 
73 
74 template<class CloudType>
76 (
78 )
79 :
81  g_(sfm.g_),
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
95 template<class CloudType>
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 template<class CloudType>
103 template<class TrackData>
105 {
106  if (!this->active())
107  {
108  return;
109  }
110 
111  // Retrieve the film model from the owner database
113  this->owner().mesh().time().objectRegistry::template lookupObject
115  (
116  "surfaceFilmProperties"
117  );
118 
119  if (!filmModel.active())
120  {
121  return;
122  }
123 
124  const labelList& filmPatches = filmModel.intCoupledPatchIDs();
125  const labelList& primaryPatches = filmModel.primaryPatchIDs();
126 
127  const fvMesh& mesh = this->owner().mesh();
128  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
129 
130  forAll(filmPatches, i)
131  {
132  const label filmPatchi = filmPatches[i];
133  const label primaryPatchi = primaryPatches[i];
134 
135  const labelList& injectorCellsPatch = pbm[primaryPatchi].faceCells();
136 
137  cacheFilmFields(filmPatchi, primaryPatchi, filmModel);
138 
139  const vectorField& Cf = mesh.C().boundaryField()[primaryPatchi];
140  const vectorField& Sf = mesh.Sf().boundaryField()[primaryPatchi];
141  const scalarField& magSf = mesh.magSf().boundaryField()[primaryPatchi];
142 
143  forAll(injectorCellsPatch, j)
144  {
145  if (diameterParcelPatch_[j] > 0)
146  {
147  const label celli = injectorCellsPatch[j];
148 
149  // The position could bein any tet of the decomposed cell,
150  // so arbitrarily choose the first face of the cell as the
151  // tetFace and the first point on the face after the base
152  // point as the tetPt. The tracking will pick the cell
153  // consistent with the motion in the first tracking step.
154  const label tetFacei = this->owner().mesh().cells()[celli][0];
155  const label tetPtI = 1;
156 
157 // const point& pos = this->owner().mesh().C()[celli];
158 
159  const scalar offset =
160  max
161  (
163  deltaFilmPatch_[primaryPatchi][j]
164  );
165  const point pos = Cf[j] - 1.1*offset*Sf[j]/magSf[j];
166 
167  // Create a new parcel
168  parcelType* pPtr =
169  new parcelType
170  (
171  this->owner().pMesh(),
172  pos,
173  celli,
174  tetFacei,
175  tetPtI
176  );
177 
178  // Check/set new parcel thermo properties
179  td.cloud().setParcelThermoProperties(*pPtr, 0.0);
180 
181  setParcelProperties(*pPtr, j);
182 
183  if (pPtr->nParticle() > 0.001)
184  {
185  // Check new parcel properties
186  // td.cloud().checkParcelProperties(*pPtr, 0.0, true);
187  td.cloud().checkParcelProperties(*pPtr, 0.0, false);
188 
189  // Add the new parcel to the cloud
190  td.cloud().addParticle(pPtr);
191 
193  }
194  else
195  {
196  // TODO: cache mass and re-distribute?
197  delete pPtr;
198  }
199  }
200  }
201  }
202 }
203 
204 
205 template<class CloudType>
207 (
208  const label filmPatchi,
209  const label primaryPatchi,
211 )
212 {
213  massParcelPatch_ = filmModel.cloudMassTrans().boundaryField()[filmPatchi];
214  filmModel.toPrimary(filmPatchi, massParcelPatch_);
215 
217  filmModel.cloudDiameterTrans().boundaryField()[filmPatchi];
218  filmModel.toPrimary(filmPatchi, diameterParcelPatch_, maxEqOp<scalar>());
219 
220  UFilmPatch_ = filmModel.Us().boundaryField()[filmPatchi];
221  filmModel.toPrimary(filmPatchi, UFilmPatch_);
222 
223  rhoFilmPatch_ = filmModel.rho().boundaryField()[filmPatchi];
224  filmModel.toPrimary(filmPatchi, rhoFilmPatch_);
225 
226  deltaFilmPatch_[primaryPatchi] =
227  filmModel.delta().boundaryField()[filmPatchi];
228  filmModel.toPrimary(filmPatchi, deltaFilmPatch_[primaryPatchi]);
229 }
230 
231 
232 template<class CloudType>
234 (
235  parcelType& p,
236  const label filmFacei
237 ) const
238 {
239  // Set parcel properties
240  scalar vol = mathematical::pi/6.0*pow3(diameterParcelPatch_[filmFacei]);
241  p.d() = diameterParcelPatch_[filmFacei];
242  p.U() = UFilmPatch_[filmFacei];
243  p.rho() = rhoFilmPatch_[filmFacei];
244 
245  p.nParticle() = massParcelPatch_[filmFacei]/p.rho()/vol;
246 
247  if (ejectedParcelType_ >= 0)
248  {
249  p.typeId() = ejectedParcelType_;
250  }
251 }
252 
253 
254 template<class CloudType>
256 {
257  label nTrans0 =
258  this->template getModelProperty<label>("nParcelsTransferred");
259 
260  label nInject0 =
261  this->template getModelProperty<label>("nParcelsInjected");
262 
263  label nTransTotal =
265 
266  label nInjectTotal =
268 
269  os << " Parcels absorbed into film = " << nTransTotal << nl
270  << " New film detached parcels = " << nInjectTotal << endl;
271 
272  if (this->writeTime())
273  {
274  this->setModelProperty("nParcelsTransferred", nTransTotal);
275  this->setModelProperty("nParcelsInjected", nInjectTotal);
277  nParcelsInjected_ = 0;
278  }
279 }
280 
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 #include "SurfaceFilmModelNew.C"
285 
286 // ************************************************************************* //
Collection of constants.
virtual void info(Ostream &os)
Write surface film info to stream.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const surfaceVectorField & Sf() const
Return cell face area vectors.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
const Switch & active() const
Return the active flag.
Definition: regionModelI.H:43
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
label nParcelsInjected_
Number of parcels injected from the film model.
Base class for cloud sub-models.
const volVectorField & C() const
Return cell centres as volVectorField.
const cellList & cells() const
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmModel &filmModel)
Cache the film fields in preparation for injection.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
virtual bool writeTime() const
Flag to indicate when to write a property.
dimensionedScalar pos(const dimensionedScalar &ds)
scalarListList deltaFilmPatch_
Film height of all film patches / patch face.
dynamicFvMesh & mesh
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
void inject(TrackData &td)
Inject parcels into the cloud.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A class for handling words, derived from string.
Definition: word.H:59
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Foam::polyBoundaryMesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
List< vector > UFilmPatch_
Film velocity / patch face.
const dimensionedVector & g
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:172
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
Templated wall surface film model class.
virtual ~SurfaceFilmModel()
Destructor.
dimensionedScalar pow3(const dimensionedScalar &ds)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
CloudType::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 volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
scalarList diameterParcelPatch_
Parcel diameter / patch face.
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
label nParcelsTransferred_
Number of parcels transferred to the film model.
const CloudType & owner() const
Return const access to the owner cloud.
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 dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
virtual bool active() const
Return the model &#39;active&#39; status - default active = true.
Definition: subModelBase.C:154
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
label ejectedParcelType_
Ejected parcel type label - id assigned to identify parcel for.
const dimensionedVector & g_
Gravitational acceleration constant.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
scalarList massParcelPatch_
Parcel mass / patch face.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].