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-2024 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 "mappedFvPatchBaseBase.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  film_.alpha.name(),
79  film_.thermo.he().name(),
80  film_.U.name()
81  };
82 }
83 
84 
86 {
87  if (ejection_.valid() && correctEjection_)
88  {
89  ejection_->correct();
90 
91  // Do not correct ejection rate until the cloud has evolved
92  // to include the last set of ejected parcels
93  correctEjection_ = false;
94  }
95 }
96 
97 
98 template<class Type>
100 inline Foam::fv::filmCloudTransfer::CloudToFilmTransferRate
101 (
102  const Field<Type>& prop,
103  const dimensionSet& dimProp
104 ) const
105 {
107  (
109  (
110  "Su",
111  mesh(),
113  )
114  );
115 
116  if (cloudFieldsTransferred_)
117  {
118  const fvMesh& cloudMesh = film_.surfacePatchMap().nbrMesh();
119  const label cloudPatchi = film_.surfacePatchMap().nbrFvPatch().index();
120 
121  UIndirectList<Type>(tSu.ref(), film_.surfacePatch().faceCells()) =
122  film_.surfacePatchMap().fromNeighbour
123  (
124  prop/cloudMesh.boundary()[cloudPatchi].magSf()
125  );
126 
127  tSu.ref().primitiveFieldRef() /= film_.VbyA;
128  tSu.ref().primitiveFieldRef() /= mesh().time().deltaTValue();
129  }
130 
131  return tSu;
132 }
133 
134 
136 (
137  const volScalarField& rho,
138  const volScalarField& alpha,
139  fvMatrix<scalar>& eqn
140 ) const
141 {
142  if (debug)
143  {
144  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
145  }
146 
147  if (&alpha == &film_.alpha && &eqn.psi() == &film_.alpha)
148  {
149  eqn += CloudToFilmTransferRate<scalar>(massFromCloud_, dimMass);
150 
151  if (ejection_.valid())
152  {
153  eqn -= fvm::Sp(ejection_->rate()*rho(), eqn.psi());
154  }
155  }
156  else
157  {
159  << "Support for field " << alpha.name() << " is not implemented"
160  << exit(FatalError);
161  }
162 }
163 
164 
166 (
167  const volScalarField& alpha,
168  const volScalarField& rho,
169  const volScalarField& he,
170  fvMatrix<scalar>& eqn
171 ) const
172 {
173  if (debug)
174  {
175  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
176  }
177 
178  if (&he == &film_.thermo.he() && &eqn.psi() == &film_.thermo.he())
179  {
180  eqn += CloudToFilmTransferRate<scalar>(energyFromCloud_, dimEnergy);
181 
182  if (ejection_.valid())
183  {
184  eqn -= fvm::Sp(alpha()*rho()*ejection_->rate(), eqn.psi());
185  }
186  }
187  else
188  {
190  << "Support for field " << he.name() << " is not implemented"
191  << exit(FatalError);
192  }
193 }
194 
195 
197 (
198  const volScalarField& alpha,
199  const volScalarField& rho,
200  const volVectorField& U,
201  fvMatrix<vector>& eqn
202 ) const
203 {
204  if (debug)
205  {
206  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
207  }
208 
209  if (&U == &film_.U && &U == &film_.U)
210  {
211  eqn += CloudToFilmTransferRate<vector>(momentumFromCloud_, dimMomentum);
212 
213  if (ejection_.valid())
214  {
215  eqn -= fvm::Sp(alpha()*rho()*ejection_->rate(), eqn.psi());
216  }
217  }
218  else
219  {
221  << "Support for field " << U.name() << " is not implemented"
222  << exit(FatalError);
223  }
224 }
225 
226 
228 {
229  const fvMesh& cloudMesh = film_.surfacePatchMap().nbrMesh();
230  const label cloudPatchi = film_.surfacePatchMap().nbrFvPatch().index();
231  const label nCloudPatchFaces = cloudMesh.boundary()[cloudPatchi].size();
232 
233  if (massFromCloud_.size() != nCloudPatchFaces)
234  {
235  massFromCloud_.setSize(nCloudPatchFaces);
236  momentumFromCloud_.setSize(nCloudPatchFaces);
237  energyFromCloud_.setSize(nCloudPatchFaces);
238  }
239 
240  massFromCloud_ = 0;
241  momentumFromCloud_ = Zero;
242  energyFromCloud_ = 0;
243 
244  cloudFieldsTransferred_ = true;
245 
246  // Enable ejection correction on next call to correct()
247  correctEjection_ = true;
248 }
249 
250 
252 (
253  const label facei,
254  const scalar mass,
255  const vector& momentum,
256  const scalar energy
257 )
258 {
259  massFromCloud_[facei] += mass;
260  momentumFromCloud_[facei] += momentum;
261  energyFromCloud_[facei] += energy;
262 }
263 
264 
265 template<class Type>
267 inline Foam::fv::filmCloudTransfer::filmToCloudTransfer
268 (
269  const VolInternalField<Type>& prop
270 ) const
271 {
272  return film_.surfacePatchMap().toNeighbour
273  (
274  Field<Type>(prop, film_.surfacePatch().faceCells())
275  );
276 }
277 
278 
280 {
281  return ejection_.valid();
282 }
283 
284 
287 {
288  return filmToCloudTransfer<scalar>
289  (
290  (
291  mesh().V()
292  *mesh().time().deltaTValue()
293  *film_.alpha()*film_.rho()*ejection_->rate()
294  )()
295  );
296 }
297 
298 
301 {
302  return filmToCloudTransfer<scalar>(ejection_->diameter());
303 }
304 
305 
308 {
309  return filmToCloudTransfer<scalar>(film_.delta);
310 }
311 
312 
315 {
316  return filmToCloudTransfer<vector>(film_.U);
317 }
318 
319 
322 {
323  return filmToCloudTransfer<scalar>(film_.rho);
324 }
325 
326 
329 {
330  return filmToCloudTransfer<scalar>(film_.thermo.T());
331 }
332 
333 
336 {
337  return filmToCloudTransfer<scalar>(film_.thermo.Cp());
338 }
339 
340 
342 {
343  // Set the cloud field state to false, will be updated by the cloud tracking
344  // If the film is evaluated before the cloud it would be better
345  // if the cloud fields were mapped
346  cloudFieldsTransferred_ = false;
347 
348  if (ejection_.valid())
349  {
350  ejection_->topoChange(map);
351  }
352 }
353 
354 
356 {
357  // Set the cloud field state to false, will be updated by the cloud tracking
358  // If the film is evaluated before the cloud it would be better
359  // if the cloud fields were mapped
360  cloudFieldsTransferred_ = false;
361 
362  if (ejection_.valid())
363  {
364  ejection_->mapMesh(map);
365  }
366 }
367 
368 
370 {
371  // Set the cloud field state to false, will be updated by the cloud tracking
372  // If the film is evaluated before the cloud it would be better
373  // if the cloud fields were mapped
374  cloudFieldsTransferred_ = false;
375 
376  if (ejection_.valid())
377  {
378  ejection_->distribute(map);
379  }
380 }
381 
382 
384 {
385  if (ejection_.valid())
386  {
387  ejection_->movePoints();
388  }
389 
390  return true;
391 }
392 
393 
394 // ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
Pre-declare SubField and related Field type.
Definition: Field.H:83
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:162
Dimension set for the base types.
Definition: dimensionSet.H:125
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:99
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
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.
virtual void addSup(const volScalarField &rho, const volScalarField &alpha, fvMatrix< scalar > &eqn) const
Add source to phase continuity equation.
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 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.
void parcelFromCloud(const label facei, const scalar mass, const vector &momentum, const scalar energy)
Transfer parcel properties from cloud to the film.
tmp< Field< scalar > > ejectedMassToCloud() const
Transfer the ejected mass to the cloud.
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:56
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:334
Calculate the matrix for implicit and explicit sources.
U
Definition: pEqn.H:72
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
const dimensionSet dimMomentum
const dimensionSet dimTime
const dimensionSet dimVolume
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
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
thermo he()
labelList fv(nPoints)
dictionary dict