|
Ifpack2 Templated Preconditioning Package
Version 1.0
|
00001 /*@HEADER 00002 // *********************************************************************** 00003 // 00004 // Ifpack2: Tempated Object-Oriented Algebraic Preconditioner Package 00005 // Copyright (2009) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // *********************************************************************** 00040 //@HEADER 00041 */ 00042 00043 #ifndef IFPACK2_BLOCKRELAXATION_DEF_HPP 00044 #define IFPACK2_BLOCKRELAXATION_DEF_HPP 00045 00046 #include "Ifpack2_BlockRelaxation_decl.hpp" 00047 #include "Ifpack2_LinearPartitioner.hpp" 00048 #include "Ifpack2_Details_UserPartitioner_decl.hpp" 00049 #include "Ifpack2_Details_UserPartitioner_def.hpp" 00050 #include <Ifpack2_Condest.hpp> 00051 #include <Ifpack2_Parameters.hpp> 00052 00053 namespace Ifpack2 { 00054 00055 template<class MatrixType,class ContainerType> 00056 void BlockRelaxation<MatrixType,ContainerType>:: 00057 setMatrix (const Teuchos::RCP<const row_matrix_type>& A) 00058 { 00059 if (A.getRawPtr () != A_.getRawPtr ()) { // it's a different matrix 00060 IsInitialized_ = false; 00061 IsComputed_ = false; 00062 Condest_ = -STM::one (); 00063 Partitioner_ = Teuchos::null; 00064 Importer_ = Teuchos::null; 00065 W_ = Teuchos::null; 00066 00067 if (! A.is_null ()) { 00068 IsParallel_ = (A->getRowMap ()->getComm ()->getSize () > 1); 00069 } 00070 00071 std::vector<Teuchos::RCP<ContainerType> > emptyVec; 00072 std::swap (Containers_, emptyVec); 00073 NumLocalBlocks_ = 0; 00074 00075 A_ = A; 00076 } 00077 } 00078 00079 template<class MatrixType,class ContainerType> 00080 BlockRelaxation<MatrixType,ContainerType>:: 00081 BlockRelaxation (const Teuchos::RCP<const row_matrix_type>& A) 00082 : A_ (A), 00083 Time_ (Teuchos::rcp (new Teuchos::Time ("Ifpack2::BlockRelaxation"))), 00084 OverlapLevel_ (0), 00085 PartitionerType_ ("linear"), 00086 NumSweeps_ (1), 00087 NumLocalBlocks_(0), 00088 PrecType_ (Ifpack2::Details::JACOBI), 00089 DampingFactor_ (STS::one ()), 00090 IsParallel_ (false), 00091 ZeroStartingSolution_ (true), 00092 DoBackwardGS_ (false), 00093 Condest_ (-STM::one ()), 00094 IsInitialized_ (false), 00095 IsComputed_ (false), 00096 NumInitialize_ (0), 00097 NumCompute_ (0), 00098 NumApply_ (0), 00099 InitializeTime_ (0.0), 00100 ComputeTime_ (0.0), 00101 ApplyTime_ (0.0), 00102 ComputeFlops_ (0.0), 00103 ApplyFlops_ (0.0), 00104 NumMyRows_ (0), 00105 NumGlobalRows_ (0), 00106 NumGlobalNonzeros_ (0) 00107 { 00108 TEUCHOS_TEST_FOR_EXCEPTION( 00109 A_.is_null (), std::invalid_argument, 00110 Teuchos::typeName(*this) << "::BlockRelaxation(): input matrix is null."); 00111 } 00112 00113 //========================================================================== 00114 template<class MatrixType,class ContainerType> 00115 BlockRelaxation<MatrixType,ContainerType>::~BlockRelaxation() {} 00116 00117 //========================================================================== 00118 template<class MatrixType,class ContainerType> 00119 void 00120 BlockRelaxation<MatrixType,ContainerType>:: 00121 setParameters (const Teuchos::ParameterList& List) 00122 { 00123 Teuchos::ParameterList validparams; 00124 Ifpack2::getValidParameters (validparams); 00125 List.validateParameters (validparams); 00126 00127 std::string PT; 00128 if (PrecType_ == Ifpack2::Details::JACOBI) { 00129 PT = "Jacobi"; 00130 } else if (PrecType_ == Ifpack2::Details::GS) { 00131 PT = "Gauss-Seidel"; 00132 } else if (PrecType_ == Ifpack2::Details::SGS) { 00133 PT = "Symmetric Gauss-Seidel"; 00134 } 00135 00136 Ifpack2::getParameter (List, "relaxation: type", PT); 00137 00138 if (PT == "Jacobi") { 00139 PrecType_ = Ifpack2::Details::JACOBI; 00140 } 00141 else if (PT == "Gauss-Seidel") { 00142 PrecType_ = Ifpack2::Details::GS; 00143 } 00144 else if (PT == "Symmetric Gauss-Seidel") { 00145 PrecType_ = Ifpack2::Details::SGS; 00146 } 00147 else { 00148 TEUCHOS_TEST_FOR_EXCEPTION( 00149 true, std::invalid_argument, "Ifpack2::BlockRelaxation::setParameters: " 00150 "Invalid parameter value \"" << PT << "\" for parameter \"relaxation: type\"."); 00151 } 00152 00153 Ifpack2::getParameter (List, "relaxation: sweeps",NumSweeps_); 00154 Ifpack2::getParameter (List, "relaxation: damping factor", DampingFactor_); 00155 Ifpack2::getParameter (List, "relaxation: zero starting solution", ZeroStartingSolution_); 00156 Ifpack2::getParameter (List, "relaxation: backward mode",DoBackwardGS_); 00157 Ifpack2::getParameter (List, "partitioner: type",PartitionerType_); 00158 Ifpack2::getParameter (List, "partitioner: local parts",NumLocalBlocks_); 00159 Ifpack2::getParameter (List, "partitioner: overlap",OverlapLevel_); 00160 00161 // check parameters 00162 if (PrecType_ != Ifpack2::Details::JACOBI) { 00163 OverlapLevel_ = 0; 00164 } 00165 if (NumLocalBlocks_ < 0) { 00166 NumLocalBlocks_ = A_->getNodeNumRows() / (-NumLocalBlocks_); 00167 } 00168 // other checks are performed in Partitioner_ 00169 00170 // NTS: Sanity check to be removed at a later date when Backward mode is enabled 00171 TEUCHOS_TEST_FOR_EXCEPTION( 00172 DoBackwardGS_, std::runtime_error, 00173 "Ifpack2::BlockRelaxation:setParameters: \"relaxation: backward mode\" == " 00174 "true is not supported yet."); 00175 00176 // copy the list as each subblock's constructor will 00177 // require it later 00178 List_ = List; 00179 } 00180 00181 //========================================================================== 00182 template<class MatrixType,class ContainerType> 00183 Teuchos::RCP<const Teuchos::Comm<int> > 00184 BlockRelaxation<MatrixType,ContainerType>::getComm() const{ 00185 return A_->getComm(); 00186 } 00187 00188 //========================================================================== 00189 template<class MatrixType,class ContainerType> 00190 Teuchos::RCP<const Tpetra::RowMatrix<typename MatrixType::scalar_type, 00191 typename MatrixType::local_ordinal_type, 00192 typename MatrixType::global_ordinal_type, 00193 typename MatrixType::node_type> > 00194 BlockRelaxation<MatrixType,ContainerType>::getMatrix() const { 00195 return(A_); 00196 } 00197 00198 //========================================================================== 00199 template<class MatrixType,class ContainerType> 00200 Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type, 00201 typename MatrixType::global_ordinal_type, 00202 typename MatrixType::node_type> > 00203 BlockRelaxation<MatrixType,ContainerType>::getDomainMap() const { 00204 return A_->getDomainMap(); 00205 } 00206 00207 //========================================================================== 00208 template<class MatrixType,class ContainerType> 00209 Teuchos::RCP<const Tpetra::Map<typename MatrixType::local_ordinal_type, 00210 typename MatrixType::global_ordinal_type, 00211 typename MatrixType::node_type> > 00212 BlockRelaxation<MatrixType,ContainerType>::getRangeMap() const { 00213 return A_->getRangeMap(); 00214 } 00215 00216 //========================================================================== 00217 template<class MatrixType,class ContainerType> 00218 bool BlockRelaxation<MatrixType,ContainerType>::hasTransposeApply() const { 00219 return true; 00220 } 00221 00222 //========================================================================== 00223 template<class MatrixType,class ContainerType> 00224 int BlockRelaxation<MatrixType,ContainerType>::getNumInitialize() const { 00225 return(NumInitialize_); 00226 } 00227 00228 //========================================================================== 00229 template<class MatrixType,class ContainerType> 00230 int BlockRelaxation<MatrixType,ContainerType>::getNumCompute() const { 00231 return(NumCompute_); 00232 } 00233 00234 //========================================================================== 00235 template<class MatrixType,class ContainerType> 00236 int BlockRelaxation<MatrixType,ContainerType>::getNumApply() const { 00237 return(NumApply_); 00238 } 00239 00240 //========================================================================== 00241 template<class MatrixType,class ContainerType> 00242 double BlockRelaxation<MatrixType,ContainerType>::getInitializeTime() const { 00243 return(InitializeTime_); 00244 } 00245 00246 //========================================================================== 00247 template<class MatrixType,class ContainerType> 00248 double BlockRelaxation<MatrixType,ContainerType>::getComputeTime() const { 00249 return(ComputeTime_); 00250 } 00251 00252 //========================================================================== 00253 template<class MatrixType,class ContainerType> 00254 double BlockRelaxation<MatrixType,ContainerType>::getApplyTime() const { 00255 return(ApplyTime_); 00256 } 00257 00258 //========================================================================== 00259 template<class MatrixType,class ContainerType> 00260 double BlockRelaxation<MatrixType,ContainerType>::getComputeFlops() const { 00261 return(ComputeFlops_); 00262 } 00263 00264 //========================================================================== 00265 template<class MatrixType,class ContainerType> 00266 double BlockRelaxation<MatrixType,ContainerType>::getApplyFlops() const { 00267 return ApplyFlops_; 00268 } 00269 00270 //========================================================================== 00271 template<class MatrixType,class ContainerType> 00272 typename BlockRelaxation<MatrixType, ContainerType>::magnitude_type 00273 BlockRelaxation<MatrixType,ContainerType>::getCondEst () const 00274 { 00275 return Condest_; 00276 } 00277 00278 //========================================================================== 00279 template<class MatrixType,class ContainerType> 00280 typename BlockRelaxation<MatrixType, ContainerType>::magnitude_type 00281 BlockRelaxation<MatrixType,ContainerType>:: 00282 computeCondEst (CondestType CT, 00283 typename MatrixType::local_ordinal_type MaxIters, 00284 magnitude_type Tol, 00285 const Teuchos::Ptr<const row_matrix_type>& matrix) 00286 { 00287 if (! isComputed ()) {// cannot compute right now 00288 return -STM::one (); 00289 } 00290 00291 // always compute it. Call Condest() with no parameters to get 00292 // the previous estimate. 00293 Condest_ = Ifpack2::Condest (*this, CT, MaxIters, Tol, matrix); 00294 00295 return Condest_; 00296 } 00297 00298 //========================================================================== 00299 template<class MatrixType,class ContainerType> 00300 void 00301 BlockRelaxation<MatrixType,ContainerType>:: 00302 apply (const Tpetra::MultiVector<typename MatrixType::scalar_type, 00303 typename MatrixType::local_ordinal_type, 00304 typename MatrixType::global_ordinal_type, 00305 typename MatrixType::node_type>& X, 00306 Tpetra::MultiVector<typename MatrixType::scalar_type, 00307 typename MatrixType::local_ordinal_type, 00308 typename MatrixType::global_ordinal_type, 00309 typename MatrixType::node_type>& Y, 00310 Teuchos::ETransp mode, 00311 scalar_type alpha, 00312 scalar_type beta) const 00313 { 00314 TEUCHOS_TEST_FOR_EXCEPTION( 00315 ! isComputed (), std::runtime_error, "Ifpack2::BlockRelaxation::apply: " 00316 "isComputed() must be true prior to calling apply."); 00317 TEUCHOS_TEST_FOR_EXCEPTION( 00318 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument, 00319 "Ifpack2::BlockRelaxation::apply: X.getNumVectors() = " 00320 << X.getNumVectors () << " != Y.getNumVectors() = " 00321 << Y.getNumVectors () << "."); 00322 TEUCHOS_TEST_FOR_EXCEPTION( 00323 mode != Teuchos::NO_TRANS, std::logic_error, "Ifpack2::BlockRelaxation::" 00324 "apply: This method currently only implements the case mode == Teuchos::" 00325 "NO_TRANS. Transposed modes are not currently supported."); 00326 TEUCHOS_TEST_FOR_EXCEPTION( 00327 alpha != Teuchos::ScalarTraits<scalar_type>::one (), std::logic_error, 00328 "Ifpack2::BlockRelaxation::apply: This method currently only implements " 00329 "the case alpha == 1. You specified alpha = " << alpha << "."); 00330 TEUCHOS_TEST_FOR_EXCEPTION( 00331 beta != Teuchos::ScalarTraits<scalar_type>::zero (), std::logic_error, 00332 "Ifpack2::BlockRelaxation::apply: This method currently only implements " 00333 "the case beta == 0. You specified beta = " << beta << "."); 00334 00335 Time_->start(true); 00336 00337 // If X and Y are pointing to the same memory location, 00338 // we need to create an auxiliary vector, Xcopy 00339 Teuchos::RCP<const MV> X_copy; 00340 if (X.getLocalMV ().getValues () == Y.getLocalMV ().getValues ()) { 00341 X_copy = Teuchos::rcp (new MV (X, Teuchos::Copy)); 00342 } else { 00343 X_copy = Teuchos::rcpFromRef (X); 00344 } 00345 00346 if (ZeroStartingSolution_) { 00347 Y.putScalar (STS::zero ()); 00348 } 00349 00350 // Flops are updated in each of the following. 00351 switch (PrecType_) { 00352 case Ifpack2::Details::JACOBI: 00353 ApplyInverseJacobi(*X_copy,Y); 00354 break; 00355 case Ifpack2::Details::GS: 00356 ApplyInverseGS(*X_copy,Y); 00357 break; 00358 case Ifpack2::Details::SGS: 00359 ApplyInverseSGS(*X_copy,Y); 00360 break; 00361 default: 00362 throw std::runtime_error("Ifpack2::BlockRelaxation::apply internal logic error."); 00363 } 00364 00365 ++NumApply_; 00366 Time_->stop(); 00367 ApplyTime_ += Time_->totalElapsedTime(); 00368 } 00369 00370 //========================================================================== 00371 template<class MatrixType,class ContainerType> 00372 void BlockRelaxation<MatrixType,ContainerType>::applyMat( 00373 const Tpetra::MultiVector<typename MatrixType::scalar_type, 00374 typename MatrixType::local_ordinal_type, 00375 typename MatrixType::global_ordinal_type, 00376 typename MatrixType::node_type>& X, 00377 Tpetra::MultiVector<typename MatrixType::scalar_type, 00378 typename MatrixType::local_ordinal_type, 00379 typename MatrixType::global_ordinal_type, 00380 typename MatrixType::node_type>& Y, 00381 Teuchos::ETransp mode) const 00382 { 00383 TEUCHOS_TEST_FOR_EXCEPTION(isComputed() == false, std::runtime_error, 00384 "Ifpack2::BlockRelaxation::applyMat() ERROR: isComputed() must be true prior to calling applyMat()."); 00385 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error, 00386 "Ifpack2::BlockRelaxation::applyMat() ERROR: X.getNumVectors() != Y.getNumVectors()."); 00387 A_->apply (X, Y, mode); 00388 } 00389 00390 //========================================================================== 00391 template<class MatrixType,class ContainerType> 00392 void BlockRelaxation<MatrixType,ContainerType>::initialize() { 00393 using Teuchos::rcp; 00394 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> 00395 row_graph_type; 00396 00397 IsInitialized_ = false; 00398 Time_->start (true); 00399 00400 NumMyRows_ = A_->getNodeNumRows (); 00401 NumGlobalRows_ = A_->getGlobalNumRows (); 00402 NumGlobalNonzeros_ = A_->getGlobalNumEntries (); 00403 00404 // NTS: Will need to add support for Zoltan2 partitions later Also, 00405 // will need a better way of handling the Graph typing issue. 00406 // Especially with ordinal types w.r.t the container. 00407 00408 if (PartitionerType_ == "linear") { 00409 Partitioner_ = 00410 rcp (new Ifpack2::LinearPartitioner<row_graph_type> (A_->getGraph ())); 00411 } else if (PartitionerType_ == "user") { 00412 Partitioner_ = 00413 rcp (new Ifpack2::Details::UserPartitioner<row_graph_type> (A_->getGraph () ) ); 00414 } else { 00415 // We should have checked for this in setParameters(), so it's a 00416 // logic_error, not an invalid_argument or runtime_error. 00417 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, 00418 "Ifpack2::BlockRelaxation::initialize, invalid partitioner type."); 00419 } 00420 00421 // need to partition the graph of A 00422 Partitioner_->setParameters (List_); 00423 Partitioner_->compute (); 00424 00425 // get actual number of partitions 00426 NumLocalBlocks_ = Partitioner_->numLocalParts (); 00427 00428 // Note: Unlike Ifpack, we'll punt on computing W_ until compute(), which is where 00429 // we assume that the type of relaxation has been chosen. 00430 00431 if (A_->getComm()->getSize() != 1) { 00432 IsParallel_ = true; 00433 } else { 00434 IsParallel_ = false; 00435 } 00436 00437 ++NumInitialize_; 00438 Time_->stop (); 00439 InitializeTime_ += Time_->totalElapsedTime (); 00440 IsInitialized_ = true; 00441 } 00442 00443 //========================================================================== 00444 template<class MatrixType,class ContainerType> 00445 void BlockRelaxation<MatrixType,ContainerType>::compute() 00446 { 00447 using Teuchos::rcp; 00448 typedef Tpetra::Vector<scalar_type, 00449 local_ordinal_type, global_ordinal_type, node_type> vector_type; 00450 typedef Tpetra::Import<local_ordinal_type, 00451 global_ordinal_type, node_type> import_type; 00452 00453 // We should have checked for this in setParameters(), so it's a 00454 // logic_error, not an invalid_argument or runtime_error. 00455 TEUCHOS_TEST_FOR_EXCEPTION(NumSweeps_ < 0, std::logic_error, 00456 "Ifpack2::BlockRelaxation::compute, NumSweeps_ must be >= 0"); 00457 00458 if (! isInitialized ()) { 00459 initialize (); 00460 } 00461 Time_->start (true); 00462 00463 // reset values 00464 IsComputed_ = false; 00465 Condest_ = -STM::one (); 00466 00467 // Extract the submatrices 00468 ExtractSubmatrices (); 00469 00470 // Compute the weight vector if we're doing overlapped Jacobi (and 00471 // only if we're doing overlapped Jacobi). 00472 if (PrecType_ == Ifpack2::Details::JACOBI && OverlapLevel_ > 0) { 00473 // weight of each vertex 00474 W_ = rcp (new vector_type (A_->getRowMap ())); 00475 W_->putScalar (STS::zero ()); 00476 Teuchos::ArrayRCP<scalar_type > w_ptr = W_->getDataNonConst(0); 00477 00478 for (local_ordinal_type i = 0 ; i < NumLocalBlocks_ ; ++i) { 00479 for (size_t j = 0 ; j < Partitioner_->numRowsInPart(i) ; ++j) { 00480 // FIXME (mfh 12 Sep 2014) Should this really be int? 00481 // Perhaps it should be local_ordinal_type instead. 00482 int LID = (*Partitioner_)(i,j); 00483 w_ptr[LID]+= STS::one(); 00484 } 00485 } 00486 W_->reciprocal (*W_); 00487 } 00488 00489 // We need to import data from external processors. Here I create a 00490 // Tpetra::Import object if needed (stealing from A_ if possible) 00491 // Marzio's comment: 00492 // Note that I am doing some strange stuff to set the components of Y 00493 // from Y2 (to save some time). 00494 // 00495 if (IsParallel_ && (PrecType_ == Ifpack2::Details::GS || 00496 PrecType_ == Ifpack2::Details::SGS)) { 00497 Importer_ = A_->getGraph ()->getImporter (); 00498 if (Importer_.is_null ()) { 00499 Importer_ = rcp (new import_type (A_->getDomainMap (), A_->getColMap ())); 00500 } 00501 } 00502 00503 ++NumCompute_; 00504 Time_->stop (); 00505 ComputeTime_ += Time_->totalElapsedTime(); 00506 IsComputed_ = true; 00507 } 00508 00509 //============================================================================== 00510 template<class MatrixType,class ContainerType> 00511 void BlockRelaxation<MatrixType,ContainerType>::ExtractSubmatrices() 00512 { 00513 typedef Tpetra::Vector<scalar_type, local_ordinal_type, 00514 global_ordinal_type, node_type> vec_type; 00515 TEUCHOS_TEST_FOR_EXCEPTION( 00516 Partitioner_.is_null (), std::runtime_error, 00517 "Ifpack2::BlockRelaxation::ExtractSubmatrices: Partitioner object is null."); 00518 00519 NumLocalBlocks_ = Partitioner_->numLocalParts (); 00520 Containers_.resize (NumLocalBlocks_); 00521 vec_type D (A_->getRowMap ()); 00522 A_->getLocalDiagCopy (D); 00523 DiagRCP = D.getData (); 00524 00525 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) { 00526 const size_t numRows = Partitioner_->numRowsInPart (i); 00527 00528 // Extract a list of the indices of each partitioner row. 00529 Teuchos::Array<local_ordinal_type> localRows (numRows); 00530 for (size_t j = 0; j < numRows; ++j) { 00531 localRows[j] = (*Partitioner_) (i,j); 00532 } 00533 if(numRows>1) { // only do for non-singletons 00534 Containers_[i] = Teuchos::rcp (new ContainerType (A_, localRows ())); 00535 Containers_[i]->setParameters (List_); 00536 Containers_[i]->initialize (); 00537 Containers_[i]->compute (); 00538 } 00539 } 00540 } 00541 00542 00543 00544 //========================================================================== 00545 template<class MatrixType,class ContainerType> 00546 void BlockRelaxation<MatrixType,ContainerType>::ApplyInverseJacobi (const MV& X, MV& Y) const 00547 { 00548 const size_t NumVectors = X.getNumVectors (); 00549 MV AY (Y.getMap (), NumVectors); 00550 00551 // Initial matvec not needed 00552 int starting_iteration = 0; 00553 if (ZeroStartingSolution_) { 00554 DoJacobi (X, Y); 00555 starting_iteration = 1; 00556 } 00557 00558 const scalar_type ONE = STS::one (); 00559 for (int j = starting_iteration; j < NumSweeps_; ++j) { 00560 applyMat (Y, AY); 00561 AY.update (ONE, X, -ONE); 00562 DoJacobi (AY, Y); 00563 00564 // Flops for matrix apply & update 00565 ApplyFlops_ += NumVectors * (2 * NumGlobalNonzeros_ + 2 * NumGlobalRows_); 00566 } 00567 00568 } 00569 00570 00571 //========================================================================== 00572 template<class MatrixType,class ContainerType> 00573 void BlockRelaxation<MatrixType,ContainerType>::DoJacobi(const MV& X, MV& Y) const 00574 { 00575 const size_t NumVectors = X.getNumVectors (); 00576 const scalar_type one = STS::one (); 00577 // Note: Flop counts copied naively from Ifpack. 00578 00579 if (OverlapLevel_ == 0) { 00580 // Non-overlapping Jacobi 00581 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) { 00582 // may happen that a partition is empty 00583 if( Partitioner_->numRowsInPart (i) != 1 ) { 00584 if(Containers_[i]->getNumRows () == 0 ) continue; 00585 Containers_[i]->apply (X, Y, Teuchos::NO_TRANS, DampingFactor_, one); 00586 ApplyFlops_ += NumVectors * 2 * NumGlobalRows_; 00587 } 00588 else { // singleton, can't access Containers_[i] as it was never filled and may be null. 00589 local_ordinal_type LRID = (*Partitioner_)(i,0); // by definition, a singleton 1 row in block. 00590 Teuchos::ArrayView< const scalar_type > Diag = DiagRCP(); 00591 scalar_type d = Diag[LRID]; 00592 for(unsigned int nv = 0;nv < NumVectors ; ++nv ) { 00593 Teuchos::ArrayRCP< const scalar_type > xRCP = X.getData(nv); 00594 scalar_type x = xRCP[LRID]; 00595 Teuchos::ArrayRCP< scalar_type > yRCP = Y.getDataNonConst(nv); 00596 00597 scalar_type newy= x/d; 00598 yRCP[LRID]= newy; 00599 } 00600 } 00601 } 00602 } 00603 else { 00604 // Overlapping Jacobi 00605 for (local_ordinal_type i = 0 ; i < NumLocalBlocks_ ; i++) { 00606 // may happen that a partition is empty 00607 if(Containers_[i]->getNumRows() == 0) continue; 00608 if ( Partitioner_->numRowsInPart (i) != 1 ) { 00609 try { 00610 Containers_[i]->weightedApply(X,Y,*W_,Teuchos::NO_TRANS,DampingFactor_,one); 00611 } catch (std::exception& e) { 00612 std::cerr << "BlockRelaxation::DoJacobi: Containers_[" << i 00613 << "]->weightedApply() threw an exception: " << e.what () 00614 << std::endl; 00615 throw; 00616 } 00617 } // end Partitioner_->numRowsInPart (i) != 1 00618 00619 // NOTE: do not count (for simplicity) the flops due to overlapping rows 00620 ApplyFlops_ += NumVectors * 4 * NumGlobalRows_; 00621 } 00622 } 00623 } 00624 00625 //========================================================================== 00626 template<class MatrixType,class ContainerType> 00627 void BlockRelaxation<MatrixType,ContainerType>:: 00628 ApplyInverseGS (const MV& X, MV& Y) const 00629 { 00630 MV Xcopy (X, Teuchos::Copy); 00631 for (int j = 0; j < NumSweeps_; ++j) { 00632 DoGaussSeidel (Xcopy, Y); 00633 if (j != NumSweeps_ - 1) { 00634 Tpetra::deep_copy (Xcopy, X); 00635 } 00636 } 00637 } 00638 00639 00640 //============================================================================== 00641 template<class MatrixType,class ContainerType> 00642 void BlockRelaxation<MatrixType,ContainerType>:: 00643 DoGaussSeidel (MV& X, MV& Y) const 00644 { 00645 using Teuchos::Array; 00646 using Teuchos::ArrayRCP; 00647 using Teuchos::ArrayView; 00648 using Teuchos::RCP; 00649 using Teuchos::rcp; 00650 using Teuchos::rcpFromRef; 00651 00652 // Note: Flop counts copied naively from Ifpack. 00653 00654 const scalar_type one = STS::one (); 00655 const size_t Length = A_->getNodeMaxNumRowEntries(); 00656 const size_t NumVectors = X.getNumVectors(); 00657 Array<scalar_type> Values; 00658 Array<local_ordinal_type> Indices; 00659 Values.resize (Length); 00660 Indices.resize (Length); 00661 00662 // an additonal vector is needed by parallel computations 00663 // (note that applications through Ifpack2_AdditiveSchwarz 00664 // are always seen are serial) 00665 RCP<MV> Y2; 00666 if (IsParallel_) { 00667 Y2 = rcp (new MV (Importer_->getTargetMap (), NumVectors)); 00668 } else { 00669 Y2 = rcpFromRef (Y); 00670 } 00671 00672 // I think I decided I need two extra vectors: 00673 // One to store the sum of the corrections (initialized to zero) 00674 // One to store the temporary residual (doesn't matter if it is zeroed or not) 00675 // My apologies for making the names clear and meaningful. (X=RHS, Y=guess?! Nice.) 00676 MV Residual (X.getMap (), NumVectors, false); 00677 00678 ArrayRCP<ArrayRCP<scalar_type> > x_ptr = X.get2dViewNonConst(); 00679 ArrayRCP<ArrayRCP<scalar_type> > y_ptr = Y.get2dViewNonConst(); 00680 ArrayRCP<ArrayRCP<scalar_type> > y2_ptr = Y2->get2dViewNonConst(); 00681 ArrayRCP<ArrayRCP<scalar_type> > residual_ptr = Residual.get2dViewNonConst(); 00682 00683 // data exchange is here, once per sweep 00684 if (IsParallel_) Y2->doImport(Y,*Importer_,Tpetra::INSERT); 00685 00686 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) { 00687 if( Partitioner_->numRowsInPart (i) != 1 ) { 00688 if (Containers_[i]->getNumRows () == 0) continue; 00689 // update from previous block 00690 ArrayView<const local_ordinal_type> localRows = 00691 Containers_[i]->getLocalRows (); 00692 const size_t localNumRows = Containers_[i]->getNumRows (); 00693 for (size_t j = 0; j < localNumRows; ++j) { 00694 const local_ordinal_type LID = localRows[j]; // Containers_[i]->ID (j); 00695 size_t NumEntries; 00696 A_->getLocalRowCopy (LID, Indices (), Values (), NumEntries); 00697 00698 for (size_t m = 0; m < NumVectors; ++m) { 00699 ArrayView<const scalar_type> x_local = (x_ptr())[m](); 00700 ArrayView<scalar_type> y2_local = (y2_ptr())[m](); 00701 ArrayView<scalar_type> r_local = (residual_ptr())[m](); 00702 00703 r_local[LID] = x_local[LID]; 00704 for (size_t k = 0; k < NumEntries; ++k) { 00705 const local_ordinal_type col = Indices[k]; 00706 r_local[LID] -= Values[k] * y2_local[col]; 00707 } 00708 } 00709 } 00710 // solve with this block 00711 // 00712 // Note: I'm abusing the ordering information, knowing that X/Y 00713 // and Y2 have the same ordering for on-proc unknowns. 00714 // 00715 // Note: Add flop counts for inverse apply 00716 Containers_[i]->apply (Residual, *Y2, Teuchos::NO_TRANS, 00717 DampingFactor_,one); 00718 00719 // operations for all getrow's 00720 ApplyFlops_ += NumVectors * (2 * NumGlobalNonzeros_ + 2 * NumGlobalRows_); 00721 } 00722 else { // singleton, can't access Containers_[i] as it was never filled and may be null. 00723 // a singleton calculation is exact, all residuals should be zero. 00724 local_ordinal_type LRID = (*Partitioner_)(i,0); // by definition, a singleton 1 row in block. 00725 Teuchos::ArrayView< const scalar_type > Diag = DiagRCP(); 00726 scalar_type d = Diag[LRID]; 00727 ArrayRCP<ArrayRCP<scalar_type> > y2_ptr2 = Y2->get2dViewNonConst(); 00728 for(unsigned int nv = 0;nv < NumVectors ; ++nv ) { 00729 Teuchos::ArrayRCP< const scalar_type > xRCP = X.getData(nv); 00730 scalar_type x = xRCP[LRID]; 00731 ArrayView<scalar_type> y2_local = (y2_ptr2())[nv](); 00732 scalar_type newy= x/d; 00733 y2_local[LRID]= newy; 00734 } 00735 } // end else 00736 } // end for NumLocalBlocks_ 00737 // Attention: this is delicate... Not all combinations 00738 // of Y2 and Y will always work (tough for ML it should be ok) 00739 if (IsParallel_) { 00740 for (size_t m = 0; m < NumVectors; ++m) { 00741 ArrayView<scalar_type> y2_local = (y2_ptr())[m](); 00742 ArrayView<scalar_type> y_local = (y_ptr())[m](); 00743 for (size_t i = 0; i < NumMyRows_; ++i) { 00744 y_local[i] = y2_local[i]; 00745 } 00746 } 00747 } 00748 } 00749 00750 //========================================================================== 00751 template<class MatrixType,class ContainerType> 00752 void 00753 BlockRelaxation<MatrixType,ContainerType>:: 00754 ApplyInverseSGS (const MV& X, MV& Y) const 00755 { 00756 MV Xcopy (X, Teuchos::Copy); 00757 for (int j = 0; j < NumSweeps_; ++j) { 00758 DoSGS (Xcopy, Y); 00759 if (j != NumSweeps_ - 1) { 00760 Tpetra::deep_copy (Xcopy, X); 00761 } 00762 } 00763 } 00764 00765 //========================================================================== 00766 template<class MatrixType,class ContainerType> 00767 void 00768 BlockRelaxation<MatrixType,ContainerType>::DoSGS (MV& X, MV& Y) const 00769 { 00770 using Teuchos::Array; 00771 using Teuchos::ArrayRCP; 00772 using Teuchos::ArrayView; 00773 using Teuchos::RCP; 00774 using Teuchos::rcp; 00775 using Teuchos::rcpFromRef; 00776 00777 const scalar_type one = STS::one (); 00778 const size_t Length = A_->getNodeMaxNumRowEntries(); 00779 const size_t NumVectors = X.getNumVectors(); 00780 Array<scalar_type> Values; 00781 Array<local_ordinal_type> Indices; 00782 Values.resize(Length); 00783 Indices.resize(Length); 00784 00785 // an additonal vector is needed by parallel computations 00786 // (note that applications through Ifpack2_AdditiveSchwarz 00787 // are always seen are serial) 00788 RCP<MV> Y2; 00789 if (IsParallel_) { 00790 Y2 = rcp (new MV (Importer_->getTargetMap (), NumVectors)); 00791 } else { 00792 Y2 = rcpFromRef (Y); 00793 } 00794 00795 // I think I decided I need two extra vectors: 00796 // One to store the sum of the corrections (initialized to zero) 00797 // One to store the temporary residual (doesn't matter if it is zeroed or not) 00798 // My apologies for making the names clear and meaningful. (X=RHS, Y=guess?! Nice.) 00799 MV Residual (X.getMap (), NumVectors, false); 00800 00801 ArrayRCP<ArrayRCP<scalar_type> > x_ptr = X.get2dViewNonConst(); 00802 ArrayRCP<ArrayRCP<scalar_type> > y_ptr = Y.get2dViewNonConst(); 00803 ArrayRCP<ArrayRCP<scalar_type> > y2_ptr = Y2->get2dViewNonConst(); 00804 ArrayRCP<ArrayRCP<scalar_type> > residual_ptr = Residual.get2dViewNonConst(); 00805 00806 // data exchange is here, once per sweep 00807 if (IsParallel_) { 00808 Y2->doImport (Y, *Importer_, Tpetra::INSERT); 00809 } 00810 00811 // Forward Sweep 00812 for (local_ordinal_type i = 0; i < NumLocalBlocks_; ++i) { 00813 if( Partitioner_->numRowsInPart (i) != 1 ) { 00814 if (Containers_[i]->getNumRows () == 0) { 00815 continue; // Skip empty partitions 00816 } 00817 // update from previous block 00818 ArrayView<const local_ordinal_type> localRows = 00819 Containers_[i]->getLocalRows (); 00820 for (size_t j = 0; j < Containers_[i]->getNumRows (); ++j) { 00821 const local_ordinal_type LID = localRows[j]; // Containers_[i]->ID (j); 00822 size_t NumEntries; 00823 A_->getLocalRowCopy (LID, Indices (), Values (), NumEntries); 00824 00825 //set tmpresid = initresid - A*correction 00826 for (size_t m = 0; m < NumVectors; ++m) { 00827 ArrayView<const scalar_type> x_local = (x_ptr())[m](); 00828 ArrayView<scalar_type> y2_local = (y2_ptr())[m](); 00829 ArrayView<scalar_type> r_local = (residual_ptr())[m](); 00830 00831 r_local[LID] = x_local[LID]; 00832 for (size_t k = 0 ; k < NumEntries ; k++) { 00833 local_ordinal_type col = Indices[k]; 00834 r_local[LID] -= Values[k] * y2_local[col]; 00835 } 00836 } 00837 } 00838 // solve with this block 00839 // 00840 // Note: I'm abusing the ordering information, knowing that X/Y 00841 // and Y2 have the same ordering for on-proc unknowns. 00842 // 00843 // Note: Add flop counts for inverse apply 00844 Containers_[i]->apply (Residual, *Y2, Teuchos::NO_TRANS, 00845 DampingFactor_, one); 00846 00847 // operations for all getrow's 00848 ApplyFlops_ += NumVectors * (2 * NumGlobalNonzeros_ + 2 * NumGlobalRows_); 00849 00850 } 00851 else { // singleton, can't access Containers_[i] as it was never filled and may be null. 00852 local_ordinal_type LRID = (*Partitioner_)(i,0); // by definition, a singleton 1 row in block. 00853 Teuchos::ArrayView< const scalar_type > Diag = DiagRCP(); 00854 scalar_type d = Diag[LRID]; 00855 for(unsigned int nv = 0;nv < NumVectors ; ++nv ) { 00856 Teuchos::ArrayRCP< const scalar_type > xRCP = X.getData(nv); 00857 scalar_type x = xRCP[LRID]; 00858 Teuchos::ArrayRCP< scalar_type > yRCP = Y.getDataNonConst(nv); 00859 00860 scalar_type newy= x/d; 00861 yRCP[LRID]= newy; 00862 } 00863 } // end else 00864 } // end forward sweep over NumLocalBlocks 00865 00866 // Reverse Sweep 00867 // 00868 // mfh 12 July 2013: The unusual iteration bounds, and the use of 00869 // i-1 rather than i in the loop body, ensure correctness even if 00870 // local_ordinal_type is unsigned. "i = NumLocalBlocks_-1; i >= 0; 00871 // i--" will loop forever if local_ordinal_type is unsigned, because 00872 // unsigned integers are (trivially) always nonnegative. 00873 for (local_ordinal_type i = NumLocalBlocks_; i > 0; --i) { 00874 if( Partitioner_->numRowsInPart (i) != 1 ) { 00875 if (Containers_[i-1]->getNumRows () == 0) continue; 00876 00877 // update from previous block 00878 ArrayView<const local_ordinal_type> localRows = 00879 Containers_[i-1]->getLocalRows (); 00880 for (size_t j = 0; j < Containers_[i-1]->getNumRows (); ++j) { 00881 const local_ordinal_type LID = localRows[j]; // Containers_[i-1]->ID (j); 00882 size_t NumEntries; 00883 A_->getLocalRowCopy (LID, Indices (), Values (), NumEntries); 00884 00885 //set tmpresid = initresid - A*correction 00886 for (size_t m = 0; m < NumVectors; ++m) { 00887 ArrayView<const scalar_type> x_local = (x_ptr())[m](); 00888 ArrayView<scalar_type> y2_local = (y2_ptr())[m](); 00889 ArrayView<scalar_type> r_local = (residual_ptr())[m](); 00890 00891 r_local [LID] = x_local[LID]; 00892 for (size_t k = 0; k < NumEntries; ++k) { 00893 local_ordinal_type col = Indices[k]; 00894 r_local[LID] -= Values[k] * y2_local[col]; 00895 } 00896 } 00897 } 00898 00899 // solve with this block 00900 // 00901 // Note: I'm abusing the ordering information, knowing that X/Y 00902 // and Y2 have the same ordering for on-proc unknowns. 00903 // 00904 // Note: Add flop counts for inverse apply 00905 Containers_[i-1]->apply (Residual, *Y2, Teuchos::NO_TRANS, 00906 DampingFactor_, one); 00907 00908 // operations for all getrow's 00909 ApplyFlops_ += NumVectors * (2 * NumGlobalNonzeros_ + 2 * NumGlobalRows_); 00910 } // end Partitioner_->numRowsInPart (i) != 1 ) { 00911 // else do nothing, as by definition with a singleton, the residuals are zero. 00912 } //end reverse sweep 00913 00914 // Attention: this is delicate... Not all combinations 00915 // of Y2 and Y will always work (though for ML it should be ok) 00916 if (IsParallel_) { 00917 for (size_t m = 0; m < NumVectors; ++m) { 00918 ArrayView<scalar_type> y_local = (y_ptr())[m](); 00919 ArrayView<scalar_type> y2_local = (y2_ptr())[m](); 00920 for (size_t i = 0 ; i < NumMyRows_ ; ++i) { 00921 y_local[i] = y2_local[i]; 00922 } 00923 } 00924 } 00925 } 00926 00927 template<class MatrixType, class ContainerType> 00928 std::string BlockRelaxation<MatrixType,ContainerType>::description () const 00929 { 00930 std::ostringstream out; 00931 00932 // Output is a valid YAML dictionary in flow style. If you don't 00933 // like everything on a single line, you should call describe() 00934 // instead. 00935 out << "\"Ifpack2::BlockRelaxation\": {"; 00936 if (this->getObjectLabel () != "") { 00937 out << "Label: \"" << this->getObjectLabel () << "\", "; 00938 } 00939 out << "Initialized: " << (isInitialized () ? "true" : "false") << ", "; 00940 out << "Computed: " << (isComputed () ? "true" : "false") << ", "; 00941 00942 if (A_.is_null ()) { 00943 out << "Matrix: null"; 00944 } 00945 else { 00946 out << "Matrix: not null" 00947 << ", Global matrix dimensions: [" 00948 << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"; 00949 } 00950 00951 // It's useful to print this instance's relaxation method. If you 00952 // want more info than that, call describe() instead. 00953 out << "\"relaxation: type\": "; 00954 if (PrecType_ == Ifpack2::Details::JACOBI) { 00955 out << "Block Jacobi"; 00956 } else if (PrecType_ == Ifpack2::Details::GS) { 00957 out << "Block Gauss-Seidel"; 00958 } else if (PrecType_ == Ifpack2::Details::SGS) { 00959 out << "Block Symmetric Gauss-Seidel"; 00960 } else { 00961 out << "INVALID"; 00962 } 00963 00964 out << "}"; 00965 return out.str(); 00966 } 00967 00968 00969 template<class MatrixType,class ContainerType> 00970 void BlockRelaxation<MatrixType,ContainerType>:: 00971 describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const 00972 { 00973 using std::endl; 00974 using std::setw; 00975 using Teuchos::TypeNameTraits; 00976 using Teuchos::VERB_DEFAULT; 00977 using Teuchos::VERB_NONE; 00978 using Teuchos::VERB_LOW; 00979 using Teuchos::VERB_MEDIUM; 00980 using Teuchos::VERB_HIGH; 00981 using Teuchos::VERB_EXTREME; 00982 00983 Teuchos::EVerbosityLevel vl = verbLevel; 00984 if (vl == VERB_DEFAULT) vl = VERB_LOW; 00985 const int myImageID = A_->getComm()->getRank(); 00986 00987 // Convention requires that describe() begin with a tab. 00988 Teuchos::OSTab tab (out); 00989 00990 // none: print nothing 00991 // low: print O(1) info from node 0 00992 // medium: 00993 // high: 00994 // extreme: 00995 if (vl != VERB_NONE && myImageID == 0) { 00996 out << "Ifpack2::BlockRelaxation<" 00997 << TypeNameTraits<MatrixType>::name () << ", " 00998 << TypeNameTraits<ContainerType>::name () << " >:"; 00999 Teuchos::OSTab tab1 (out); 01000 01001 if (this->getObjectLabel () != "") { 01002 out << "label: \"" << this->getObjectLabel () << "\"" << endl; 01003 } 01004 out << "initialized: " << (isInitialized () ? "true" : "false") << endl 01005 << "computed: " << (isComputed () ? "true" : "false") << endl; 01006 01007 std::string precType; 01008 if (PrecType_ == Ifpack2::Details::JACOBI) { 01009 precType = "Block Jacobi"; 01010 } else if (PrecType_ == Ifpack2::Details::GS) { 01011 precType = "Block Gauss-Seidel"; 01012 } else if (PrecType_ == Ifpack2::Details::SGS) { 01013 precType = "Block Symmetric Gauss-Seidel"; 01014 } 01015 out << "type: " << precType << endl 01016 << "global number of rows: " << A_->getGlobalNumRows () << endl 01017 << "global number of columns: " << A_->getGlobalNumCols () << endl 01018 << "number of sweeps: " << NumSweeps_ << endl 01019 << "damping factor: " << DampingFactor_ << endl 01020 << "backwards mode: " 01021 << ((PrecType_ == Ifpack2::Details::GS && DoBackwardGS_) ? "true" : "false") 01022 << endl 01023 << "zero starting solution: " 01024 << (ZeroStartingSolution_ ? "true" : "false") << endl 01025 << "condition number estimate: " << Condest_ << endl; 01026 01027 out << "===============================================================================" << endl; 01028 out << "Phase # calls Total Time (s) Total MFlops MFlops/s " << endl; 01029 out << "------------ ------- --------------- --------------- ---------------" << endl; 01030 out << setw(12) << "initialize()" << setw(5) << getNumInitialize() << " " << setw(15) << getInitializeTime() << endl; 01031 out << setw(12) << "compute()" << setw(5) << getNumCompute() << " " << setw(15) << getComputeTime() << " " 01032 << setw(15) << getComputeFlops() << " " 01033 << setw(15) << (getComputeTime() != 0.0 ? getComputeFlops() / getComputeTime() * 1.0e-6 : 0.0) << endl; 01034 out << setw(12) << "apply()" << setw(5) << getNumApply() << " " << setw(15) << getApplyTime() << " " 01035 << setw(15) << getApplyFlops() << " " 01036 << setw(15) << (getApplyTime() != 0.0 ? getApplyFlops() / getApplyTime() * 1.0e-6 : 0.0) << endl; 01037 out << "===============================================================================" << endl; 01038 out << endl; 01039 } 01040 } 01041 01042 }//namespace Ifpack2 01043 01044 01045 #ifdef HAVE_IFPACK2_EXPLICIT_INSTANTIATION 01046 01047 // For ETI 01048 #include "Ifpack2_DenseContainer_decl.hpp" 01049 #include "Ifpack2_SparseContainer_decl.hpp" 01050 #include "Ifpack2_ILUT_decl.hpp" 01051 01052 // FIXME (mfh 16 Sep 2014) We should really only use RowMatrix here! 01053 // There's no need to instantiate for CrsMatrix too. All Ifpack2 01054 // preconditioners can and should do dynamic casts if they need a type 01055 // more specific than RowMatrix. 01056 01057 #define IFPACK2_BLOCKRELAXATION_INSTANT(S,LO,GO,N) \ 01058 template \ 01059 class Ifpack2::BlockRelaxation< \ 01060 Tpetra::RowMatrix<S, LO, GO, N>, \ 01061 Ifpack2::SparseContainer< \ 01062 Tpetra::RowMatrix<S, LO, GO, N>, \ 01063 Ifpack2::ILUT< ::Tpetra::RowMatrix<S,LO,GO,N> > > >; \ 01064 template \ 01065 class Ifpack2::BlockRelaxation< \ 01066 Tpetra::RowMatrix<S, LO, GO, N>, \ 01067 Ifpack2::DenseContainer< \ 01068 Tpetra::RowMatrix<S, LO, GO, N>, \ 01069 S > >; \ 01070 template \ 01071 class Ifpack2::BlockRelaxation< \ 01072 Tpetra::CrsMatrix<S, LO, GO, N>, \ 01073 Ifpack2::SparseContainer< \ 01074 Tpetra::CrsMatrix<S, LO, GO, N>, \ 01075 Ifpack2::ILUT< ::Tpetra::CrsMatrix<S,LO,GO,N> > > >; \ 01076 template \ 01077 class Ifpack2::BlockRelaxation< \ 01078 Tpetra::CrsMatrix<S, LO, GO, N>, \ 01079 Ifpack2::DenseContainer< \ 01080 Tpetra::CrsMatrix<S, LO, GO, N>, \ 01081 S > >; 01082 01083 01084 01085 #endif // HAVE_IFPACK2_EXPLICIT_INSTANTIATION 01086 01087 #endif // IFPACK2_BLOCKRELAXATION_DEF_HPP
1.7.6.1