icoFoam.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-2023 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  icoFoam
26 
27 Description
28  Transient solver for incompressible, laminar flow of Newtonian fluids.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "pisoControl.H"
34 #include "pressureReference.H"
35 #include "findRefCell.H"
36 #include "constrainPressure.H"
37 #include "constrainHbyA.H"
38 #include "adjustPhi.H"
39 
40 #include "fvcDdt.H"
41 #include "fvcGrad.H"
42 #include "fvcFlux.H"
43 
44 #include "fvmDdt.H"
45 #include "fvmDiv.H"
46 #include "fvmLaplacian.H"
47 
48 using namespace Foam;
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  #include "setRootCase.H"
55  #include "createTime.H"
56  #include "createMesh.H"
57 
58  pisoControl piso(mesh);
59 
60  #include "createFields.H"
61  #include "initContinuityErrs.H"
62 
63  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65  Info<< "\nStarting time loop\n" << endl;
66 
67  while (runTime.loop())
68  {
69  Info<< "Time = " << runTime.userTimeName() << nl << endl;
70 
71  #include "CourantNo.H"
72 
73  // Momentum predictor
74 
76  (
77  fvm::ddt(U)
78  + fvm::div(phi, U)
79  - fvm::laplacian(nu, U)
80  );
81 
82  if (piso.momentumPredictor())
83  {
84  solve(UEqn == -fvc::grad(p));
85  }
86 
87  // --- PISO loop
88  while (piso.correct())
89  {
90  volScalarField rAU(1.0/UEqn.A());
93  (
94  "phiHbyA",
97  );
98 
99  adjustPhi(phiHbyA, U, p);
100 
101  // Update the pressure BCs to ensure flux consistency
103 
104  // Non-orthogonal pressure corrector loop
105  while (piso.correctNonOrthogonal())
106  {
107  // Pressure corrector
108 
109  fvScalarMatrix pEqn
110  (
112  );
113 
114  pEqn.setReference(pRefCell, pRefValue);
115 
116  pEqn.solve();
117 
119  {
120  phi = phiHbyA - pEqn.flux();
121  }
122  }
123 
124  #include "continuityErrs.H"
125 
126  U = HbyA - rAU*fvc::grad(p);
127  U.correctBoundaryConditions();
128  }
129 
130  runTime.write();
131 
132  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133  << " ClockTime = " << runTime.elapsedClockTime() << " s"
134  << nl << endl;
135  }
136 
137  Info<< "End\n" << endl;
138 
139  return 0;
140 }
141 
142 
143 // ************************************************************************* //
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Generic GeometricField class.
bool momentumPredictor() const
Flag to indicate to solve for momentum.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:809
tmp< VolField< Type > > H() const
Return the H operation source.
Definition: fvMatrix.C:867
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:963
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:563
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
Piso control class. Provides time-loop and piso-loop control methods. No convergence checking is done...
Definition: pisoControl.H:55
bool correct(const bool finalIter=false)
Piso loop within outer loop.
Definition: pisoControl.C:77
bool correctNonOrthogonal(const bool finalIter=true)
Non-orthogonal corrector loop.
Definition: pisoControl.C:98
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculates and prints the continuity errors.
pisoControl piso(mesh)
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Calculate the first temporal derivative.
Calculate the face-flux of the given field.
Calculate the gradient of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
int main(int argc, char *argv[])
Definition: icoFoam.C:49
Declare and initialise the cumulative continuity error.
volScalarField rAU(1.0/UEqn.A())
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho) *fvc::flux(HbyA))
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:35
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
static const char nl
Definition: Ostream.H:260
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField & p