Loading...
Searching...
No Matches
potentialFoam.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 potentialFoam
29
30Group
31 grpBasicSolvers
32
33Description
34 Potential flow solver which solves for the velocity potential, to
35 calculate the flux-field, from which the velocity field is obtained by
36 reconstructing the flux.
37
38 \heading Solver details
39 The potential flow solution is typically employed to generate initial fields
40 for full Navier-Stokes codes. The flow is evolved using the equation:
41
42 \f[
43 \laplacian \Phi = \div(\vec{U})
44 \f]
45
46 Where:
47 \vartable
48 \Phi | Velocity potential [m2/s]
49 \vec{U} | Velocity [m/s]
50 \endvartable
51
52 The corresponding pressure field could be calculated from the divergence
53 of the Euler equation:
54
55 \f[
56 \laplacian p + \div(\div(\vec{U}\otimes\vec{U})) = 0
57 \f]
58
59 but this generates excessive pressure variation in regions of large
60 velocity gradient normal to the flow direction. A better option is to
61 calculate the pressure field corresponding to velocity variation along the
62 stream-lines:
63
64 \f[
65 \laplacian p + \div(\vec{F}\cdot\div(\vec{U}\otimes\vec{U})) = 0
66 \f]
67 where the flow direction tensor \f$\vec{F}\f$ is obtained from
68 \f[
69 \vec{F} = \hat{\vec{U}}\otimes\hat{\vec{U}}
70 \f]
71
72 \heading Required fields
73 \plaintable
74 U | Velocity [m/s]
75 \endplaintable
76
77 \heading Optional fields
78 \plaintable
79 p | Kinematic pressure [m2/s2]
80 Phi | Velocity potential [m2/s]
81 | Generated from p (if present) or U if not present
82 \endplaintable
83
84 \heading Options
85 \plaintable
86 -writep | write the Euler pressure
87 -writephi | Write the final volumetric flux
88 -writePhi | Write the final velocity potential
89 -initialiseUBCs | Update the velocity boundaries before solving for Phi
90 \endplaintable
91
92\*---------------------------------------------------------------------------*/
93
94#include "fvCFD.H"
95#include "regionProperties.H"
96#include "pisoControl.H"
97
98// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99
100int main(int argc, char *argv[])
101{
102 argList::addNote
103 (
104 "Potential flow solver which solves for the velocity potential"
105 );
106
107 argList::addOption
108 (
109 "pName",
110 "pName",
111 "Name of the pressure field"
112 );
113
114 argList::addBoolOption
115 (
116 "initialiseUBCs",
117 "Initialise U boundary conditions"
118 );
119
120 argList::addBoolOption
121 (
122 "writephi",
123 "Write the final volumetric flux field"
124 );
125
126 argList::addBoolOption
127 (
128 "writePhi",
129 "Write the final velocity potential field"
130 );
131
132 argList::addBoolOption
133 (
134 "writep",
135 "Calculate and write the Euler pressure field"
136 );
137
138 argList::addBoolOption
139 (
140 "withFunctionObjects",
141 "Execute functionObjects"
142 );
143
144 // Prevent volume BCs from triggering finite-area
145 regionModels::allowFaModels(false);
146
147 #include "addRegionOption.H"
148 #include "addCheckCaseOptions.H"
149 #include "setRootCaseLists.H"
150 #include "createTime.H"
151 #include "createMesh.H"
152
153 pisoControl potentialFlow(mesh, "potentialFlow");
154
155 #include "createFields.H"
156
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159 Info<< nl << "Calculating potential flow" << endl;
160
161 // Since solver contains no time loop it would never execute
162 // function objects so do it ourselves
163 runTime.functionObjects().start();
164
165 MRF.makeRelative(phi);
166 adjustPhi(phi, U, p);
167
168 // Non-orthogonal velocity potential corrector loop
169 while (potentialFlow.correctNonOrthogonal())
170 {
171 fvScalarMatrix PhiEqn
172 (
173 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
174 ==
175 fvc::div(phi)
176 );
177
178 PhiEqn.setReference(PhiRefCell, PhiRefValue);
179 PhiEqn.solve();
180
181 if (potentialFlow.finalNonOrthogonalIter())
182 {
183 phi -= PhiEqn.flux();
184 }
185 }
186
187 MRF.makeAbsolute(phi);
188
189 Info<< "Continuity error = "
190 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
191 << endl;
192
193 U = fvc::reconstruct(phi);
194 U.correctBoundaryConditions();
195
196 Info<< "Interpolated velocity error = "
197 << (sqrt(sum(sqr(fvc::flux(U) - phi)))/sum(mesh.magSf())).value()
198 << endl;
199
200 // Write U
201 U.write();
202
203 // Optionally write the volumetric flux, phi
204 if (args.found("writephi"))
205 {
206 phi.write();
207 }
208
209 // Optionally write velocity potential, Phi
210 if (args.found("writePhi"))
211 {
212 Phi.write();
213 }
214
215 // Calculate the pressure field from the Euler equation
216 if (args.found("writep"))
217 {
218 Info<< nl << "Calculating approximate pressure field" << endl;
219
220 label pRefCell = 0;
221 scalar pRefValue = 0.0;
223 (
224 p,
225 potentialFlow.dict(),
226 pRefCell,
228 );
229
230 // Calculate the flow-direction filter tensor
231 volScalarField magSqrU(magSqr(U));
232 volSymmTensorField F(sqr(U)/(magSqrU + SMALL*average(magSqrU)));
233
234 // Calculate the divergence of the flow-direction filtered div(U*U)
235 // Filtering with the flow-direction generates a more reasonable
236 // pressure distribution in regions of high velocity gradient in the
237 // direction of the flow
238 volScalarField divDivUU
239 (
240 fvc::div
241 (
242 F & fvc::div(phi, U),
243 "div(div(phi,U))"
244 )
245 );
246
247 // Solve a Poisson equation for the approximate pressure
248 while (potentialFlow.correctNonOrthogonal())
249 {
250 fvScalarMatrix pEqn
251 (
252 fvm::laplacian(p) + divDivUU
253 );
254
255 pEqn.setReference(pRefCell, pRefValue);
256 pEqn.solve();
257 }
258
259 p.write();
260 }
261
262 runTime.functionObjects().end();
263
264 runTime.printExecutionTime(Info);
265
266 Info<< "End\n" << endl;
267
268 return 0;
269}
270
271
272// ************************************************************************* //
Required Classes.
Required Classes.
const dictionary & potentialFlow(mesh.solutionDict().subDict("potentialFlow"))
const scalar pRefValue
const label pRefCell
IOMRFZoneList & MRF
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
adjustPhi(phiHbyA, U, p_rgh)
volVectorField F(fluid.F())
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition facAverage.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Type weightedAverage(const UList< scalar > &weights, const UList< Type > &fld)
The local weighted average of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
setRefCell(p, pimple.dict(), pRefCell, pRefValue)