HermesCommon 1.0
precond_ifpack.cpp
Go to the documentation of this file.
00001 // This file is part of HermesCommon
00002 //
00003 // Hermes2D is free software: you can redistribute it and/or modify
00004 // it under the terms of the GNU General Public License as published by
00005 // the Free Software Foundation, either version 2 of the License, or
00006 // (at your option) any later version.
00007 //
00008 // Hermes2D is distributed in the hope that it will be useful,
00009 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00010 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00011 // GNU General Public License for more details.
00012 //
00013 // You should have received a copy of the GNU General Public License
00014 // along with Hermes2D.  If not, see <http://www.gnu.org/licenses/>.
00018 #include "config.h"
00019 #ifdef HAVE_IFPACK
00020 #include "precond_ifpack.h"
00021 #include "Ifpack_PointRelaxation.h"
00022 #include "Ifpack_BlockRelaxation.h"
00023 #include "Ifpack_DenseContainer.h"
00024 #include "Ifpack_AdditiveSchwarz.h"
00025 #include "Ifpack_ILU.h"
00026 #include "Ifpack_ILUT.h"
00027 #include "Ifpack_IC.h"
00028 #include "Ifpack_ICT.h"
00029 #include "Ifpack_Graph_Epetra_CrsGraph.h"
00030 
00031 namespace Hermes
00032 {
00033   namespace Preconditioners
00034 {
00035     template<typename Scalar>
00036     IfpackPrecond<Scalar>::IfpackPrecond(const char *cls, const char *type)
00037     {
00038       this->prec = NULL;
00039       this->owner = true;
00040       this->mat = NULL;
00041 
00042       this->cls = cls;
00043       this->type = type;
00044     }
00045 
00046     template<typename Scalar>
00047     IfpackPrecond<Scalar>::IfpackPrecond(const char *cls, const char *type, int overlap)
00048     {
00049       this->prec = NULL;
00050       this->owner = true;
00051       this->mat = NULL;
00052 
00053       this->cls = cls;
00054       this->type = type;
00055       this->overlap = overlap;
00056     }
00057 
00058     template<typename Scalar>
00059     IfpackPrecond<Scalar>::IfpackPrecond(Ifpack_Preconditioner *ipc)
00060     {
00061       this->prec = ipc;
00062       this->owner = false;
00063       this->mat = NULL;    // FIXME: take the matrix from ipc
00064     }
00065 
00066     template<typename Scalar>
00067     IfpackPrecond<Scalar>::~IfpackPrecond()
00068     {
00069       if(owner) delete prec;
00070     }
00071 
00072     template<typename Scalar>
00073     void IfpackPrecond<Scalar>::set_param(const char *name, const char *value)
00074     {
00075       ilist.set(name, value);
00076     }
00077 
00078     template<typename Scalar>
00079     void IfpackPrecond<Scalar>::set_param(const char *name, int value)
00080     {
00081       ilist.set(name, value);
00082     }
00083 
00084     template<typename Scalar>
00085     void IfpackPrecond<Scalar>::set_param(const char *name, double value)
00086     {
00087       ilist.set(name, value);
00088     }
00089 
00090     template<typename Scalar>
00091     void IfpackPrecond<Scalar>::create(Matrix<Scalar> *m)
00092     {
00093       EpetraMatrix<Scalar> *mt = static_cast<EpetraMatrix<Scalar> *>(m);
00094       assert(mt != NULL);
00095       mat = mt;
00096       if(strcmp(cls, "point-relax") == 0)
00097       {
00098         create_point_relax(mat, type);
00099         apply_params();
00100         initialize();
00101       }
00102       else if(strcmp(cls, "block-relax") == 0)
00103       {
00104         create_block_relax(mat, type);
00105         apply_params();
00106       }
00107       else if(strcmp(cls, "add-schwartz") == 0)
00108       {
00109         create_add_schwartz(mat, type, overlap);
00110         apply_params();
00111         initialize();
00112       }
00113     }
00114 
00115     template<typename Scalar>
00116     void IfpackPrecond<Scalar>::create_point_relax(EpetraMatrix<Scalar> *a, const char *name)
00117     {
00118       prec = new Ifpack_PointRelaxation(a->mat);
00119       ilist.set("relaxation: type", name);
00120     }
00121 
00122     template<typename Scalar>
00123     void IfpackPrecond<Scalar>::create_block_relax(EpetraMatrix<Scalar> *a, const char *name)
00124     {
00125       Teuchos::RCP<const Epetra_CrsGraph> rgraph = Teuchos::rcp(&(a->mat->Graph()));
00126       Ifpack_Graph *graph = new Ifpack_Graph_Epetra_CrsGraph(rgraph);
00127 
00128       Ifpack_Partitioner *partitioner = new Ifpack_GreedyPartitioner(graph);
00129 
00130       Teuchos::ParameterList list;
00131       list.set("partitioner: local parts", 1000);  //\todo parametrize me
00132       partitioner->SetParameters(list);
00133       partitioner->Compute();
00134 
00135       prec = new Ifpack_BlockRelaxation<Ifpack_DenseContainer>(a->mat);
00136       ilist.set("relaxation: type", name);
00137 
00138       rgraph.release();
00139     }
00140 
00141     template<typename Scalar>
00142     void IfpackPrecond<Scalar>::create_add_schwartz(EpetraMatrix<Scalar> *a, const char *name, int overlap)
00143     {
00144       if(strcasecmp(name, "ilu") == 0)
00145       {
00146         prec = new Ifpack_AdditiveSchwarz<Ifpack_ILU>(a->mat, overlap);
00147       }
00148       else if(strcasecmp(name, "ilut") == 0)
00149       {
00150         prec = new Ifpack_AdditiveSchwarz<Ifpack_ILUT>(a->mat, overlap);
00151       }
00152       else if(strcasecmp(name, "ic") == 0)
00153       {
00154         prec = new Ifpack_AdditiveSchwarz<Ifpack_IC>(a->mat, overlap);
00155       }
00156       else if(strcasecmp(name, "ict") == 0)
00157       {
00158         prec = new Ifpack_AdditiveSchwarz<Ifpack_ICT>(a->mat, overlap);
00159       }
00160       else
00161         prec = NULL;
00162     }
00163 
00164     template<typename Scalar>
00165     int IfpackPrecond<Scalar>::initialize()
00166     {
00167       assert(prec != NULL);
00168       return prec->Initialize();
00169     }
00170 
00171     template<typename Scalar>
00172     void IfpackPrecond<Scalar>::compute()
00173     {
00174       assert(prec != NULL);
00175       prec->Compute();
00176     }
00177 
00178     template<typename Scalar>
00179     void IfpackPrecond<Scalar>::apply_params()
00180     {
00181       prec->SetParameters(ilist);
00182     }
00183 
00184     template<typename Scalar>
00185     int IfpackPrecond<Scalar>::ApplyInverse(const Epetra_MultiVector &r, Epetra_MultiVector &z) const
00186     {
00187       assert(prec != NULL);
00188       return prec->ApplyInverse(r, z);
00189     }
00190 
00191     template<typename Scalar>
00192     const Epetra_Comm &IfpackPrecond<Scalar>::Comm() const
00193     {
00194       return mat->mat->Comm();
00195     }
00196 
00197     template<typename Scalar>
00198     const Epetra_Map &IfpackPrecond<Scalar>::OperatorDomainMap() const
00199     {
00200       return mat->mat->OperatorDomainMap();
00201     }
00202 
00203     template<typename Scalar>
00204     const Epetra_Map &IfpackPrecond<Scalar>::OperatorRangeMap() const
00205     {
00206       return mat->mat->OperatorRangeMap();
00207     }
00208     template class HERMES_API IfpackPrecond<double>;
00209     template class HERMES_API IfpackPrecond<std::complex<double> >;
00210   }
00211 }
00212 #endif