/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2019 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . Class Foam::TDILUPreconditioner Description Simplified diagonal-based incomplete LU preconditioner for asymmetric matrices. The inverse (reciprocal for scalar) of the preconditioned diagonal is calculated and stored. SourceFiles TDILUPreconditioner.C \*---------------------------------------------------------------------------*/ #ifndef TDILUPreconditioner_H #define TDILUPreconditioner_H #include "LduMatrix.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { /*---------------------------------------------------------------------------*\ Class TDILUPreconditioner Declaration \*---------------------------------------------------------------------------*/ template class TDILUPreconditioner : public LduMatrix::preconditioner { // Private data //- The inverse (reciprocal for scalar) preconditioned diagonal Field rD_; public: //- Runtime type information TypeName("DILU"); // Constructors //- Construct from matrix components and preconditioner data dictionary TDILUPreconditioner ( const typename LduMatrix::solver& sol, const dictionary& preconditionerDict ); // Destructor virtual ~TDILUPreconditioner() = default; // Member Functions //- Calculate the reciprocal of the preconditioned diagonal static void calcInvD ( Field& rD, const LduMatrix& matrix ); //- Return wA the preconditioned form of residual rA virtual void precondition ( Field& wA, const Field& rA ) const; //- Return wT the transpose-matrix preconditioned form of // residual rT. virtual void preconditionT ( Field& wT, const Field& rT ) const; }; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace Foam // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #ifdef NoRepository #include "TDILUPreconditioner.C" #endif // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #endif // ************************************************************************* //