filmCloudTransfer.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) 2023 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 "filmCloudTransfer.H"
27 #include "mappedPatchBase.H"
28 #include "fvmSup.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
38 
40  (
41  fvModel,
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
59  fvModel(sourceName, modelType, mesh, dict),
60  film_(mesh.lookupObject<solvers::isothermalFilm>(solver::typeName)),
61  cloudFieldsTransferred_(false),
62  correctEjection_(false),
63  ejection_
64  (
65  dict.found("ejection")
66  ? ejectionModel::New(dict.subDict("ejection"), film_)
67  : autoPtr<ejectionModel>(nullptr)
68  )
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
76  return wordList
77  {
78  "pi",
79  film_.alpha.name(),
80  film_.thermo.he().name(),
81  film_.U.name()
82  };
83 }
84 
85 
87 {
88  if (ejection_.valid() && correctEjection_)
89  {
90  ejection_->correct();
91 
92  // Do not correct ejection rate until the cloud has evolved
93  // to include the last set of ejected parcels
94  correctEjection_ = false;
95  }
96 }
97 
98 
99 template<class Type>
101 inline Foam::fv::filmCloudTransfer::CloudToFilmTransferRate
102 (
103  const Field<Type>& prop,
104  const dimensionSet& dimProp
105 ) const
106 {
108  (
110  (
111  "Su",
112  mesh(),
114  )
115  );
116 
117  if (cloudFieldsTransferred_)
118  {
119  const fvMesh& cloudMesh =
120  refCast<const fvMesh>(film_.surfacePatchMap().nbrMesh());
121 
122  const label cloudPatchi =
123  film_.surfacePatchMap().nbrPolyPatch().index();
124 
125  UIndirectList<Type>(tSu.ref(), film_.surfacePatch().faceCells()) =
126  film_.surfacePatchMap().fromNeighbour
127  (
128  prop/cloudMesh.boundary()[cloudPatchi].magSf()
129  );
130 
131  tSu.ref().field() /= film_.VbyA;
132  tSu.ref().field() /= mesh().time().deltaTValue();
133  }
134 
135  return tSu;
136 }
137 
138 
140 (
141  fvMatrix<scalar>& eqn,
142  const word& fieldName
143 ) const
144 {
145  if (debug)
146  {
147  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
148  }
149 
150  // Droplet impingement pressure
151  if (fieldName == "pi")
152  {
153  eqn +=
154  CloudToFilmTransferRate<scalar>
155  (
156  pressureFromCloud_,
158  );
159  }
160  else
161  {
163  << "Support for field " << fieldName << " is not implemented"
164  << exit(FatalError);
165  }
166 }
167 
168 
170 (
171  const volScalarField& rho,
172  fvMatrix<scalar>& eqn,
173  const word& fieldName
174 ) const
175 {
176  if (debug)
177  {
178  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
179  }
180 
181  if (fieldName == film_.alpha.name())
182  {
183  eqn += CloudToFilmTransferRate<scalar>(massFromCloud_, dimMass);
184 
185  if (ejection_.valid())
186  {
187  eqn -= fvm::Sp(ejection_->rate()*rho(), eqn.psi());
188  }
189  }
190  else
191  {
193  << "Support for field " << fieldName << " is not implemented"
194  << exit(FatalError);
195  }
196 }
197 
198 
200 (
201  const volScalarField& alpha,
202  const volScalarField& rho,
203  fvMatrix<scalar>& eqn,
204  const word& fieldName
205 ) const
206 {
207  if (debug)
208  {
209  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
210  }
211 
212  if (fieldName == film_.thermo.he().name())
213  {
214  eqn += CloudToFilmTransferRate<scalar>(energyFromCloud_, dimEnergy);
215 
216  if (ejection_.valid())
217  {
218  eqn -= fvm::Sp(alpha()*rho()*ejection_->rate(), eqn.psi());
219  }
220  }
221  else
222  {
224  << "Support for field " << fieldName << " is not implemented"
225  << exit(FatalError);
226  }
227 }
228 
229 
231 (
232  const volScalarField& alpha,
233  const volScalarField& rho,
234  fvMatrix<vector>& eqn,
235  const word& fieldName
236 ) const
237 {
238  if (debug)
239  {
240  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
241  }
242 
243  eqn += CloudToFilmTransferRate<vector>(momentumFromCloud_, dimMomentum);
244 
245  if (ejection_.valid())
246  {
247  eqn -= fvm::Sp(alpha()*rho()*ejection_->rate(), eqn.psi());
248  }
249 }
250 
251 
253 {
254  const fvMesh& cloudMesh =
255  refCast<const fvMesh>(film_.surfacePatchMap().nbrMesh());
256  const label cloudPatchi = film_.surfacePatchMap().nbrPolyPatch().index();
257  const label nCloudPatchFaces = cloudMesh.boundary()[cloudPatchi].size();
258 
259  if (massFromCloud_.size() != nCloudPatchFaces)
260  {
261  massFromCloud_.setSize(nCloudPatchFaces);
262  momentumFromCloud_.setSize(nCloudPatchFaces);
263  pressureFromCloud_.setSize(nCloudPatchFaces);
264  energyFromCloud_.setSize(nCloudPatchFaces);
265  }
266 
267  massFromCloud_ = 0;
268  momentumFromCloud_ = Zero;
269  pressureFromCloud_ = 0;
270  energyFromCloud_ = 0;
271 
272  cloudFieldsTransferred_ = true;
273 
274  // Enable ejection correction on next call to correct()
275  correctEjection_ = true;
276 }
277 
278 
280 (
281  const label facei,
282  const scalar mass,
283  const vector& momentum,
284  const scalar pressure,
285  const scalar energy
286 )
287 {
288  massFromCloud_[facei] += mass;
289  momentumFromCloud_[facei] += momentum;
290  pressureFromCloud_[facei] += pressure;
291  energyFromCloud_[facei] += energy;
292 }
293 
294 
295 template<class Type>
297 inline Foam::fv::filmCloudTransfer::filmToCloudTransfer
298 (
299  const VolInternalField<Type>& prop
300 ) const
301 {
302  return film_.surfacePatchMap().toNeighbour
303  (
304  Field<Type>(prop, film_.surfacePatch().faceCells())
305  );
306 }
307 
308 
310 {
311  return ejection_.valid();
312 }
313 
314 
317 {
318  return filmToCloudTransfer<scalar>
319  (
320  (
321  mesh().V()
322  *mesh().time().deltaTValue()
323  *film_.alpha()*film_.rho()*ejection_->rate()
324  )()
325  );
326 }
327 
328 
331 {
332  return filmToCloudTransfer<scalar>(ejection_->diameter());
333 }
334 
335 
338 {
339  return filmToCloudTransfer<scalar>(film_.delta);
340 }
341 
342 
345 {
346  return filmToCloudTransfer<vector>(film_.U);
347 }
348 
349 
352 {
353  return filmToCloudTransfer<scalar>(film_.rho);
354 }
355 
356 
359 {
360  return filmToCloudTransfer<scalar>(film_.thermo.T());
361 }
362 
363 
366 {
367  return filmToCloudTransfer<scalar>(film_.thermo.Cp());
368 }
369 
370 
372 {
373  // Set the cloud field state to false, will be updated by the cloud tracking
374  // If the film is evaluated before the cloud it would be better
375  // if the cloud fields were mapped
376  cloudFieldsTransferred_ = false;
377 
378  if (ejection_.valid())
379  {
380  ejection_->topoChange(map);
381  }
382 }
383 
384 
386 {
387  // Set the cloud field state to false, will be updated by the cloud tracking
388  // If the film is evaluated before the cloud it would be better
389  // if the cloud fields were mapped
390  cloudFieldsTransferred_ = false;
391 
392  if (ejection_.valid())
393  {
394  ejection_->mapMesh(map);
395  }
396 }
397 
398 
400 {
401  // Set the cloud field state to false, will be updated by the cloud tracking
402  // If the film is evaluated before the cloud it would be better
403  // if the cloud fields were mapped
404  cloudFieldsTransferred_ = false;
405 
406  if (ejection_.valid())
407  {
408  ejection_->distribute(map);
409  }
410 }
411 
412 
414 {
415  if (ejection_.valid())
416  {
417  ejection_->movePoints();
418  }
419 
420  return true;
421 }
422 
423 
424 // ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricField class.
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
Finite volume model abstract base class.
Definition: fvModel.H:59
Film<->cloud transfer model.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
tmp< Field< scalar > > deltaToCloud() const
Transfer the film delta field to the cloud.
void resetFromCloudFields()
Reset the fields accumulated cloud transfer fields.
tmp< Field< scalar > > TToCloud() const
Transfer the film temperature field to the cloud.
virtual void correct()
Solve the film and update the sources.
tmp< Field< scalar > > CpToCloud() const
Transfer the film heat capacity field to the cloud.
tmp< Field< scalar > > ejectedDiameterToCloud() const
Transfer the ejected droplet diameter to the cloud.
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit droplet impingement contribution to pressure field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
tmp< Field< vector > > UToCloud() const
Transfer the film velocity field to the cloud.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
filmCloudTransfer(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
tmp< Field< scalar > > ejectedMassToCloud() const
Transfer the ejected mass to the cloud.
void parcelFromCloud(const label facei, const scalar mass, const vector &momentum, const scalar pressure, const scalar energy)
Transfer parcel properties from cloud to the film.
bool ejecting() const
Return true if the film is ejecting to the cloud.
tmp< Field< scalar > > rhoToCloud() const
Transfer the film density field to the cloud.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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 dimensionSet dimEnergy
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
const dimensionSet dimMomentum
const dimensionSet dimTime
const dimensionSet dimVolume
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
const dimensionSet dimMass
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict