Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #include "TSFIfpackOperator.hpp"
00030 #include "Teuchos_Array.hpp"
00031 #include "Teuchos_MPIComm.hpp"
00032 #include "TSFEpetraVector.hpp"
00033
00034
00035 #ifndef HAVE_TEUCHOS_EXPLICIT_INSTANTIATION
00036 #include "TSFVectorImpl.hpp"
00037 #include "TSFLinearOperatorImpl.hpp"
00038 #endif
00039
00040 using namespace TSFExtended;
00041 using namespace Teuchos;
00042
00043 IfpackOperator::IfpackOperator(const EpetraMatrix* A,
00044 int fillLevels,
00045 int overlapFill,
00046 double relaxationValue,
00047 double relativeThreshold,
00048 double absoluteThreshold)
00049 : precondGraph_(),
00050 precond_(),
00051 domain_(A->domain()),
00052 range_(A->range())
00053 {
00054 const Epetra_CrsMatrix* matrix = A->crsMatrix();
00055
00056 const Epetra_CrsGraph& matrixGraph = matrix->Graph();
00057
00058 Ifpack_IlukGraph* precondGraph
00059 = new Ifpack_IlukGraph(matrixGraph, fillLevels, overlapFill);
00060 precondGraph_ = rcp(precondGraph);
00061
00062 int ierr = precondGraph->ConstructFilledGraph();
00063
00064 TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00065 "IfpackOperator ctor: "
00066 "precondGraph->ConstructFilledGraph() failed with ierr="
00067 << ierr);
00068
00069 Ifpack_CrsRiluk* precond = new Ifpack_CrsRiluk(*precondGraph);
00070 precond_ = rcp(precond);
00071
00072 precond->SetRelaxValue(relaxationValue);
00073 precond->SetRelativeThreshold(relativeThreshold);
00074 precond->SetAbsoluteThreshold(absoluteThreshold);
00075
00076 ierr = precond->InitValues(*matrix);
00077
00078 TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00079 "IfpackOperator ctor: "
00080 "precond->InitValues() failed with ierr="
00081 << ierr);
00082
00083 ierr = precond->Factor();
00084
00085 TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00086 "IfpackOperator ctor: "
00087 "precond->Factor() failed with ierr="
00088 << ierr);
00089 }
00090
00091 bool IfpackOperator::opSupportedImpl(Thyra::EOpTransp M_trans) const
00092 {
00093 return (M_trans == Thyra::NOTRANS);
00094 }
00095
00096 void IfpackOperator::generalApply(const Thyra::EOpTransp M_trans,
00097 const Thyra::VectorBase<double>& x,
00098 Thyra::VectorBase<double>* y,
00099 const double alpha,
00100 const double beta) const
00101 {
00102
00103
00104 const EpetraVector* epIn = dynamic_cast<const EpetraVector*>(&x);
00105 TEST_FOR_EXCEPTION(epIn == 0, std::runtime_error,
00106 "IfpackOperator apply: input vector is "
00107 "not an EpetraVector");
00108
00109 const Epetra_Vector* in = epIn->epetraVec().get();
00110
00111 EpetraVector* epy = dynamic_cast<EpetraVector*>(y);
00112 TEST_FOR_EXCEPTION(epy == 0, std::runtime_error,
00113 "IfpackOperator apply: output vector is "
00114 "not an EpetraVector");
00115
00116 Epetra_Vector* yy = epy->epetraVec().get();
00117
00118
00119
00120 RCP<Thyra::VectorBase<double> > tsfOut;
00121 Epetra_Vector* tmp;
00122
00123 if (beta!=0.0)
00124 {
00125 tsfOut = createMember(y->space());
00126 EpetraVector* ep = dynamic_cast<EpetraVector*>(tsfOut.get());
00127 tmp = ep->epetraVec().get();
00128 }
00129 else
00130 {
00131 tmp = yy;
00132 }
00133
00134
00135
00136 Ifpack_CrsRiluk* p = const_cast<Ifpack_CrsRiluk*>(precond_.get());
00137
00138 int ierr;
00139
00140
00141 if (M_trans==NOTRANS)
00142 {
00143 ierr = p->Solve(false, *in, *tmp);
00144 }
00145 else
00146 {
00147 ierr = p->Solve(true, *in, *tmp);
00148 }
00149
00150
00151 if (beta != 0.0)
00152 {
00153
00154 yy->Update(alpha, *tmp, beta);
00155 }
00156 else if (alpha != 1.0)
00157 {
00158
00159
00160
00161 tmp->Scale(alpha);
00162 }
00163
00164
00165
00166 }
00167