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 | Copyright (C) 2011-2015 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 Application
25  potentialFoam
26 
27 Group
28  grpBasicSolvers
29 
30 Description
31  Potential flow solver.
32 
33  \heading Solver details
34  The potential flow solution is typically employed to generate initial fields
35  for full Navier-Stokes codes. The flow is evolved using the equation:
36 
37  \f[
38  \laplacian \Phi = \div(\vec{U})
39  \f]
40 
41  Where:
42  \vartable
43  \Phi | Velocity potential [m2/s]
44  \vec{U} | Velocity [m/s]
45  \endvartable
46 
47  The corresponding pressure field could be calculated from the divergence
48  of the Euler equation:
49 
50  \f[
51  \laplacian p + \div(\div(\vec{U}\otimes\vec{U})) = 0
52  \f]
53 
54  but this generates excessive pressure variation in regions of large
55  velocity gradient normal to the flow direction. A better option is to
56  calculate the pressure field corresponding to velocity variation along the
57  stream-lines:
58 
59  \f[
60  \laplacian p + \div(\vec{F}\cdot\div(\vec{U}\otimes\vec{U})) = 0
61  \f]
62  where the flow direction tensor \f$\vec{F}\f$ is obtained from
63  \f[
64  \vec{F} = \hat{\vec{U}}\otimes\hat{\vec{U}}
65  \f]
66 
67  \heading Required fields
68  \plaintable
69  U | Velocity [m/s]
70  \endplaintable
71 
72  \heading Optional fields
73  \plaintable
74  p | Kinematic pressure [m2/s2]
75  Phi | Velocity potential [m2/s]
76  | Generated from p (if present) or U if not present
77  \endplaintable
78 
79  \heading Options
80  \plaintable
81  -writep | write the Euler pressure
82  -writePhi | Write the final velocity potential
83  -initialiseUBCs | Update the velocity boundaries before solving for Phi
84  \endplaintable
85 
86 \*---------------------------------------------------------------------------*/
87 
88 #include "fvCFD.H"
89 #include "pisoControl.H"
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 int main(int argc, char *argv[])
94 {
95  argList::addOption
96  (
97  "pName",
98  "pName",
99  "Name of the pressure field"
100  );
101 
102  argList::addBoolOption
103  (
104  "initialiseUBCs",
105  "Initialise U boundary conditions"
106  );
107 
108  argList::addBoolOption
109  (
110  "writePhi",
111  "Write the final velocity potential field"
112  );
113 
114  argList::addBoolOption
115  (
116  "writep",
117  "Calculate and write the Euler pressure field"
118  );
119 
120  argList::addBoolOption
121  (
122  "withFunctionObjects",
123  "execute functionObjects"
124  );
125 
126  #include "setRootCase.H"
127  #include "createTime.H"
128  #include "createMesh.H"
129 
130  pisoControl potentialFlow(mesh, "potentialFlow");
131 
132  #include "createFields.H"
133  #include "createMRF.H"
134 
135  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137  Info<< nl << "Calculating potential flow" << endl;
138 
139  // Since solver contains no time loop it would never execute
140  // function objects so do it ourselves
141  runTime.functionObjects().start();
142 
143  MRF.makeRelative(phi);
144  adjustPhi(phi, U, p);
145 
146  // Non-orthogonal velocity potential corrector loop
147  while (potentialFlow.correctNonOrthogonal())
148  {
149  fvScalarMatrix PhiEqn
150  (
152  ==
153  fvc::div(phi)
154  );
155 
156  PhiEqn.setReference(PhiRefCell, PhiRefValue);
157  PhiEqn.solve();
158 
159  if (potentialFlow.finalNonOrthogonalIter())
160  {
161  phi -= PhiEqn.flux();
162  }
163  }
164 
165  MRF.makeAbsolute(phi);
166 
167  Info<< "Continuity error = "
168  << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
169  << endl;
170 
172  U.correctBoundaryConditions();
173 
174  Info<< "Interpolated velocity error = "
175  << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
176  /sum(mesh.magSf())).value()
177  << endl;
178 
179  // Write U and phi
180  U.write();
181  phi.write();
182 
183  // Optionally write Phi
184  if (args.optionFound("writePhi"))
185  {
186  Phi.write();
187  }
188 
189  // Calculate the pressure field from the Euler equation
190  if (args.optionFound("writep"))
191  {
192  Info<< nl << "Calculating approximate pressure field" << endl;
193 
194  label pRefCell = 0;
195  scalar pRefValue = 0.0;
196  setRefCell
197  (
198  p,
199  potentialFlow.dict(),
200  pRefCell,
201  pRefValue
202  );
203 
204  // Calculate the flow-direction filter tensor
205  volScalarField magSqrU(magSqr(U));
206  volSymmTensorField F(sqr(U)/(magSqrU + SMALL*average(magSqrU)));
207 
208  // Calculate the divergence of the flow-direction filtered div(U*U)
209  // Filtering with the flow-direction generates a more reasonable
210  // pressure distribution in regions of high velocity gradient in the
211  // direction of the flow
212  volScalarField divDivUU
213  (
214  fvc::div
215  (
216  F & fvc::div(phi, U),
217  "div(div(phi,U))"
218  )
219  );
220 
221  // Solve a Poisson equation for the approximate pressure
222  while (potentialFlow.correctNonOrthogonal())
223  {
224  fvScalarMatrix pEqn
225  (
226  fvm::laplacian(p) + divDivUU
227  );
228 
229  pEqn.setReference(pRefCell, pRefValue);
230  pEqn.solve();
231  }
232 
233  p.write();
234  }
235 
236  runTime.functionObjects().end();
237 
238  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
239  << " ClockTime = " << runTime.elapsedClockTime() << " s"
240  << nl << endl;
241 
242  Info<< "End\n" << endl;
243 
244  return 0;
245 }
246 
247 
248 // ************************************************************************* //
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:54
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:27
pisoControl.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
p
p
Definition: pEqn.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvc::average
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
Foam::constant::physicoChemical::F
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
U
U
Definition: pEqn.H:46
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::label
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
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::setRefCell
void setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:31
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:15
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
createMesh.H
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
fvCFD.H
main
int main(int argc, char *argv[])
Definition: potentialFoam.C:126
args
Foam::argList args(argc, argv)
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)