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
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #include "PlayaEpetraMatrix.hpp"
00043 #include "PlayaEpetraVector.hpp"
00044 #include "PlayaVectorSpaceDecl.hpp"
00045 #include "PlayaVectorDecl.hpp"
00046 #include "PlayaLinearOperatorDecl.hpp"
00047 #include "PlayaIfpackILUOperator.hpp"
00048 #include "PlayaIfpackICCOperator.hpp"
00049 #include "PlayaPreconditioner.hpp"
00050 #include "PlayaGenericLeftPreconditioner.hpp"
00051 #include "PlayaGenericRightPreconditioner.hpp"
00052
00053 #include "PlayaGenericTwoSidedPreconditioner.hpp"
00054
00055
00056
00057 #ifndef HAVE_TEUCHOS_EXPLICIT_INSTANTIATION
00058 #include "PlayaVectorImpl.hpp"
00059 #include "PlayaLinearOperatorImpl.hpp"
00060 #include "PlayaSimpleTransposedOpImpl.hpp"
00061 #endif
00062
00063 using namespace Playa;
00064 using namespace Teuchos;
00065
00066
00067 EpetraMatrix::EpetraMatrix(const Epetra_CrsGraph& graph,
00068 const VectorSpace<double>& domain,
00069 const VectorSpace<double>& range)
00070 : LinearOpWithSpaces<double>(domain, range),
00071 matrix_(rcp(new Epetra_CrsMatrix(Copy, graph)))
00072 {}
00073
00074 EpetraMatrix::EpetraMatrix(const RCP<Epetra_CrsMatrix>& mat,
00075 const VectorSpace<double>& domain,
00076 const VectorSpace<double>& range)
00077 : LinearOpWithSpaces<double>(domain, range),
00078 matrix_(mat)
00079 {}
00080
00081
00082 void EpetraMatrix::apply(
00083 Teuchos::ETransp applyType,
00084 const Vector<double>& in,
00085 Vector<double> out) const
00086 {
00087 const Epetra_Vector& epIn = EpetraVector::getConcrete(in);
00088 Epetra_Vector& epOut = EpetraVector::getConcrete(out);
00089 using Teuchos::rcp_dynamic_cast;
00090
00091 bool trans = applyType == TRANS;
00092
00093 int ierr = matrix_->Multiply(trans, epIn, epOut);
00094 TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00095 "EpetraMatrix::generalApply() detected ierr="
00096 << ierr << " in matrix_->Multiply()");
00097 }
00098
00099
00100
00101 void EpetraMatrix::addToRow(int globalRowIndex,
00102 int nElemsToInsert,
00103 const int* globalColumnIndices,
00104 const double* elementValues)
00105 {
00106 int ierr = crsMatrix()->SumIntoGlobalValues(globalRowIndex,
00107 nElemsToInsert,
00108 (double*) elementValues,
00109 (int*) globalColumnIndices);
00110
00111 TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00112 "failed to add to row " << globalRowIndex
00113 << " in EpetraMatrix::addToRow() with nnz="
00114 << nElemsToInsert
00115 << ". Error code was " << ierr);
00116 }
00117
00118
00119
00120 void EpetraMatrix::addToElementBatch(int numRows,
00121 int rowBlockSize,
00122 const int* globalRowIndices,
00123 int numColumnsPerRow,
00124 const int* globalColumnIndices,
00125 const double* values,
00126 const int* skipRow)
00127 {
00128 Epetra_CrsMatrix* crs = crsMatrix();
00129
00130 int numRowBlocks = numRows/rowBlockSize;
00131 int row = 0;
00132
00133 for (int rb=0; rb<numRowBlocks; rb++)
00134 {
00135 const int* cols = globalColumnIndices + rb*numColumnsPerRow;
00136 for (int r=0; r<rowBlockSize; r++, row++)
00137 {
00138 if (skipRow[row]) continue;
00139 const double* rowVals = values + row*numColumnsPerRow;
00140 int ierr=crs->SumIntoGlobalValues(globalRowIndices[row],
00141 numColumnsPerRow,
00142 (double*) rowVals,
00143 (int*) cols);
00144 TEUCHOS_TEST_FOR_EXCEPTION(ierr < 0, std::runtime_error,
00145 "failed to add to row " << globalRowIndices[row]
00146 << " in EpetraMatrix::addToRow() with nnz="
00147 << numColumnsPerRow
00148 << ". Error code was " << ierr);
00149 }
00150 }
00151 }
00152
00153
00154 void EpetraMatrix::zero()
00155 {
00156 crsMatrix()->PutScalar(0.0);
00157 }
00158
00159
00160
00161 void EpetraMatrix::getILUKPreconditioner(int fillLevels,
00162 int overlapFill,
00163 double relaxationValue,
00164 double relativeThreshold,
00165 double absoluteThreshold,
00166 LeftOrRight leftOrRight,
00167 Preconditioner<double>& rtn) const
00168 {
00169 RCP<LinearOperatorBase<double> > a = rcp(new IfpackILUOperator(this,
00170 fillLevels,
00171 overlapFill,
00172 relaxationValue,
00173 relativeThreshold,
00174 absoluteThreshold));
00175
00176 LinearOperator<double> ilu = a;
00177
00178
00179 if (leftOrRight == Left)
00180 {
00181 rtn = new GenericLeftPreconditioner<double>(ilu);
00182 }
00183 else
00184 {
00185 rtn = new GenericRightPreconditioner<double>(ilu);
00186 }
00187 }
00188
00189 void EpetraMatrix::getICCPreconditioner(int fillLevels,
00190 int overlapFill,
00191 double dropTolerance,
00192 double relaxationValue,
00193 double relativeThreshold,
00194 double absoluteThreshold,
00195 Preconditioner<double>& rtn) const
00196 {
00197 RCP<LinearOperatorBase<double> > ROp
00198 = rcp(new IfpackICCOperator(this,
00199 fillLevels,
00200 overlapFill,
00201 dropTolerance,
00202 relaxationValue,
00203 relativeThreshold,
00204 absoluteThreshold));
00205
00206 LinearOperator<double> R = ROp;
00207 LinearOperator<double> Rt= R.transpose();
00208 rtn = new GenericTwoSidedPreconditioner<double>(Rt, R);
00209 }
00210
00211
00212 void EpetraMatrix::print(std::ostream& os) const
00213 {
00214 crsMatrix()->Print(os);
00215 }
00216
00217
00218 string EpetraMatrix::description() const
00219 {
00220 std::string rtn = "EpetraMatrix[nRow="
00221 + Teuchos::toString(crsMatrix()->NumGlobalRows())
00222 + ", nCol=" + Teuchos::toString(crsMatrix()->NumGlobalCols())
00223 + "]";
00224 return rtn;
00225 }
00226
00227
00228 Epetra_CrsMatrix* EpetraMatrix::crsMatrix()
00229 {
00230 return matrix_.get();
00231 }
00232
00233
00234 const Epetra_CrsMatrix* EpetraMatrix::crsMatrix() const
00235 {
00236 return matrix_.get();
00237 }
00238
00239
00240 Epetra_CrsMatrix& EpetraMatrix::getConcrete(const LinearOperator<double>& A)
00241 {
00242 return *Teuchos::dyn_cast<EpetraMatrix>(*A.ptr()).crsMatrix();
00243 }
00244
00245
00246 RCP<const Epetra_CrsMatrix>
00247 EpetraMatrix::getConcretePtr(const LinearOperator<double>& A)
00248 {
00249 return Teuchos::rcp_dynamic_cast<EpetraMatrix>(A.ptr())->matrix_;
00250 }
00251
00252
00253 void EpetraMatrix::getRow(const int& row,
00254 Teuchos::Array<int>& indices,
00255 Teuchos::Array<double>& values) const
00256 {
00257 const Epetra_CrsMatrix* crs = crsMatrix();
00258
00259 int numEntries;
00260 int* epIndices;
00261 double* epValues;
00262
00263 int info = crs->ExtractGlobalRowView(row, numEntries, epValues, epIndices);
00264 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::runtime_error,
00265 "call to ExtractGlobalRowView not successful");
00266
00267 indices.resize(numEntries);
00268 values.resize(numEntries);
00269 for (int i = 0; i < numEntries; i++)
00270 {
00271 indices[i] = *epIndices;
00272 values[i] = *epValues;
00273 epIndices++;
00274 epValues++;
00275 }
00276 }
00277
00278
00279