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