isothermalFilm.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 "isothermalFilm.H"
27 #include "filmWallPolyPatch.H"
28 #include "filmSurfacePolyPatch.H"
29 #include "mappedFvPatchBaseBase.H"
32 #include "constantSurfaceTension.H"
33 #include "fvcVolumeIntegrate.H"
34 #include "fvcDdt.H"
35 #include "fvcDiv.H"
36 #include "fvcFlux.H"
37 #include "fvcSurfaceIntegrate.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace solvers
45 {
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 bool Foam::solvers::isothermalFilm::initFilmMesh()
55 {
56  // Search for film wall patches
57 
58  label nWallFaces = 0;
59  DynamicList<label> wallPatchIndices_;
60 
61  const polyBoundaryMesh& bm = mesh.boundaryMesh();
62 
63  forAll(bm, patchi)
64  {
65  const polyPatch& p = bm[patchi];
66 
67  if (isA<filmWallPolyPatch>(p))
68  {
69  wallPatchIndices_.append(patchi);
70  nWallFaces += p.faceCells().size();
71  }
72  }
73 
74  if (nWallFaces != mesh.nCells())
75  {
77  << "The number of film wall faces in the mesh "
78  << nWallFaces
79  << " is not equal to the number of cells "
80  << mesh.nCells()
81  << exit(FatalError);
82  }
83 
84  if (returnReduce(nWallFaces, sumOp<label>()) == 0)
85  {
87  << "There are no filmWall faces in the mesh"
88  << exit(FatalError);
89  }
90 
91  wallPatchIDs.transfer(wallPatchIndices_);
92 
93 
94  // Search for film surface patch
95  surfacePatchID = -1;
96 
97  forAll(bm, patchi)
98  {
99  const polyPatch& p = bm[patchi];
100 
101  if (isA<filmSurfacePolyPatch>(p))
102  {
103  if (surfacePatchID == -1)
104  {
106  }
107  else
108  {
110  << "More than one filmSurface patch defined: "
111  << surfacePatchID << " and " << patchi
112  << exit(FatalError);
113  }
114  }
115  }
116 
117  if (surfacePatchID == -1)
118  {
119  Info<< "The filmSurface patch is not defined"
120  << endl;
121  }
122 
123 
124  // Calculate film specific mesh geometry from the film wall patches
125 
126  forAll(wallPatchIDs, i)
127  {
128  const label patchi = wallPatchIDs[i];
129  const polyPatch& wallp = bm[patchi];
130  const labelList& fCells = wallp.faceCells();
131 
132  UIndirectList<vector>(nHat_, fCells) = wallp.faceNormals();
133  UIndirectList<scalar>(magSf_, fCells) = wallp.magFaceAreas();
134  }
135 
137 
140 
141  return true;
142 }
143 
144 
145 Foam::wordList Foam::solvers::isothermalFilm::alphaTypes() const
146 {
147  wordList alphaTypes(delta_.boundaryField().types());
148 
149  forAll(delta_.boundaryField(), patchi)
150  {
151  if (!delta_.boundaryField()[patchi].assignable())
152  {
153  alphaTypes[patchi] = fixedValueFvPatchScalarField::typeName;
154  }
155  }
156 
157  forAll(wallPatchIDs, i)
158  {
159  alphaTypes[wallPatchIDs[i]] = alphaOneFvPatchScalarField::typeName;
160  }
161 
162  if (surfacePatchID != -1)
163  {
164  alphaTypes[surfacePatchID] = alphaOneFvPatchScalarField::typeName;
165  }
166 
167  return alphaTypes;
168 }
169 
170 
171 void Foam::solvers::isothermalFilm::correctCoNum()
172 {
173  const scalarField sumPhi(fvc::surfaceSum(mag(phi))().primitiveField());
174 
175  CoNum = 0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue();
176 
177  const scalar meanCoNum =
178  0.5
179  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
180  *runTime.deltaTValue();
181 
182  Info<< "Courant Number mean: " << meanCoNum
183  << " max: " << CoNum << endl;
184 }
185 
186 
187 void Foam::solvers::isothermalFilm::continuityErrors()
188 {
189  const dimensionedScalar mass = fvc::domainIntegrate(rho()*delta()*magSf);
190 
191  correctContinuityError();
192 
193  if (mass.value() > small)
194  {
195  const volScalarField::Internal massContErr
196  (
197  runTime.deltaT()*magSf*contErr()
198  );
199 
200  const scalar sumLocalContErr =
201  (fvc::domainIntegrate(mag(massContErr))/mass).value();
202 
203  const scalar globalContErr =
204  (fvc::domainIntegrate(massContErr)/mass).value();
205 
206  Info<< "time step continuity errors : sum local = " << sumLocalContErr
207  << ", global = " << globalContErr;
208 
209  if (pimple.finalPisoIter() && pimple.finalIter())
210  {
212 
213  Info<< ", cumulative = " << cumulativeContErr;
214  }
215 
216  Info<< endl;
217  }
218 }
219 
220 
221 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
222 
224 {
225  return runTime.controlDict().modified();
226 }
227 
228 
230 {
231  solver::read();
232 
233  maxCo =
234  runTime.controlDict().lookupOrDefault<scalar>("maxCo", vGreat);
235 
236  maxDeltaT_ =
237  runTime.controlDict().found("maxDeltaT")
238  ? runTime.controlDict().lookup<scalar>("maxDeltaT", runTime.userUnits())
239  : vGreat;
240 
241  return true;
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
248 (
249  fvMesh& mesh,
250  autoPtr<rhoFluidThermo> thermoPtr
251 )
252 :
253  solver(mesh),
254 
255  CoNum(0),
257 
258  thermoPtr_(thermoPtr),
259  thermo_(thermoPtr_()),
260 
261  p(thermo_.p()),
262 
263  nHat_
264  (
265  IOobject
266  (
267  "nHat",
268  runTime.name(),
269  mesh
270  ),
271  mesh,
274  ),
275 
276  magSf_
277  (
278  IOobject
279  (
280  "magSf",
281  runTime.name(),
282  mesh
283  ),
284  mesh,
286  ),
287 
288  VbyA_
289  (
290  IOobject
291  (
292  "VbyA",
293  runTime.name(),
294  mesh
295  ),
296  mesh,
299  ),
300 
301  initialised_(initFilmMesh()),
302 
303  delta_
304  (
305  IOobject
306  (
307  "delta",
308  runTime.name(),
309  mesh,
310  IOobject::MUST_READ,
311  IOobject::AUTO_WRITE
312  ),
313  mesh
314  ),
315 
316  alpha_
317  (
318  IOobject
319  (
320  "alpha",
321  runTime.name(),
322  mesh,
323  IOobject::READ_IF_PRESENT,
324  IOobject::AUTO_WRITE
325  ),
326  delta_/VbyA_,
327  alphaTypes()
328  ),
329 
330  deltaWet("deltaWet", dimLength, thermo_.properties()),
331 
332  U_
333  (
334  IOobject
335  (
336  "U",
337  runTime.name(),
338  mesh,
339  IOobject::MUST_READ,
340  IOobject::AUTO_WRITE
341  ),
342  mesh
343  ),
344 
345  alphaRhoPhi_
346  (
347  IOobject
348  (
349  "alphaRhoPhi",
350  runTime.name(),
351  mesh,
352  IOobject::READ_IF_PRESENT,
353  IOobject::AUTO_WRITE
354  ),
355  fvc::flux(alpha_*thermo_.rho()*U_)
356  ),
357 
358  phi_
359  (
360  IOobject
361  (
362  "phi",
363  runTime.name(),
364  mesh,
365  IOobject::READ_IF_PRESENT,
366  IOobject::AUTO_WRITE
367  ),
368  fvc::flux(U_)
369  ),
370 
371  surfaceTension(surfaceTensionModel::New(thermo_.properties(), mesh)),
372  thermocapillary(!isType<surfaceTensionModels::constant>(surfaceTension())),
373 
374  g
375  (
376  IOobject
377  (
378  "g",
379  runTime.constant(),
380  mesh,
381  IOobject::MUST_READ,
382  IOobject::NO_WRITE
383  )
384  ),
385 
386  nHat(nHat_),
387  magSf(magSf_),
388  VbyA(VbyA_),
389  delta(delta_),
390  alpha(alpha_),
391  thermo(thermo_),
392  rho(thermo_.rho()),
393  U(U_),
394  alphaRhoPhi(alphaRhoPhi_),
395  phi(phi_),
396 
397  momentumTransport
398  (
399  filmCompressible::momentumTransportModel::New
400  (
401  alpha,
402  thermo.rho(),
403  U,
404  alphaRhoPhi,
405  phi,
406  thermo
407  )
408  )
409 {
410  // Read the controls
411  read();
412 
414  momentumTransport->validate();
415 }
416 
417 
419 :
420  isothermalFilm(mesh, rhoFluidThermo::New(mesh))
421 {}
422 
423 
424 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
425 
427 {}
428 
429 
430 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
431 
433 {
434  return mesh.boundary()[surfacePatchID];
435 }
436 
437 
440 {
441  return refCast<const mappedFvPatchBaseBase>(surfacePatch());
442 }
443 
444 
446 {
447  scalar deltaT = min(fvModels().maxDeltaT(), maxDeltaT_);
448 
449  if (CoNum > small)
450  {
451  deltaT = min(deltaT, maxCo/CoNum*runTime.deltaTValue());
452  }
453 
454  return deltaT;
455 }
456 
457 
459 {
460  correctCoNum();
461 }
462 
463 
465 {}
466 
467 
469 {}
470 
471 
473 {
474  thermo_.correct();
475 }
476 
477 
479 {
480  correctAlpha();
481 }
482 
483 
485 {
486  if (pimple.correctTransport())
487  {
488  momentumTransport->correct();
489  }
490 }
491 
492 
494 {}
495 
496 
497 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static const char *const typeName
Definition: Field.H:106
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
void correctBoundaryConditions()
Correct boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
Base class for fv patches that provide mapping between two fv patches.
Abstract base class for turbulence models (RAS, LES and laminar).
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
label nCells() const
Base-class for fluid thermodynamic properties based on density.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
virtual bool read()
Read controls.
Definition: solver.C:52
Solver module for flow of compressible isothermal liquid films.
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
volScalarField & p
The thermodynamic pressure field.
autoPtr< filmCompressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
virtual bool dependenciesModified() const
Return true if the solver's dependencies have been modified.
const volScalarField & alpha
Film volume fraction in the cell layer.
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
virtual void motionCorrector()
Corrections that follow mesh motion.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
isothermalFilm(fvMesh &mesh, autoPtr< rhoFluidThermo >)
Construct from region mesh and thermophysical properties.
volScalarField VbyA_
Film cell volume/wall face area.
const mappedFvPatchBaseBase & surfacePatchMap() const
Return the film surface patch region-region map.
const fvPatch & surfacePatch() const
Return the film surface patch.
volScalarField::Internal magSf_
Film cell cross-sectional area magnitude.
volVectorField nHat_
Film wall normal.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
label surfacePatchID
Film surface patch ID.
virtual ~isothermalFilm()
Destructor.
virtual bool read()
Read controls.
labelList wallPatchIDs
List of film wall patch IDs.
'Patch' on surface as subset of triSurface.
Definition: surfacePatch.H:61
Abstract base-class for surface tension models which return the surface tension coefficient field.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
scalar CoNum
scalar meanCoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
scalar maxDeltaT
scalar maxCo
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Volume integrate volField creating a volField.
label patchi
scalar cumulativeContErr
label nWallFaces(0)
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
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
List< word > wordList
A List of words.
Definition: fileName.H:54
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:158
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
const dimensionSet dimArea
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
scalar sumLocalContErr
Definition: volContinuity.H:9
scalar globalContErr
Definition: volContinuity.H:12