nonNewtonianIcoFoam.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) 2011-2016 OpenFOAM Foundation
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  nonNewtonianIcoFoam
28 
29 Group
30  grpIncompressibleSolvers
31 
32 Description
33  Transient solver for incompressible, laminar flow of non-Newtonian fluids.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
39 #include "pisoControl.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  argList::addNote
46  (
47  "Transient solver for incompressible laminar flow"
48  " of non-Newtonian fluids."
49  );
50 
51  #include "postProcess.H"
52 
53  #include "addCheckCaseOptions.H"
54  #include "setRootCaseLists.H"
55  #include "createTime.H"
56  #include "createMeshNoClear.H"
57  #include "createControl.H"
58  #include "createFields.H"
59  #include "initContinuityErrs.H"
60 
61  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63  Info<< "\nStarting time loop\n" << endl;
64 
65  while (runTime.loop())
66  {
67  Info<< "Time = " << runTime.timeName() << nl << endl;
68 
69  #include "CourantNo.H"
70 
71  fluid.correct();
72 
73  // Momentum predictor
74 
76  (
77  fvm::ddt(U)
78  + fvm::div(phi, U)
79  - fvm::laplacian(fluid.nu(), U)
80  - (fvc::grad(U) & fvc::grad(fluid.nu()))
81  );
82 
83  if (piso.momentumPredictor())
84  {
85  solve(UEqn == -fvc::grad(p));
86  }
87 
88  // --- PISO loop
89  while (piso.correct())
90  {
91  volScalarField rAU(1.0/UEqn.A());
94  (
95  "phiHbyA",
98  );
99 
100  adjustPhi(phiHbyA, U, p);
101 
102  // Update the pressure BCs to ensure flux consistency
104 
105  // Non-orthogonal pressure corrector loop
106  while (piso.correctNonOrthogonal())
107  {
108  // Pressure corrector
109 
110  fvScalarMatrix pEqn
111  (
113  );
114 
115  pEqn.setReference(pRefCell, pRefValue);
116 
117  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
118 
119  if (piso.finalNonOrthogonalIter())
120  {
121  phi = phiHbyA - pEqn.flux();
122  }
123  }
124 
125  #include "continuityErrs.H"
126 
127  U = HbyA - rAU*fvc::grad(p);
128  U.correctBoundaryConditions();
129  }
130 
131  runTime.write();
132 
133  runTime.printExecutionTime(Info);
134  }
135 
136  Info<< "End\n" << endl;
137 
138  return 0;
139 }
140 
141 
142 // ************************************************************************* //
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:39
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:35
pisoControl.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
createMeshNoClear.H
Required Variables.
singlePhaseTransportModel.H
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
HbyA
HbyA
Definition: pcEqn.H:74
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
solve
CEqn solve()
piso
pisoControl piso(mesh)
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:47
phiHbyA
phiHbyA
Definition: pcEqn.H:73
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
ddtCorr
ddtCorr
Definition: readControls.H:9
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
createTime.H
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvCFD.H
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:36
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:37