Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / incompressible / simpleSRFFoam / simpleSRFFoam.C
blob4e8768b3096450169972598ba4bacb4be0ffa32a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     simpleSRFFoam
28 Description
29     Steady-state solver for incompressible, turbulent flow of non-Newtonian
30     fluids with single rotating frame.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
36 #include "incompressible/RAS/RASModel/RASModel.H"
37 #include "SRFModel.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
44 #   include "setRootCase.H"
46 #   include "createTime.H"
47 #   include "createMesh.H"
48 #   include "createFields.H"
49 #   include "initContinuityErrs.H"
51     //mesh.clearPrimitives();
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55     Info<< "\nStarting time loop\n" << endl;
57     for (runTime++; !runTime.end(); runTime++)
58     {
59         Info<< "Time = " << runTime.timeName() << nl << endl;
61 #       include "readSIMPLEControls.H"
63         p.storePrevIter();
65         // Pressure-velocity SIMPLE corrector
66         {
67             // Momentum predictor
68             tmp<fvVectorMatrix> UrelEqn
69             (
70                 fvm::div(phi, Urel)
71               + turbulence->divDevReff(Urel)
72               + SRF->Su()
73             );
75             UrelEqn().relax();
77             solve(UrelEqn() == -fvc::grad(p));
79             p.boundaryField().updateCoeffs();
80             volScalarField AUrel = UrelEqn().A();
81             Urel = UrelEqn().H()/AUrel;
82             UrelEqn.clear();
83             phi = fvc::interpolate(Urel) & mesh.Sf();
84             adjustPhi(phi, Urel, p);
86             // Non-orthogonal pressure corrector loop
87             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
88             {
89                 fvScalarMatrix pEqn
90                 (
91                     fvm::laplacian(1.0/AUrel, p) == fvc::div(phi)
92                 );
94                 pEqn.setReference(pRefCell, pRefValue);
95                 pEqn.solve();
97                 if (nonOrth == nNonOrthCorr)
98                 {
99                     phi -= pEqn.flux();
100                 }
101             }
103 #           include "continuityErrs.H"
105             // Explicitly relax pressure for momentum corrector
106             p.relax();
108             // Momentum corrector
109             Urel -= fvc::grad(p)/AUrel;
110             Urel.correctBoundaryConditions();
111         }
113         turbulence->correct();
115         // Recalculate Uabs
116         Uabs = Urel + SRF->U();
118         runTime.write();
120         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
121             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
122             << nl << endl;
123     }
125     Info<< "End\n" << endl;
127     return(0);
131 // ************************************************************************* //