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 | Copyright (C) 2011-2015 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  nonNewtonianIcoFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Transient solver for incompressible, laminar flow of non-Newtonian fluids.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
37 #include "pisoControl.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 int main(int argc, char *argv[])
42 {
43  #include "setRootCase.H"
44  #include "createTime.H"
45  #include "createMeshNoClear.H"
46 
47  pisoControl piso(mesh);
48 
49  #include "createFields.H"
50  #include "initContinuityErrs.H"
51 
52  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54  Info<< "\nStarting time loop\n" << endl;
55 
56  while (runTime.loop())
57  {
58  Info<< "Time = " << runTime.timeName() << nl << endl;
59 
60  #include "CourantNo.H"
61 
62  fluid.correct();
63 
64  // Momentum predictor
65 
67  (
68  fvm::ddt(U)
69  + fvm::div(phi, U)
70  - fvm::laplacian(fluid.nu(), U)
71  - (fvc::grad(U) & fvc::grad(fluid.nu()))
72  );
73 
74  if (piso.momentumPredictor())
75  {
76  solve(UEqn == -fvc::grad(p));
77  }
78 
79  // --- PISO loop
80  while (piso.correct())
81  {
82  volScalarField rAU(1.0/UEqn.A());
83 
84  volVectorField HbyA("HbyA", U);
85  HbyA = rAU*UEqn.H();
87  (
88  "phiHbyA",
89  (fvc::interpolate(HbyA) & mesh.Sf())
91  );
92 
93  adjustPhi(phiHbyA, U, p);
94 
95  // Non-orthogonal pressure corrector loop
96  while (piso.correctNonOrthogonal())
97  {
98  // Pressure corrector
99 
100  fvScalarMatrix pEqn
101  (
103  );
104 
105  pEqn.setReference(pRefCell, pRefValue);
106 
107  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
108 
109  if (piso.finalNonOrthogonalIter())
110  {
111  phi = phiHbyA - pEqn.flux();
112  }
113  }
114 
115  #include "continuityErrs.H"
116 
117  U = HbyA - rAU*fvc::grad(p);
118  U.correctBoundaryConditions();
119  }
120 
121  runTime.write();
122 
123  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
124  << " ClockTime = " << runTime.elapsedClockTime() << " s"
125  << nl << endl;
126  }
127 
128  Info<< "End\n" << endl;
129 
130  return 0;
131 }
132 
133 
134 // ************************************************************************* //
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:27
pisoControl.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
createMeshNoClear.H
singlePhaseTransportModel.H
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
U
U
Definition: pEqn.H:46
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
setRootCase.H
HbyA
HbyA
Definition: pEqn.H:4
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createTime.H
fvCFD.H
phiHbyA
phiHbyA
Definition: pEqn.H:21
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28