initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / test / LduMatrix / LduMatrixTest.C
blob91a6dfb68d2dbf41dff4b12ef34cc2d7e9f8fc6e
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 \*---------------------------------------------------------------------------*/
27 #include "argList.H"
28 #include "Time.H"
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "LduMatrix.H"
32 #include "diagTensorField.H"
33 #include "TPCG.H"
34 #include "TPBiCG.H"
35 #include "NoPreconditioner.H"
37 using namespace Foam;
39 typedef Foam::LduMatrix<vector, diagTensor, scalar>
40     lduVectorMatrix;
41 defineNamedTemplateTypeNameAndDebug(lduVectorMatrix, 0);
43 typedef Foam::DiagonalSolver<vector, diagTensor, scalar>
44     lduVectorDiagonalSolver;
45 defineNamedTemplateTypeNameAndDebug(lduVectorDiagonalSolver, 0);
47 template<>
48 const vector lduVectorMatrix::great_(1e15, 1e15, 1e15);
50 template<>
51 const vector lduVectorMatrix::small_(1e-15, 1e-15, 1e-15);
53 namespace Foam
55     typedef LduMatrix<vector, diagTensor, scalar>::preconditioner
56         lduVectorPreconditioner;
57     defineTemplateRunTimeSelectionTable(lduVectorPreconditioner, symMatrix);
58     defineTemplateRunTimeSelectionTable(lduVectorPreconditioner, asymMatrix);
60     typedef LduMatrix<vector, diagTensor, scalar>::smoother
61         lduVectorSmoother;
62     defineTemplateRunTimeSelectionTable(lduVectorSmoother, symMatrix);
63     defineTemplateRunTimeSelectionTable(lduVectorSmoother, asymMatrix);
65     typedef LduMatrix<vector, diagTensor, scalar>::solver
66         lduVectorSolver;
67     defineTemplateRunTimeSelectionTable(lduVectorSolver, symMatrix);
68     defineTemplateRunTimeSelectionTable(lduVectorSolver, asymMatrix);
70     typedef TPCG<vector, diagTensor, scalar> TPCGVector;
71     defineNamedTemplateTypeNameAndDebug(TPCGVector, 0);
73     LduMatrix<vector, diagTensor, scalar>::solver::
74         addsymMatrixConstructorToTable<TPCGVector>
75         addTPCGSymMatrixConstructorToTable_;
77     typedef TPBiCG<vector, diagTensor, scalar> TPBiCGVector;
78     defineNamedTemplateTypeNameAndDebug(TPBiCGVector, 0);
80     LduMatrix<vector, diagTensor, scalar>::solver::
81         addasymMatrixConstructorToTable<TPBiCGVector>
82         addTPBiCGSymMatrixConstructorToTable_;
84     typedef NoPreconditioner<vector, diagTensor, scalar> NoPreconditionerVector;
85     defineNamedTemplateTypeNameAndDebug(NoPreconditionerVector, 0);
87     LduMatrix<vector, diagTensor, scalar>::preconditioner::
88         addsymMatrixConstructorToTable<NoPreconditionerVector>
89         addNoPreconditionerSymMatrixConstructorToTable_;
91     LduMatrix<vector, diagTensor, scalar>::preconditioner::
92         addasymMatrixConstructorToTable<NoPreconditionerVector>
93         addNoPreconditionerAsymMatrixConstructorToTable_;
97 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98 // Main program:
100 int main(int argc, char *argv[])
103 #   include "setRootCase.H"
105 #   include "createTime.H"
106 #   include "createMesh.H"
108     volVectorField psi
109     (
110         IOobject
111         (
112             "U",
113             runTime.timeName(),
114             mesh,
115             IOobject::MUST_READ,
116             IOobject::AUTO_WRITE
117         ),
118         mesh
119     );
121     lduVectorMatrix testMatrix(mesh);
122     testMatrix.diag() = 2*pTraits<diagTensor>::one;
123     testMatrix.source() = pTraits<vector>::one;
124     testMatrix.upper() = 0.1;
125     testMatrix.lower() = -0.1;
127     Info<< testMatrix << endl;
129     FieldField<Field, scalar> boundaryCoeffs(0);
130     FieldField<Field, scalar> internalCoeffs(0);
132     autoPtr<lduVectorMatrix::solver> testMatrixSolver =
133     lduVectorMatrix::solver::New
134     (
135         psi.name(),
136         testMatrix,
137         boundaryCoeffs,
138         internalCoeffs,
139         psi.boundaryField().interfaces(),
140         IStringStream
141         (
142             "PBiCG"
143             "{"
144             "    preconditioner   none;"
145             "    tolerance        (1e-05 1e-05 1e-05);"
146             "    relTol           (0 0 0);"
147             "}"
148         )()
149     );
151     lduVectorMatrix::solverPerformance solverPerf =
152         testMatrixSolver->solve(psi);
154     solverPerf.print();
156     Info<< psi << endl;
158     Info << "End\n" << endl;
160     return 0;
164 // ************************************************************************* //