FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 Usage
31  - nonNewtonianIcoFoam [OPTION]
32 
33  @param -case <dir> \n
34  Specify the case directory
35 
36  @param -parallel \n
37  Run the case in parallel
38 
39  @param -help \n
40  Display short usage message
41 
42  @param -doc \n
43  Display Doxygen documentation page
44 
45  @param -srcDoc \n
46  Display source code
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #include <finiteVolume/fvCFD.H>
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #include <OpenFOAM/setRootCase.H>
58 
59  #include <OpenFOAM/createTime.H>
61  #include "createFields.H"
63 
64  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66  Info<< "\nStarting time loop\n" << endl;
67 
68  while (runTime.loop())
69  {
70  Info<< "Time = " << runTime.timeName() << nl << endl;
71 
73  #include <finiteVolume/CourantNo.H>
74 
75  fluid.correct();
76 
78  (
79  fvm::ddt(U)
80  + fvm::div(phi, U)
81  - fvm::laplacian(fluid.nu(), U)
82  );
83 
84  solve(UEqn == -fvc::grad(p));
85 
86  // --- PISO loop
87 
88  for (int corr=0; corr<nCorr; corr++)
89  {
90  volScalarField rUA = 1.0/UEqn.A();
91 
92  U = rUA*UEqn.H();
93  phi = (fvc::interpolate(U) & mesh.Sf())
94  + fvc::ddtPhiCorr(rUA, U, phi);
95 
96  adjustPhi(phi, U, p);
97 
98  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
99  {
100  fvScalarMatrix pEqn
101  (
102  fvm::laplacian(rUA, p) == fvc::div(phi)
103  );
104 
105  pEqn.setReference(pRefCell, pRefValue);
106  pEqn.solve();
107 
108  if (nonOrth == nNonOrthCorr)
109  {
110  phi -= pEqn.flux();
111  }
112  }
113 
115 
116  U -= rUA*fvc::grad(p);
117  U.correctBoundaryConditions();
118  }
119 
120  runTime.write();
121 
122  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
123  << " ClockTime = " << runTime.elapsedClockTime() << " s"
124  << nl << endl;
125  }
126 
127  Info<< "End\n" << endl;
128 
129  return 0;
130 }
131 
132 
133 // ************************ vim: set sw=4 sts=4 et: ************************ //