nonNewtonianIcoFoam.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-2018 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 Description
28  Transient solver for incompressible, laminar flow of non-Newtonian fluids.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
34 #include "pisoControl.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 int main(int argc, char *argv[])
39 {
40  #include "postProcess.H"
41 
42  #include "setRootCaseLists.H"
43  #include "createTime.H"
44  #include "createMeshNoClear.H"
45  #include "createControl.H"
46  #include "createFields.H"
47  #include "initContinuityErrs.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< "\nStarting time loop\n" << endl;
52 
53  while (runTime.loop())
54  {
55  Info<< "Time = " << runTime.timeName() << nl << endl;
56 
57  #include "CourantNo.H"
58 
59  fluid.correct();
60 
61  // Momentum predictor
62 
64  (
65  fvm::ddt(U)
66  + fvm::div(phi, U)
67  - fvm::laplacian(fluid.nu(), U)
68  - (fvc::grad(U) & fvc::grad(fluid.nu()))
69  );
70 
71  if (piso.momentumPredictor())
72  {
73  solve(UEqn == -fvc::grad(p));
74  }
75 
76  // --- PISO loop
77  while (piso.correct())
78  {
79  volScalarField rAU(1.0/UEqn.A());
82  (
83  "phiHbyA",
86  );
87 
88  adjustPhi(phiHbyA, U, p);
89 
90  // Update the pressure BCs to ensure flux consistency
92 
93  // Non-orthogonal pressure corrector loop
94  while (piso.correctNonOrthogonal())
95  {
96  // Pressure corrector
97 
98  fvScalarMatrix pEqn
99  (
101  );
102 
103  pEqn.setReference(pRefCell, pRefValue);
104 
105  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
106 
107  if (piso.finalNonOrthogonalIter())
108  {
109  phi = phiHbyA - pEqn.flux();
110  }
111  }
112 
113  #include "continuityErrs.H"
114 
115  U = HbyA - rAU*fvc::grad(p);
116  U.correctBoundaryConditions();
117  }
118 
119  runTime.write();
120 
121  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
122  << " ClockTime = " << runTime.elapsedClockTime() << " s"
123  << nl << endl;
124  }
125 
126  Info<< "End\n" << endl;
127 
128  return 0;
129 }
130 
131 
132 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
surfaceScalarField & phi
multiphaseSystem & fluid
Definition: createFields.H:11
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
pisoControl piso(mesh)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
HbyA
Definition: pcEqn.H:74
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
phiHbyA
Definition: pcEqn.H:73
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
rhoEqn solve()
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
static const char nl
Definition: Ostream.H:265
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
U
Definition: pEqn.H:72
label pRefCell
Definition: createFields.H:106
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
scalar pRefValue
Definition: createFields.H:107