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-------------------------------------------------------------------------------
10License
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
26Application
27 nonNewtonianIcoFoam
28
29Group
30 grpIncompressibleSolvers
31
32Description
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
43int 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",
96 fvc::flux(HbyA)
97 + fvc::interpolate(rAU)*fvc::ddtCorr(U, phi)
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 (
112 fvm::laplacian(rAU, p) == fvc::div(phiHbyA)
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// ************************************************************************* //
Required Classes.
surfaceScalarField & phi
const scalar pRefValue
const label pRefCell
twoPhaseSystem & fluid
U
Definition: pEqn.H:72
volScalarField & p
fvVectorMatrix & UEqn
Definition: UEqn.H:13
phiHbyA
Definition: pcEqn.H:73
HbyA
Definition: pcEqn.H:74
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
pisoControl piso(mesh)
adjustPhi(phiHbyA, U, p_rgh)
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:46
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
CEqn solve()