PlayaEpetraMatrix.cpp
Go to the documentation of this file.
00001 /* @HEADER@ */
00002 //   
00003  /* @HEADER@ */
00004 
00005 #include "PlayaEpetraMatrix.hpp"
00006 #include "PlayaEpetraVector.hpp"
00007 #include "PlayaVectorSpaceDecl.hpp"  // changed from Impl
00008 #include "PlayaVectorDecl.hpp"
00009 #include "PlayaLinearOperatorDecl.hpp"  // changed from Impl
00010 #include "PlayaIfpackILUOperator.hpp"
00011 #include "PlayaIfpackICCOperator.hpp"
00012 #include "PlayaPreconditioner.hpp"
00013 #include "PlayaGenericLeftPreconditioner.hpp"
00014 #include "PlayaGenericRightPreconditioner.hpp"
00015 
00016 #include "PlayaGenericTwoSidedPreconditioner.hpp"
00017 
00018 
00019 
00020 #ifndef HAVE_TEUCHOS_EXPLICIT_INSTANTIATION
00021 #include "PlayaVectorImpl.hpp"
00022 #include "PlayaLinearOperatorImpl.hpp"
00023 #endif
00024 
00025 using namespace Playa;
00026 using namespace Teuchos;
00027 
00028 
00029 EpetraMatrix::EpetraMatrix(const Epetra_CrsGraph& graph,
00030   const VectorSpace<double>& domain,
00031   const VectorSpace<double>& range)
00032   : LinearOpWithSpaces<double>(domain, range),
00033     matrix_(rcp(new Epetra_CrsMatrix(Copy, graph)))
00034 {}
00035 
00036 EpetraMatrix::EpetraMatrix(const RCP<Epetra_CrsMatrix>& mat,
00037   const VectorSpace<double>& domain,
00038   const VectorSpace<double>& range)
00039   : LinearOpWithSpaces<double>(domain, range),
00040     matrix_(mat)
00041 {}
00042 
00043 
00044 void EpetraMatrix::apply(
00045   Teuchos::ETransp applyType,
00046   const Vector<double>& in,
00047   Vector<double> out) const
00048 {
00049   const Epetra_Vector& epIn = EpetraVector::getConcrete(in);
00050   Epetra_Vector& epOut = EpetraVector::getConcrete(out);
00051   using Teuchos::rcp_dynamic_cast;
00052 
00053   bool trans = applyType == TRANS;
00054   
00055   int ierr = matrix_->Multiply(trans, epIn, epOut);
00056   TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error, 
00057     "EpetraMatrix::generalApply() detected ierr="
00058     << ierr << " in matrix_->Multiply()");
00059 }
00060 
00061 
00062 
00063 void EpetraMatrix::addToRow(int globalRowIndex,
00064   int nElemsToInsert,
00065   const int* globalColumnIndices,
00066   const double* elementValues)
00067 {
00068   int ierr = crsMatrix()->SumIntoGlobalValues(globalRowIndex,
00069     nElemsToInsert,
00070     (double*) elementValues,
00071     (int*) globalColumnIndices);
00072 
00073   TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error, 
00074     "failed to add to row " << globalRowIndex
00075     << " in EpetraMatrix::addToRow() with nnz="
00076     << nElemsToInsert 
00077     << ". Error code was " << ierr);
00078 }
00079 
00080 
00081 
00082 void EpetraMatrix::addToElementBatch(int numRows, 
00083   int rowBlockSize,
00084   const int* globalRowIndices,
00085   int numColumnsPerRow,
00086   const int* globalColumnIndices,
00087   const double* values,
00088   const int* skipRow)
00089 {
00090   Epetra_CrsMatrix* crs = crsMatrix();
00091 
00092   int numRowBlocks = numRows/rowBlockSize;
00093   int row = 0;
00094 
00095   for (int rb=0; rb<numRowBlocks; rb++)
00096   {
00097     const int* cols = globalColumnIndices + rb*numColumnsPerRow;
00098     for (int r=0; r<rowBlockSize; r++, row++)
00099     {
00100       if (skipRow[row]) continue;
00101       const double* rowVals = values + row*numColumnsPerRow;
00102       int ierr=crs->SumIntoGlobalValues(globalRowIndices[row], 
00103         numColumnsPerRow,
00104         (double*) rowVals,
00105         (int*) cols);
00106       TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error, 
00107         "failed to add to row " << globalRowIndices[row]
00108         << " in EpetraMatrix::addToRow() with nnz="
00109         << numColumnsPerRow
00110         << ". Error code was " << ierr);
00111     }
00112   }
00113 }
00114 
00115 
00116 void EpetraMatrix::zero()
00117 {
00118   crsMatrix()->PutScalar(0.0);
00119 }
00120 
00121 
00122 
00123 void EpetraMatrix::getILUKPreconditioner(int fillLevels,
00124   int overlapFill,
00125   double relaxationValue,
00126   double relativeThreshold,
00127   double absoluteThreshold,
00128   LeftOrRight leftOrRight,
00129   Preconditioner<double>& rtn) const
00130 {
00131   RCP<LinearOperatorBase<double> > a = rcp(new IfpackILUOperator(this, 
00132       fillLevels,
00133       overlapFill,
00134       relaxationValue,
00135       relativeThreshold,
00136       absoluteThreshold));
00137 
00138   LinearOperator<double> ilu = a;
00139 
00140 
00141   if (leftOrRight == Left)
00142   {
00143     rtn = new GenericLeftPreconditioner<double>(ilu);
00144   }
00145   else
00146   {
00147     rtn = new GenericRightPreconditioner<double>(ilu);
00148   }
00149 }
00150 
00151 void EpetraMatrix::getICCPreconditioner(int fillLevels,
00152   int overlapFill,
00153   double relaxationValue,
00154   double relativeThreshold,
00155   double absoluteThreshold,
00156   Preconditioner<double>& rtn) const
00157 {
00158   RCP<LinearOperatorBase<double> > a = rcp(new IfpackICCOperator(this,
00159       fillLevels,
00160       overlapFill,
00161       relaxationValue,
00162       relativeThreshold,
00163       absoluteThreshold));
00164 
00165   LinearOperator<double> icc = a;
00166   LinearOperator<double> iccTrans=icc.transpose();
00167   //Transpose call and then put in for second argument for below
00168 
00169     rtn = new GenericTwoSidedPreconditioner<double>(icc,iccTrans);
00170 }
00171 
00172 
00173 void EpetraMatrix::print(std::ostream& os) const 
00174 {
00175   crsMatrix()->Print(os);
00176 }
00177 
00178 
00179 string EpetraMatrix::description() const 
00180 {
00181   std::string rtn = "EpetraMatrix[nRow=" 
00182     + Teuchos::toString(crsMatrix()->NumGlobalRows())
00183     + ", nCol=" + Teuchos::toString(crsMatrix()->NumGlobalCols())
00184     + "]";
00185   return rtn;
00186 }
00187 
00188 
00189 Epetra_CrsMatrix* EpetraMatrix::crsMatrix()
00190 {
00191   return matrix_.get();
00192 }
00193 
00194 
00195 const Epetra_CrsMatrix* EpetraMatrix::crsMatrix() const 
00196 {
00197   return matrix_.get();
00198 }
00199 
00200 
00201 Epetra_CrsMatrix& EpetraMatrix::getConcrete(const LinearOperator<double>& A)
00202 {
00203   return *Teuchos::dyn_cast<EpetraMatrix>(*A.ptr()).crsMatrix();
00204 }
00205 
00206 
00207 RCP<const Epetra_CrsMatrix>
00208 EpetraMatrix::getConcretePtr(const LinearOperator<double>& A)
00209 {
00210   return Teuchos::rcp_dynamic_cast<EpetraMatrix>(A.ptr())->matrix_;
00211 }
00212 
00213 
00214 void EpetraMatrix::getRow(const int& row, 
00215   Teuchos::Array<int>& indices, 
00216   Teuchos::Array<double>& values) const
00217 {
00218   const Epetra_CrsMatrix* crs = crsMatrix();
00219 
00220   int numEntries;
00221   int* epIndices;
00222   double* epValues;
00223 
00224   int info = crs->ExtractGlobalRowView(row, numEntries, epValues, epIndices);
00225   TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::runtime_error,
00226     "call to ExtractGlobalRowView not successful");
00227 
00228   indices.resize(numEntries);
00229   values.resize(numEntries);
00230   for (int i = 0; i < numEntries; i++)
00231   {
00232     indices[i] = *epIndices;
00233     values[i] = *epValues;
00234     epIndices++;
00235     epValues++;
00236   }
00237 }
00238 
00239 
00240 

Site Contact