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 volumetric flux
87  -writePhi | Write the final velocity potential
88  -initialiseUBCs | Update the velocity boundaries before solving for Phi
89  \endplaintable
90 
91 \*---------------------------------------------------------------------------*/
92 
93 #include "fvCFD.H"
94 #include "pisoControl.H"
95 #include "dynamicFvMesh.H"
96 #include "cellCellStencilObject.H"
97 #include "localMin.H"
98 
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101 int main(int argc, char *argv[])
102 {
103  argList::addNote
104  (
105  "Overset potential flow solver which solves for the velocity potential"
106  );
107 
108  argList::addOption
109  (
110  "pName",
111  "pName",
112  "Name of the pressure field"
113  );
114 
115  argList::addBoolOption
116  (
117  "initialiseUBCs",
118  "Initialise U boundary conditions"
119  );
120 
121  argList::addBoolOption
122  (
123  "writephi",
124  "Write the final volumetric flux field"
125  );
126 
127  argList::addBoolOption
128  (
129  "writePhi",
130  "Write the final velocity potential field"
131  );
132 
133  argList::addBoolOption
134  (
135  "writep",
136  "Calculate and write the Euler pressure field"
137  );
138 
139  argList::addBoolOption
140  (
141  "withFunctionObjects",
142  "Execute functionObjects"
143  );
144 
145  #include "addRegionOption.H"
146  #include "addCheckCaseOptions.H"
147  #include "setRootCaseLists.H"
148  #include "createTime.H"
149  #include "createNamedDynamicFvMesh.H"
150 
151  pisoControl potentialFlow(mesh, "potentialFlow");
152 
153  #include "createFields.H"
154 
155  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157  Info<< nl << "Calculating potential flow" << endl;
158 
159  mesh.update();
160 
161 
162  // Since solver contains no time loop it would never execute
163  // function objects so do it ourselves
164  runTime.functionObjects().start();
165 
166  MRF.makeRelative(phi);
167  adjustPhi(phi, U, p);
168 
169  // Non-orthogonal velocity potential corrector loop
170  while (potentialFlow.correct())
171  {
172  phi = fvc::flux(U);
173 
174  while (potentialFlow.correctNonOrthogonal())
175  {
176  fvScalarMatrix PhiEqn
177  (
179  ==
180  fvc::div(phi)
181  );
182 
183  PhiEqn.setReference(PhiRefCell, PhiRefValue);
184  PhiEqn.solve();
185 
186  if (potentialFlow.finalNonOrthogonalIter())
187  {
188  phi -= PhiEqn.flux();
189  }
190  }
191 
192  MRF.makeAbsolute(phi);
193 
194  Info<< "Continuity error = "
195  << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
196  << endl;
197 
199  U.correctBoundaryConditions();
200 
201  Info<< "Interpolated velocity error = "
202  << (sqrt(sum(sqr(fvc::flux(U) - phi)))/sum(mesh.magSf())).value()
203  << endl;
204  }
205 
206  // Write U
207  U.write();
208 
209  // Optionally write the volumetric flux, phi
210  if (args.found("writephi"))
211  {
212  phi.write();
213  }
214 
215  // Optionally write velocity potential, Phi
216  if (args.found("writePhi"))
217  {
218  Phi.write();
219  }
220 
221  // Calculate the pressure field from the Euler equation
222  if (args.found("writep"))
223  {
224  Info<< nl << "Calculating approximate pressure field" << endl;
225 
226  label pRefCell = 0;
227  scalar pRefValue = 0.0;
228  setRefCell
229  (
230  p,
231  potentialFlow.dict(),
232  pRefCell,
233  pRefValue
234  );
235 
236  // Calculate the flow-direction filter tensor
237  volScalarField magSqrU(magSqr(U));
238  volSymmTensorField F(sqr(U)/(magSqrU + SMALL*average(magSqrU)));
239 
240  // Calculate the divergence of the flow-direction filtered div(U*U)
241  // Filtering with the flow-direction generates a more reasonable
242  // pressure distribution in regions of high velocity gradient in the
243  // direction of the flow
244  volScalarField divDivUU
245  (
246  fvc::div
247  (
248  F & fvc::div(phi, U),
249  "div(div(phi,U))"
250  )
251  );
252 
253  // Solve a Poisson equation for the approximate pressure
254  while (potentialFlow.correctNonOrthogonal())
255  {
256  fvScalarMatrix pEqn
257  (
258  fvm::laplacian(p) + divDivUU
259  );
260 
261  pEqn.setReference(pRefCell, pRefValue);
262  pEqn.solve();
263  }
264 
265  p.write();
266  }
267 
268  runTime.functionObjects().end();
269 
270  runTime.printExecutionTime(Info);
271 
272  Info<< "End\n" << endl;
273 
274  return 0;
275 }
276 
277 
278 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:84
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
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:50
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:531
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
Required Classes.
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 Classes.
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)
Required Classes.
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