$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 #include "libmesh/libmesh_common.h" 00019 00020 #ifdef LIBMESH_HAVE_TRILINOS 00021 00022 // C++ includes 00023 00024 // Local Includes 00025 #include "libmesh/trilinos_preconditioner.h" 00026 #include "libmesh/trilinos_epetra_matrix.h" 00027 #include "libmesh/trilinos_epetra_vector.h" 00028 #include "libmesh/libmesh_common.h" 00029 00030 #include "Ifpack.h" 00031 #include "Ifpack_DiagPreconditioner.h" 00032 #include "Ifpack_AdditiveSchwarz.h" 00033 #include "Ifpack_ILU.h" 00034 #include "Ifpack_ILUT.h" 00035 #include "Ifpack_IC.h" 00036 #include "Ifpack_ICT.h" 00037 00038 #ifdef LIBMESH_HAVE_ML 00039 #include "ml_MultiLevelPreconditioner.h" 00040 #endif 00041 00042 namespace libMesh 00043 { 00044 00045 template <typename T> 00046 void TrilinosPreconditioner<T>::apply(const NumericVector<T> & /* x */, 00047 NumericVector<T> & /* y */ ) 00048 { 00049 } 00050 00051 00052 00053 00054 template <typename T> 00055 void TrilinosPreconditioner<T>::init () 00056 { 00057 if(!this->_matrix) 00058 libmesh_error_msg("ERROR: No matrix set for PetscPreconditioner, but init() called"); 00059 00060 // Clear the preconditioner in case it has been created in the past 00061 if (!this->_is_initialized) 00062 { 00063 EpetraMatrix<T> * matrix = cast_ptr<EpetraMatrix<T>*, SparseMatrix<T> >(this->_matrix); 00064 _mat = matrix->mat(); 00065 } 00066 00067 set_preconditioner_type(this->_preconditioner_type); 00068 00069 this->_is_initialized = true; 00070 } 00071 00072 00073 template <typename T> 00074 void 00075 TrilinosPreconditioner<T>::set_params(Teuchos::ParameterList & list) 00076 { 00077 _param_list = list; 00078 } 00079 00080 00081 template <typename T> 00082 void 00083 TrilinosPreconditioner<T>::compute() 00084 { 00085 Ifpack_Preconditioner * ifpack = NULL; 00086 #ifdef LIBMESH_HAVE_ML 00087 ML_Epetra::MultiLevelPreconditioner * ml = NULL; 00088 #endif 00089 00090 switch (this->_preconditioner_type) 00091 { 00092 // IFPACK preconditioners 00093 case ILU_PRECOND: 00094 case SOR_PRECOND: 00095 ifpack = dynamic_cast<Ifpack_Preconditioner *>(_prec); 00096 ifpack->Compute(); 00097 break; 00098 00099 #ifdef LIBMESH_HAVE_ML 00100 // ML preconditioners 00101 case AMG_PRECOND: 00102 ml = dynamic_cast<ML_Epetra::MultiLevelPreconditioner *>(_prec); 00103 ml->ComputePreconditioner(); 00104 break; 00105 #endif 00106 00107 default: 00108 // no nothing here 00109 break; 00110 } 00111 } 00112 00113 00114 template <typename T> 00115 void 00116 TrilinosPreconditioner<T>::set_preconditioner_type (const PreconditionerType & preconditioner_type) 00117 { 00118 Ifpack_Preconditioner * pc = NULL; 00119 #ifdef LIBMESH_HAVE_ML 00120 ML_Epetra::MultiLevelPreconditioner * ml = NULL; 00121 #endif 00122 00123 switch (preconditioner_type) 00124 { 00125 case IDENTITY_PRECOND: 00126 // pc = new Ifpack_DiagPreconditioner(); 00127 break; 00128 00129 case CHOLESKY_PRECOND: 00130 break; 00131 00132 case ICC_PRECOND: 00133 break; 00134 00135 case ILU_PRECOND: 00136 pc = new Ifpack_ILU(_mat); 00137 pc->SetParameters(_param_list); 00138 pc->Initialize(); 00139 _prec = pc; 00140 break; 00141 00142 case LU_PRECOND: 00143 break; 00144 00145 case ASM_PRECOND: 00146 break; 00147 00148 case JACOBI_PRECOND: 00149 break; 00150 00151 case BLOCK_JACOBI_PRECOND: 00152 break; 00153 00154 case SOR_PRECOND: 00155 break; 00156 00157 case EISENSTAT_PRECOND: 00158 break; 00159 00160 #ifdef LIBMESH_HAVE_ML 00161 case AMG_PRECOND: 00162 ml = new ML_Epetra::MultiLevelPreconditioner(*_mat, _param_list, false);; 00163 _prec = ml; 00164 break; 00165 #endif 00166 00167 default: 00168 libMesh::err << "ERROR: Unsupported Trilinos Preconditioner: " 00169 << preconditioner_type << std::endl 00170 << "Continuing with Trilinos defaults" << std::endl; 00171 } 00172 00173 } 00174 00175 00176 template <typename T> 00177 int 00178 TrilinosPreconditioner<T>::SetUseTranspose(bool UseTranspose) 00179 { 00180 return _prec->SetUseTranspose(UseTranspose); 00181 } 00182 00183 template <typename T> 00184 int 00185 TrilinosPreconditioner<T>::Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const 00186 { 00187 return _prec->Apply(X, Y); 00188 } 00189 00190 template <typename T> 00191 int 00192 TrilinosPreconditioner<T>::ApplyInverse(const Epetra_MultiVector &r, Epetra_MultiVector &z) const 00193 { 00194 return _prec->ApplyInverse(r, z); 00195 } 00196 00197 template <typename T> 00198 double 00199 TrilinosPreconditioner<T>::NormInf() const 00200 { 00201 return _prec->NormInf(); 00202 } 00203 00204 template <typename T> 00205 const char * 00206 TrilinosPreconditioner<T>::Label() const 00207 { 00208 return _prec->Label(); 00209 } 00210 00211 template <typename T> 00212 bool 00213 TrilinosPreconditioner<T>::UseTranspose() const 00214 { 00215 return _prec->UseTranspose(); 00216 } 00217 00218 template <typename T> 00219 bool 00220 TrilinosPreconditioner<T>::HasNormInf() const 00221 { 00222 return _prec->HasNormInf(); 00223 } 00224 00225 template <typename T> 00226 const Epetra_Comm & 00227 TrilinosPreconditioner<T>::Comm() const 00228 { 00229 return _prec->Comm(); 00230 } 00231 00232 template <typename T> 00233 const Epetra_Map & 00234 TrilinosPreconditioner<T>::OperatorDomainMap() const 00235 { 00236 return _prec->OperatorDomainMap(); 00237 } 00238 00239 template <typename T> 00240 const Epetra_Map & 00241 TrilinosPreconditioner<T>::OperatorRangeMap() const 00242 { 00243 return _prec->OperatorRangeMap(); 00244 } 00245 00246 //------------------------------------------------------------------ 00247 // Explicit instantiations 00248 template class TrilinosPreconditioner<Number>; 00249 00250 } // namespace libMesh 00251 00252 #endif // #ifdef LIBMESH_HAVE_TRILINOS