/
Utilities.C
193 lines (166 loc) · 3.95 KB
/
Utilities.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
/*
* Utilities.C
*
* Created on: May 4, 2018
* Author: timo
*/
#include "Utilities.H"
Utilities::Utilities()
{
}
Utilities::Utilities(const Foam::Utilities& aUtilClass)
{
}
Utilities::~Utilities()
{
}
Foam::scalar Utilities::readAdjustTimeStep(Foam::Time& runTime)
{
Foam::scalar adjustTimeStep = runTime.controlDict().lookupOrDefault(
"adjustTimeStep",
false);
return adjustTimeStep;
}
Foam::scalar Utilities::readMaxCo(Foam::Time& runTime)
{
Foam::scalar maxCo = runTime.controlDict().lookupOrDefault("maxCo", 1.0);
return maxCo;
}
Foam::scalar Utilities::readMaxDeltaT(Foam::Time& runTime)
{
scalar maxDeltaT = runTime.controlDict().lookupOrDefault < scalar
> ("maxDeltaT", GREAT);
return maxDeltaT;
}
scalar Utilities::readMaxAcousticCo(Foam::Time& runTime)
{
scalar maxAcousticCo = readScalar(
runTime.controlDict().lookup("maxAcousticCo"));
return maxAcousticCo;
}
scalar Utilities::calculateCoNumber(
Time& runTime,
fvMesh& mesh,
surfaceScalarField& phi)
{
scalar CoNum = 0.0;
if (mesh.nInternalFaces())
{
scalarField sumPhi(fvc::surfaceSum(mag(phi))().primitiveField());
CoNum = 0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue();
}
return CoNum;
}
scalar Utilities::calculateMeanCoNumber(
Time& runTime,
fvMesh& mesh,
surfaceScalarField& phi) {
scalar meanCoNum = 0.0;
if (mesh.nInternalFaces())
{
scalarField sumPhi(fvc::surfaceSum(mag(phi))().primitiveField());
meanCoNum = 0.5
* (gSum(sumPhi) / gSum(mesh.V().field()))
* runTime.deltaTValue();
}
return meanCoNum;
}
scalar Utilities::calculateAcoustincCoNumber(
Time& runTime,
multiphaseCavitationMixture& mixture,
fvMesh& mesh,
surfaceScalarField& phi)
{
scalar acousticCoNum = 0.0;
tmp<volScalarField> psi = mixture.mixturePsi();
if (mesh.nInternalFaces())
{
scalarField sumPhi
(
fvc::surfaceSum(mag(phi))().primitiveField()
);
acousticCoNum = 0.5*gMax
(
fvc::surfaceSum
(
fvc::interpolate(scalar(1)/sqrt(psi.ref()))*mesh.magSf()
)().primitiveField()/mesh.V().field()
)*runTime.deltaTValue();
}
return acousticCoNum;
}
void Utilities::printCoNumbers(
scalar CoNumber,
scalar meanCoNumber,
scalar acousticCoNumber)
{
Info<< "phi Courant Number mean: " << meanCoNumber
<< " max: "
<< CoNumber
<< " acoustic max: "
<< acousticCoNumber
<< endl;
}
void Utilities::setInitialExtendedDeltaT(
bool adjustTimeStep,
scalar CoNum,
scalar maxCo,
scalar maxAcousticCo,
scalar acousticCoNum,
Time& runTime,
scalar maxDeltaT)
{
if (adjustTimeStep)
{
if (CoNum > SMALL)
{
scalar maxDeltaTFact = min(
maxCo / (CoNum + SMALL),
maxAcousticCo / (acousticCoNum + SMALL));
runTime.setDeltaT(
min(maxDeltaTFact * runTime.deltaTValue(), maxDeltaT));
}
}
}
scalar Utilities::calculateAlphaCoNum(
fvMesh& mesh,
Foam::multiphaseCavitationMixture& mixture,
surfaceScalarField& phi,
Time& runTime)
{
scalar alphaCoNum = 0.0;
if (mesh.nInternalFaces())
{
scalarField sumPhi(
mixture.nearInterface()().primitiveField()
* fvc::surfaceSum(mag(phi))().primitiveField());
alphaCoNum = 0.5 * gMax(sumPhi / mesh.V().field())
* runTime.deltaTValue();
}
return alphaCoNum;
}
scalar Utilities::calculateMeanAlphaCoNum(
fvMesh& mesh,
Foam::multiphaseCavitationMixture& mixture,
surfaceScalarField& phi,
Time& runTime)
{
scalar meanAlphaCoNum = 0.0;
if (mesh.nInternalFaces())
{
scalarField sumPhi(
mixture.nearInterface()().primitiveField()
* fvc::surfaceSum(mag(phi))().primitiveField());
meanAlphaCoNum = 0.5
* (gSum(sumPhi) / gSum(mesh.V().field()))
* runTime.deltaTValue();
}
return meanAlphaCoNum;
}
void Utilities::printAlphaCoNumbers(scalar meanAlphaCoNum, scalar alphaCoNum)
{
Info<< "Interface Courant Number mean: " << meanAlphaCoNum
<< " max: "
<< alphaCoNum
<< endl;
}