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-2024 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 "argList.H"
41 #include "findRefCell.H"
42 #include "IOMRFZoneList.H"
43 #include "adjustPhi.H"
44 
45 #include "fvcFlux.H"
46 #include "fvcReconstruct.H"
47 
48 #include "fvmLaplacian.H"
49 
50 using namespace Foam;
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
57  (
58  "pName",
59  "pName",
60  "Name of the pressure field"
61  );
62 
64  (
65  "initialiseUBCs",
66  "Initialise U boundary conditions"
67  );
68 
70  (
71  "writePhi",
72  "Write the velocity potential field"
73  );
74 
76  (
77  "writep",
78  "Calculate and write the pressure field"
79  );
80 
82 
83  #include "setRootCase.H"
84  #include "createTime.H"
85  #include "createMesh.H"
86 
87  nonOrthogonalSolutionControl potentialFlow(mesh, "potentialFlow");
88 
89  #include "createFields.H"
90 
91  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93  Info<< nl << "Calculating potential flow" << endl;
94 
95  // Since solver contains no time loop it would never execute
96  // function objects so do it ourselves
97  runTime.functionObjects().start();
98 
99  MRF.makeRelative(phi);
100  adjustPhi(phi, U, p);
101 
102  // Non-orthogonal velocity potential corrector loop
103  while (potentialFlow.correctNonOrthogonal())
104  {
105  fvScalarMatrix PhiEqn
106  (
108  ==
109  fvc::div(phi)
110  );
111 
112  PhiEqn.setReference(PhiRefCell, PhiRefValue);
113  PhiEqn.solve();
114 
115  if (potentialFlow.finalNonOrthogonalIter())
116  {
117  phi -= PhiEqn.flux();
118  }
119  }
120 
121  Info<< "Continuity error = "
122  << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
123  << endl;
124 
125  U = fvc::reconstruct(MRF.absolute(phi));
126  U.correctBoundaryConditions();
127 
128  Info<< "Interpolated velocity error = "
129  << (
130  sqrt(sum(sqr(fvc::flux(U) - MRF.absolute(phi))))
131  /sum(mesh.magSf())
132  ).value()
133  << endl;
134 
135  // Write U and phi
136  U.write();
137  phi.write();
138 
139  // Optionally write Phi
140  if (args.optionFound("writePhi"))
141  {
142  Phi.write();
143  }
144 
145  // Calculate the pressure field
146  if (args.optionFound("writep"))
147  {
148  Info<< nl << "Calculating approximate pressure field" << endl;
149 
150  label pRefCell = 0;
151  scalar pRefValue = 0.0;
152  setRefCell
153  (
154  p,
155  potentialFlow.dict(),
156  pRefCell,
157  pRefValue
158  );
159 
160  // Calculate the flow-direction filter tensor
161  volScalarField magSqrU(magSqr(U));
162  volSymmTensorField F(sqr(U)/(magSqrU + small*average(magSqrU)));
163 
164  // Calculate the divergence of the flow-direction filtered div(U*U)
165  // Filtering with the flow-direction generates a more reasonable
166  // pressure distribution in regions of high velocity gradient in the
167  // direction of the flow
168  volScalarField divDivUU
169  (
170  fvc::div
171  (
172  F & fvc::div(phi, U),
173  "div(div(phi,U))"
174  )
175  );
176 
177  // Solve a Poisson equation for the approximate pressure
178  while (potentialFlow.correctNonOrthogonal())
179  {
180  fvScalarMatrix pEqn
181  (
182  fvm::laplacian(p) + divDivUU
183  );
184 
185  pEqn.setReference(pRefCell, pRefValue);
186  pEqn.solve();
187  }
188 
189  p.write();
190  }
191 
192  runTime.functionObjects().end();
193 
194  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
195  << " ClockTime = " << runTime.elapsedClockTime() << " s"
196  << nl << endl;
197 
198  Info<< "End\n" << endl;
199 
200  return 0;
201 }
202 
203 
204 // ************************************************************************* //
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Generic GeometricField class.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:291
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:157
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Non-orthogonal solution control class. Provides non-orthogonal-loop control methods.
virtual bool write(const bool write=true) const
Write using setting from DB.
IOMRFZoneList MRF(mesh)
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Calculate the face-flux of the given field.
Reconstruct volField from a face flux field.
Calculate the matrix for the laplacian of the field.
U
Definition: pEqn.H:72
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:35
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimless
messageStream Info
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
const dictionary & potentialFlow(mesh.solution().dict().subDict("potentialFlow"))
Foam::argList args(argc, argv)
volScalarField & p