initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / test / LduMatrix / LduMatrixTest3.C
blob4e9ab70b0c6deee130b6a8997417539eaaa2eedc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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     icoFoam
28 Description
29     Transient solver for incompressible, laminar flow of Newtonian fluids.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "LduMatrix.H"
35 #include "diagTensorField.H"
37 typedef LduMatrix<vector, scalar, scalar> lduVectorMatrix;
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 int main(int argc, char *argv[])
45 #   include "setRootCase.H"
47 #   include "createTime.H"
48 #   include "createMesh.H"
49 #   include "createFields.H"
50 #   include "initContinuityErrs.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54     Info<< "\nStarting time loop\n" << endl;
56     while (runTime.loop())
57     {
58         Info<< "Time = " << runTime.timeName() << nl << endl;
60 #       include "readPISOControls.H"
61 #       include "CourantNo.H"
63         fvVectorMatrix UEqn
64         (
65             fvm::ddt(U)
66           + fvm::div(phi, U)
67           - fvm::laplacian(nu, U)
68         );
70         fvVectorMatrix UEqnp(UEqn == -fvc::grad(p));
72         lduVectorMatrix U3Eqnp(mesh);
73         U3Eqnp.diag() = UEqnp.diag();
74         U3Eqnp.upper() = UEqnp.upper();
75         U3Eqnp.lower() = UEqnp.lower();
76         U3Eqnp.source() = UEqnp.source();
78         UEqnp.addBoundaryDiag(U3Eqnp.diag(), 0);
79         UEqnp.addBoundarySource(U3Eqnp.source(), false);
81         autoPtr<lduVectorMatrix::solver> U3EqnpSolver =
82             lduVectorMatrix::solver::New
83             (
84                 U.name(),
85                 U3Eqnp,
86                 dictionary
87                 (
88                     IStringStream
89                     (
90                         "{"
91                         "    solver           PBiCG;"
92                         "    preconditioner   DILU;"
93                         "    tolerance        (1e-13 1e-13 1e-13);"
94                         "    relTol           (0 0 0);"
95                         "}"
96                     )()
97                 )
98             );
100         U3EqnpSolver->solve(U).print(Info);
102         // --- PISO loop
104         for (int corr=0; corr<nCorr; corr++)
105         {
106             volScalarField rUA = 1.0/UEqn.A();
108             U = rUA*UEqn.H();
109             phi = (fvc::interpolate(U) & mesh.Sf()) 
110                 + fvc::ddtPhiCorr(rUA, U, phi);
112             adjustPhi(phi, U, p);
114             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
115             {
116                 fvScalarMatrix pEqn
117                 (
118                     fvm::laplacian(rUA, p) == fvc::div(phi)
119                 );
121                 pEqn.setReference(pRefCell, pRefValue);
122                 pEqn.solve();
124                 if (nonOrth == nNonOrthCorr)
125                 {
126                     phi -= pEqn.flux();
127                 }
128             }
130 #           include "continuityErrs.H"
132             U -= rUA*fvc::grad(p);
133             U.correctBoundaryConditions();
134         }
136         runTime.write();
138         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
139             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
140             << nl << endl;
141     }
143     Info<< "End\n" << endl;
145     return 0;
149 // ************************************************************************* //