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-2015 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 :
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 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  (
162  diameterParcelPatch_[j],
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 
192  nParcelsInjected_++;
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 
216  diameterParcelPatch_ =
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 =
264  nTrans0 + returnReduce(nParcelsTransferred_, sumOp<label>());
265 
266  label nInjectTotal =
267  nInject0 + returnReduce(nParcelsInjected_, sumOp<label>());
268 
269  os << " Surface film:" << nl
270  << " - parcels absorbed = " << nTransTotal << nl
271  << " - parcels ejected = " << nInjectTotal << endl;
272 
273  if (this->outputTime())
274  {
275  this->setModelProperty("nParcelsTransferred", nTransTotal);
276  this->setModelProperty("nParcelsInjected", nInjectTotal);
277  nParcelsTransferred_ = 0;
278  nParcelsInjected_ = 0;
279  }
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 
285 #include "SurfaceFilmModelNew.C"
286 
287 // ************************************************************************* //
Foam::regionModels::regionModel::active
const Switch & active() const
Return the active flag.
Definition: regionModelI.H:43
mathematicalConstants.H
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::SurfaceFilmModel::g_
const dimensionedVector & g_
Gravitational acceleration constant.
Definition: SurfaceFilmModel.H:74
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::SurfaceFilmModel::deltaFilmPatch_
scalarListList deltaFilmPatch_
Film height of all film patches / patch face.
Definition: SurfaceFilmModel.H:97
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::SurfaceFilmModel::massParcelPatch_
scalarList massParcelPatch_
Parcel mass / patch face.
Definition: SurfaceFilmModel.H:85
Foam::maxEqOp
Definition: ops.H:77
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::CloudSubModelBase
Base class for cloud sub-models.
Definition: CloudSubModelBase.H:49
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::rho
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
Foam::SurfaceFilmModel::nParcelsTransferred_
label nParcelsTransferred_
Number of parcels transferred to the film model.
Definition: SurfaceFilmModel.H:103
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::delta
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Foam::SurfaceFilmModel::inject
void inject(TrackData &td)
Inject parcels into the cloud.
Definition: SurfaceFilmModel.C:104
Foam::SurfaceFilmModel::~SurfaceFilmModel
virtual ~SurfaceFilmModel()
Destructor.
Definition: SurfaceFilmModel.C:96
Foam::SurfaceFilmModel::nParcelsInjected_
label nParcelsInjected_
Number of parcels injected from the film model.
Definition: SurfaceFilmModel.H:106
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::Us
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
Foam::SurfaceFilmModel< CloudType >::parcelType
CloudType::parcelType parcelType
Convenience typedef to the cloud's parcel type.
Definition: SurfaceFilmModel.H:71
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::SurfaceFilmModel< CloudType >
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
Foam::nl
static const char nl
Definition: Ostream.H:260
surfaceFilmModel.H
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::cloudMassTrans
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::SurfaceFilmModel::diameterParcelPatch_
scalarList diameterParcelPatch_
Parcel diameter / patch face.
Definition: SurfaceFilmModel.H:88
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::cloudDiameterTrans
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::regionModels::regionModel::primaryPatchIDs
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:172
Foam::SurfaceFilmModel::rhoFilmPatch_
scalarList rhoFilmPatch_
Film density / patch face.
Definition: SurfaceFilmModel.H:94
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::sumOp
Definition: ops.H:162
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::SurfaceFilmModel::cacheFilmFields
virtual void cacheFilmFields(const label filmPatchI, const label primaryPatchI, const regionModels::surfaceFilmModels::surfaceFilmModel &filmModel)
Cache the film fields in preparation for injection.
Definition: SurfaceFilmModel.C:207
SurfaceFilmModelNew.C
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
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:161
Foam::SurfaceFilmModel::info
virtual void info(Ostream &os)
Write surface film info to stream.
Definition: SurfaceFilmModel.C:255
Foam::SurfaceFilmModel::ejectedParcelType_
label ejectedParcelType_
Ejected parcel type label - id assigned to identify parcel for.
Definition: SurfaceFilmModel.H:79
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
SurfaceFilmModel.H
Foam::SurfaceFilmModel::UFilmPatch_
List< vector > UFilmPatch_
Film velocity / patch face.
Definition: SurfaceFilmModel.H:91
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::SurfaceFilmModel::setParcelProperties
virtual void setParcelProperties(parcelType &p, const label filmFaceI) const
Set the individual parcel properties.
Definition: SurfaceFilmModel.C:234
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::SurfaceFilmModel::SurfaceFilmModel
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
Definition: SurfaceFilmModel.C:35
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190