potentialFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 Description
28  Potential flow solver which solves for the velocity potential, to
29  calculate the flux-field, from which the velocity field is obtained by
30  reconstructing the flux.
31 
32  This application is particularly useful to generate starting fields for
33  Navier-Stokes codes.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "pisoControl.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
44  argList::addOption
45  (
46  "pName",
47  "pName",
48  "Name of the pressure field"
49  );
50 
51  argList::addBoolOption
52  (
53  "initialiseUBCs",
54  "Initialise U boundary conditions"
55  );
56 
57  argList::addBoolOption
58  (
59  "writePhi",
60  "Write the velocity potential field"
61  );
62 
63  argList::addBoolOption
64  (
65  "writep",
66  "Calculate and write the pressure field"
67  );
68 
69  argList::addBoolOption
70  (
71  "withFunctionObjects",
72  "execute functionObjects"
73  );
74 
75  #include "setRootCaseLists.H"
76  #include "createTime.H"
77  #include "createMesh.H"
78 
79  pisoControl potentialFlow(mesh, "potentialFlow");
80 
81  #include "createFields.H"
82 
83  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85  Info<< nl << "Calculating potential flow" << endl;
86 
87  // Since solver contains no time loop it would never execute
88  // function objects so do it ourselves
89  runTime.functionObjects().start();
90 
91  MRF.makeRelative(phi);
92  adjustPhi(phi, U, p);
93 
94  // Non-orthogonal velocity potential corrector loop
95  while (potentialFlow.correctNonOrthogonal())
96  {
97  fvScalarMatrix PhiEqn
98  (
100  ==
101  fvc::div(phi)
102  );
103 
104  PhiEqn.setReference(PhiRefCell, PhiRefValue);
105  PhiEqn.solve();
106 
107  if (potentialFlow.finalNonOrthogonalIter())
108  {
109  phi -= PhiEqn.flux();
110  }
111  }
112 
113  MRF.makeAbsolute(phi);
114 
115  Info<< "Continuity error = "
116  << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
117  << endl;
118 
120  U.correctBoundaryConditions();
121 
122  Info<< "Interpolated velocity error = "
123  << (sqrt(sum(sqr(fvc::flux(U) - phi)))/sum(mesh.magSf())).value()
124  << endl;
125 
126  // Write U and phi
127  U.write();
128  phi.write();
129 
130  // Optionally write Phi
131  if (args.optionFound("writePhi"))
132  {
133  Phi.write();
134  }
135 
136  // Calculate the pressure field
137  if (args.optionFound("writep"))
138  {
139  Info<< nl << "Calculating approximate pressure field" << endl;
140 
141  label pRefCell = 0;
142  scalar pRefValue = 0.0;
143  setRefCell
144  (
145  p,
146  potentialFlow.dict(),
147  pRefCell,
148  pRefValue
149  );
150 
151  // Calculate the flow-direction filter tensor
152  volScalarField magSqrU(magSqr(U));
153  volSymmTensorField F(sqr(U)/(magSqrU + small*average(magSqrU)));
154 
155  // Calculate the divergence of the flow-direction filtered div(U*U)
156  // Filtering with the flow-direction generates a more reasonable
157  // pressure distribution in regions of high velocity gradient in the
158  // direction of the flow
159  volScalarField divDivUU
160  (
161  fvc::div
162  (
163  F & fvc::div(phi, U),
164  "div(div(phi,U))"
165  )
166  );
167 
168  // Solve a Poisson equation for the approximate pressure
169  while (potentialFlow.correctNonOrthogonal())
170  {
171  fvScalarMatrix pEqn
172  (
173  fvm::laplacian(p) + divDivUU
174  );
175 
176  pEqn.setReference(pRefCell, pRefValue);
177  pEqn.solve();
178  }
179 
180  p.write();
181  }
182 
183  runTime.functionObjects().end();
184 
185  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
186  << " ClockTime = " << runTime.elapsedClockTime() << " s"
187  << nl << endl;
188 
189  Info<< "End\n" << endl;
190 
191  return 0;
192 }
193 
194 
195 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
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
surfaceScalarField & phi
IOMRFZoneList & MRF
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
const dictionary & potentialFlow(mesh.solutionDict().subDict("potentialFlow"))
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
volVectorField F(fluid.F())
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:265
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label pRefCell
Definition: createFields.H:106
messageStream Info
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
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:31
scalar pRefValue
Definition: createFields.H:107