|
Tpetra Matrix/Vector Services
Version of the Day
|
00001 /* 00002 // @HEADER 00003 // *********************************************************************** 00004 // 00005 // Tpetra: Templated Linear Algebra Services Package 00006 // Copyright (2008) Sandia Corporation 00007 // 00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00009 // the U.S. Government retains certain rights in this software. 00010 // 00011 // Redistribution and use in source and binary forms, with or without 00012 // modification, are permitted provided that the following conditions are 00013 // met: 00014 // 00015 // 1. Redistributions of source code must retain the above copyright 00016 // notice, this list of conditions and the following disclaimer. 00017 // 00018 // 2. Redistributions in binary form must reproduce the above copyright 00019 // notice, this list of conditions and the following disclaimer in the 00020 // documentation and/or other materials provided with the distribution. 00021 // 00022 // 3. Neither the name of the Corporation nor the names of the 00023 // contributors may be used to endorse or promote products derived from 00024 // this software without specific prior written permission. 00025 // 00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00037 // 00038 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00039 // 00040 // ************************************************************************ 00041 // @HEADER 00042 */ 00043 00044 #ifndef MULTIPREC_DRIVER_HPP 00045 #define MULTIPREC_DRIVER_HPP 00046 00047 #include <Teuchos_TypeNameTraits.hpp> 00048 #include <Teuchos_TimeMonitor.hpp> 00049 #include <Teuchos_ParameterList.hpp> 00050 #include <Teuchos_FancyOStream.hpp> 00051 00052 #include <Tpetra_Version.hpp> 00053 #include <Tpetra_MatrixIO.hpp> 00054 #include <Tpetra_RTI.hpp> 00055 #include <Tpetra_CrsMatrix.hpp> 00056 #include <Tpetra_Vector.hpp> 00057 00058 #include "MultiPrecCG.hpp" 00059 00060 template <class MPStack> 00061 class MultiPrecDriver { 00062 public: 00063 // input 00064 Teuchos::RCP<Teuchos::FancyOStream> out; 00065 Teuchos::RCP<Teuchos::ParameterList> params; 00066 std::string matrixFile; 00067 bool unfusedTest; 00068 // output 00069 bool testPassed; 00070 00071 template <class Node> 00072 void run(Teuchos::ParameterList &myMachPL, const Teuchos::RCP<const Teuchos::Comm<int> > &comm, const Teuchos::RCP<Node> &node) 00073 { 00074 using std::pair; 00075 using std::make_pair; 00076 using std::plus; 00077 using std::endl; 00078 using Teuchos::null; 00079 using Teuchos::RCP; 00080 using Teuchos::ParameterList; 00081 using Teuchos::parameterList; 00082 using TpetraExamples::make_pair_op; 00083 using Tpetra::RTI::reductionGlob; 00084 using Tpetra::RTI::ZeroOp; 00085 using Tpetra::RTI::binary_pre_transform_reduce; 00086 using Tpetra::RTI::binary_transform; 00087 00088 // Static types 00089 typedef typename MPStack::type S; 00090 typedef int LO; 00091 typedef int GO; 00092 typedef Tpetra::Map<LO,GO,Node> Map; 00093 typedef Tpetra::CrsMatrix<S,LO,GO,Node> CrsMatrix; 00094 typedef Tpetra::Vector<S,LO,GO,Node> Vector; 00095 00096 *out << "Running test with Node==" << Teuchos::typeName(*node) << " on rank " << comm->getRank() << "/" << comm->getSize() << std::endl; 00097 00098 // read the matrix 00099 RCP<CrsMatrix> A; 00100 RCP<const Map> rowMap = null; 00101 RCP<ParameterList> fillParams = parameterList(); 00102 fillParams->set("Preserve Local Graph",true); 00103 // must preserve the local graph in order to do convert() calls later 00104 Tpetra::Utils::readHBMatrix(matrixFile,comm,node,A,rowMap,fillParams); 00105 rowMap = A->getRowMap(); 00106 00107 // init the solver stack 00108 TpetraExamples::RFPCGInit<S,LO,GO,Node> init(A); 00109 RCP<ParameterList> db = Tpetra::Ext::initStackDB<MPStack>(*params,init); 00110 00111 testPassed = true; 00112 00113 // choose a solution, compute a right-hand-side 00114 auto x = Tpetra::createVector<S>(rowMap), 00115 b = Tpetra::createVector<S>(rowMap); 00116 x->randomize(); 00117 A->apply(*x,*b); 00118 { 00119 // init the rhs 00120 auto bx = db->get<RCP<Vector>>("bx"); 00121 binary_transform( *bx, *b, [](S, S bi) {return bi;}); // bx = b 00122 } 00123 00124 // call the solve 00125 TpetraExamples::recursiveFPCG<MPStack,LO,GO,Node>(out,*db); 00126 00127 // check that residual is as requested 00128 { 00129 auto xhat = db->get<RCP<Vector>>("bx"), 00130 bhat = Tpetra::createVector<S>(rowMap); 00131 A->apply(*xhat,*bhat); 00132 // compute bhat-b, while simultaneously computing |bhat-b|^2 and |b|^2 00133 auto nrms = binary_pre_transform_reduce(*bhat, *b, 00134 reductionGlob<ZeroOp<pair<S,S>>>( 00135 [](S bhati, S bi){ return bi-bhati;}, // bhati = bi-bhat 00136 [](S bhati, S bi){ return make_pair(bhati*bhati, bi*bi); }, 00137 make_pair_op<S,S>(plus<S>())) ); 00138 const S enrm = Teuchos::ScalarTraits<S>::squareroot(nrms.first), 00139 bnrm = Teuchos::ScalarTraits<S>::squareroot(nrms.second); 00140 // check that residual is as requested 00141 *out << "|b - A*x|/|b|: " << enrm / bnrm << endl; 00142 const double tolerance = db->get<double>("tolerance"); 00143 if (MPStack::bottom) { 00144 // give a little slack 00145 if (enrm / bnrm > 5*tolerance) testPassed = false; 00146 } 00147 else { 00148 if (enrm / bnrm > tolerance) testPassed = false; 00149 } 00150 } 00151 00152 // 00153 // solve again, with the unfused version, just for timings purposes 00154 if (unfusedTest) 00155 { 00156 // init the rhs 00157 auto bx = db->get<RCP<Vector>>("bx"); 00158 binary_transform( *bx, *b, [](S, S bi) {return bi;}); // bx = b 00159 // call the solve 00160 TpetraExamples::recursiveFPCGUnfused<MPStack,LO,GO,Node>(out,*db); 00161 // 00162 // test the result 00163 auto xhat = db->get<RCP<Vector>>("bx"), 00164 bhat = Tpetra::createVector<S>(rowMap); 00165 A->apply(*xhat,*bhat); 00166 // compute bhat-b, while simultaneously computing |bhat-b|^2 and |b|^2 00167 auto nrms = binary_pre_transform_reduce(*bhat, *b, 00168 reductionGlob<ZeroOp<pair<S,S>>>( 00169 [](S bhati, S bi){ return bi-bhati;}, // bhati = bi-bhat 00170 [](S bhati, S bi){ return make_pair(bhati*bhati, bi*bi); }, 00171 make_pair_op<S,S>(plus<S>())) ); 00172 const S enrm = Teuchos::ScalarTraits<S>::squareroot(nrms.first), 00173 bnrm = Teuchos::ScalarTraits<S>::squareroot(nrms.second); 00174 // check that residual is as requested 00175 *out << "|b - A*x|/|b|: " << enrm / bnrm << endl; 00176 const double tolerance = db->get<double>("tolerance"); 00177 if (MPStack::bottom) { 00178 // give a little slack 00179 if (enrm / bnrm > 5*tolerance) testPassed = false; 00180 } 00181 else { 00182 if (enrm / bnrm > tolerance) testPassed = false; 00183 } 00184 } 00185 00186 00187 // print timings 00188 Teuchos::TimeMonitor::summarize( *out ); 00189 } 00190 }; 00191 00192 #endif // MULTIPREC_DRIVER_HPP
1.7.6.1