setInterfaceRDeltaT.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 "twoPhaseVoFSolver.H"
27 #include "fvcSmooth.H"
28 #include "fvcSurfaceIntegrate.H"
29 #include "fvcAverage.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
34 (
35  volScalarField& rDeltaT
36 )
37 {
38  const dictionary& pimpleDict = pimple.dict();
39 
40  const scalar maxCo
41  (
42  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.9)
43  );
44 
45  const scalar maxAlphaCo
46  (
47  pimpleDict.lookupOrDefault<scalar>("maxAlphaCo", 0.2)
48  );
49 
50  const scalar rDeltaTSmoothingCoeff
51  (
52  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
53  );
54 
55  const label nAlphaSpreadIter
56  (
57  pimpleDict.lookupOrDefault<label>("nAlphaSpreadIter", 1)
58  );
59 
60  const scalar alphaSpreadDiff
61  (
62  pimpleDict.lookupOrDefault<scalar>("alphaSpreadDiff", 0.2)
63  );
64 
65  const scalar alphaSpreadMax
66  (
67  pimpleDict.lookupOrDefault<scalar>("alphaSpreadMax", 0.99)
68  );
69 
70  const scalar alphaSpreadMin
71  (
72  pimpleDict.lookupOrDefault<scalar>("alphaSpreadMin", 0.01)
73  );
74 
75  const label nAlphaSweepIter
76  (
77  pimpleDict.lookupOrDefault<label>("nAlphaSweepIter", 5)
78  );
79 
80  if (maxAlphaCo < maxCo)
81  {
82  // Further limit the reciprocal time-step
83  // in the vicinity of the interface
84 
86 
87  rDeltaT.ref() = max
88  (
89  rDeltaT(),
90  pos0(alpha1Bar() - alphaSpreadMin)
91  *pos0(alphaSpreadMax - alpha1Bar())
92  *fvc::surfaceSum(mag(phi))()()
93  /((2*maxAlphaCo)*mesh.V())
94  );
95  }
96 
97  // Update the boundary values of the reciprocal time-step
98  rDeltaT.correctBoundaryConditions();
99 
100  Info<< "Flow and interface time scale min/max = "
101  << gMin(1/rDeltaT.primitiveField())
102  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
103 
104  if (rDeltaTSmoothingCoeff < 1.0)
105  {
106  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
107  }
108 
109  if (nAlphaSpreadIter > 0)
110  {
112  (
113  rDeltaT,
114  alpha1,
115  nAlphaSpreadIter,
116  alphaSpreadDiff,
117  alphaSpreadMax,
118  alphaSpreadMin
119  );
120  }
121 
122  if (nAlphaSweepIter > 0)
123  {
124  fvc::sweep(rDeltaT, alpha1, nAlphaSweepIter, alphaSpreadDiff);
125  }
126 
127  Info<< "Smoothed flow time scale min/max = "
128  << gMin(1/rDeltaT.primitiveField())
129  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
130 }
131 
132 
133 // ************************************************************************* //
Generic GeometricField class.
Internal & ref()
Return a reference to the dimensioned internal field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual const dictionary & dict() const
Return the solution dictionary.
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
volScalarField & alpha1
Reference to the phase1-fraction.
virtual void setInterfaceRDeltaT(volScalarField &rDeltaT)
Adjust the rDeltaT in the vicinity of the interface.
scalar maxCo
Area-weighted average a surfaceField creating a volField.
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
void sweep(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2)
Definition: fvcSmooth.C:240
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
void spread(volScalarField &field, const volScalarField &alpha, const label nLayers, const scalar alphaDiff=0.2, const scalar alphaMax=0.99, const scalar alphaMin=0.01)
Definition: fvcSmooth.C:134
dimensionedScalar pos0(const dimensionedScalar &ds)
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:251
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)