|
Stratimikos Package Browser (Single Doxygen Collection)
Version of the Day
|
00001 /* 00002 // @HEADER 00003 // *********************************************************************** 00004 // 00005 // Stratimikos: Thyra-based strategies for linear solvers 00006 // Copyright (2006) Sandia Corporation 00007 // 00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00009 // license for use of this work by or on behalf of the U.S. Government. 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 Roscoe A. Bartlett (rabartl@sandia.gov) 00039 // 00040 // *********************************************************************** 00041 // @HEADER 00042 */ 00043 00044 00045 #ifndef THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP 00046 #define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP 00047 00048 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp" 00049 #include "Thyra_VectorSpaceBase.hpp" 00050 #include "Thyra_MultiVectorBase.hpp" 00051 #include "Thyra_VectorBase.hpp" 00052 #include "Thyra_MultiVectorStdOps.hpp" 00053 #include "Thyra_VectorStdOps.hpp" 00054 #include "BelosThyraAdapter.hpp" 00055 #include "Teuchos_DebugDefaultAsserts.hpp" 00056 #include "Teuchos_VerbosityLevel.hpp" 00057 #include "Teuchos_as.hpp" 00058 00059 00060 namespace Thyra { 00061 00062 00063 // Constructors/initializers/accessors 00064 00065 00066 template<class Scalar> 00067 GeneralSolveCriteriaBelosStatusTest<Scalar>::GeneralSolveCriteriaBelosStatusTest() 00068 :convergenceTestFrequency_(-1), 00069 compute_x_(false), 00070 compute_r_(false), 00071 lastCurrIter_(-1), 00072 lastRtnStatus_(Belos::Undefined) 00073 { 00074 GeneralSolveCriteriaBelosStatusTest<Scalar>::reset(); 00075 } 00076 00077 00078 template<class Scalar> 00079 void GeneralSolveCriteriaBelosStatusTest<Scalar>::setSolveCriteria( 00080 const SolveCriteria<Scalar> &solveCriteria, 00081 const int convergenceTestFrequency 00082 ) 00083 { 00084 00085 using Teuchos::as; 00086 typedef ScalarTraits<ScalarMag> SMT; 00087 00088 // Make sure the solve criteria is valid 00089 00090 TEUCHOS_ASSERT_INEQUALITY(solveCriteria.requestedTol, >=, SMT::zero()); 00091 TEUCHOS_ASSERT_INEQUALITY(solveCriteria.solveMeasureType.numerator, !=, SOLVE_MEASURE_ONE); 00092 TEUCHOS_ASSERT(nonnull(solveCriteria.numeratorReductionFunc) || 00093 nonnull(solveCriteria.denominatorReductionFunc) ); 00094 00095 // Remember the solve criteria sett 00096 00097 solveCriteria_ = solveCriteria; 00098 convergenceTestFrequency_ = convergenceTestFrequency; 00099 00100 // Determine what needs to be computed on each check 00101 00102 compute_r_ = solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_RESIDUAL); 00103 00104 compute_x_ = (compute_r_ || 00105 solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_SOLUTION)); 00106 00107 } 00108 00109 00110 template<class Scalar> 00111 ArrayView<const typename ScalarTraits<Scalar>::magnitudeType> 00112 GeneralSolveCriteriaBelosStatusTest<Scalar>::achievedTol() const 00113 { 00114 return lastAchievedTol_; 00115 } 00116 00117 00118 // Overridden public functions from Belos::StatusTest 00119 00120 00121 template <class Scalar> 00122 Belos::StatusType 00123 GeneralSolveCriteriaBelosStatusTest<Scalar>::checkStatus( 00124 Belos::Iteration<Scalar,MV,OP> *iSolver 00125 ) 00126 { 00127 00128 using Teuchos::null; 00129 00130 #ifdef THYRA_DEBUG 00131 TEUCHOS_ASSERT(iSolver); 00132 #endif 00133 00134 const int currIter = iSolver->getNumIters(); 00135 00136 if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) { 00137 // We will skip this check! 00138 return Belos::Undefined; 00139 } 00140 00141 const RCP<FancyOStream> out = this->getOStream(); 00142 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel(); 00143 00144 const Belos::LinearProblem<Scalar,MV,OP>& lp = iSolver->getProblem(); 00145 const int numRhs = lp.getRHS()->domain()->dim(); 00146 00147 // Compute the rhs norm if requested and not already computed 00148 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_RHS)) { 00149 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||b||"); 00150 } 00151 00152 // Compute the initial residual norm if requested and not already computed 00153 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_INIT_RESIDUAL)) { 00154 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||r0||"); 00155 } 00156 00157 // Compute X if requested 00158 RCP<MultiVectorBase<Scalar> > X; 00159 if (compute_x_) { 00160 RCP<MV> X_update = iSolver->getCurrentUpdate(); 00161 X = lp.updateSolution(X_update); 00162 } 00163 00164 // Compute R if requested 00165 RCP<MultiVectorBase<Scalar> > R; 00166 if (compute_r_) { 00167 R = createMembers(lp.getOperator()->range(), X->domain()); 00168 lp.computeCurrResVec(&*R, &*X); 00169 } 00170 00171 // Form numerators and denominaotors gN(vN)/gD(vD) for each system 00172 00173 lastNumerator_.resize(numRhs); 00174 lastDenominator_.resize(numRhs); 00175 00176 for (int j = 0; j < numRhs; ++j) { 00177 const RCP<const VectorBase<Scalar> > x_j = (nonnull(X) ? X->col(j) : null); 00178 const RCP<const VectorBase<Scalar> > r_j = (nonnull(R) ? R->col(j) : null); 00179 lastNumerator_[j] = computeReductionFunctional( 00180 solveCriteria_.solveMeasureType.numerator, 00181 solveCriteria_.numeratorReductionFunc.ptr(), 00182 x_j.ptr(), r_j.ptr() ); 00183 lastDenominator_[j] = computeReductionFunctional( 00184 solveCriteria_.solveMeasureType.denominator, 00185 solveCriteria_.denominatorReductionFunc.ptr(), 00186 x_j.ptr(), r_j.ptr() ); 00187 } 00188 00189 // Determine if convRatio <= requestedTol 00190 00191 bool systemsAreConverged = true; 00192 lastAchievedTol_.resize(numRhs); 00193 00194 for (int j = 0; j < numRhs; ++j) { 00195 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j]; 00196 lastAchievedTol_[j] = convRatio; 00197 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol); 00198 if (includesVerbLevel(verbLevel, Teuchos::VERB_MEDIUM)) { 00199 printRhsStatus(currIter, j, *out); 00200 } 00201 if (!sys_converged_j) { 00202 systemsAreConverged = false; 00203 } 00204 } 00205 00206 lastRtnStatus_ = (systemsAreConverged ? Belos::Passed : Belos::Failed); 00207 lastCurrIter_ = currIter; 00208 00209 return lastRtnStatus_; 00210 00211 } 00212 00213 template <class Scalar> 00214 Belos::StatusType 00215 GeneralSolveCriteriaBelosStatusTest<Scalar>::getStatus() const 00216 { 00217 return lastRtnStatus_; 00218 } 00219 00220 00221 template <class Scalar> 00222 void GeneralSolveCriteriaBelosStatusTest<Scalar>::reset() 00223 { 00224 r0_nrm_.clear(); 00225 b_nrm_.clear(); 00226 lastNumerator_.clear(); 00227 lastDenominator_.clear(); 00228 lastAchievedTol_.clear(); 00229 lastCurrIter_ = -1; 00230 lastRtnStatus_ = Belos::Undefined; 00231 } 00232 00233 00234 template <class Scalar> 00235 void GeneralSolveCriteriaBelosStatusTest<Scalar>::print( 00236 std::ostream& os, int indent 00237 ) const 00238 { 00239 const int numRhs = lastNumerator_.size(); 00240 for (int j = 0; j < numRhs; ++j) { 00241 printRhsStatus(lastCurrIter_, j, os, indent); 00242 } 00243 } 00244 00245 00246 // private 00247 00248 00249 template <class Scalar> 00250 typename GeneralSolveCriteriaBelosStatusTest<Scalar>::ScalarMag 00251 GeneralSolveCriteriaBelosStatusTest<Scalar>::computeReductionFunctional( 00252 ESolveMeasureNormType measureType, 00253 const Ptr<const ReductionFunctional<Scalar> > &reductFunc, 00254 const Ptr<const VectorBase<Scalar> > &x, 00255 const Ptr<const VectorBase<Scalar> > &r 00256 ) const 00257 { 00258 typedef ScalarTraits<ScalarMag> SMT; 00259 ScalarMag rtn = -SMT::one(); 00260 Ptr<const VectorBase<Scalar> > v; 00261 switch(measureType) { 00262 case SOLVE_MEASURE_ONE: 00263 rtn = SMT::one(); 00264 break; 00265 case SOLVE_MEASURE_NORM_RESIDUAL: 00266 v = r; 00267 break; 00268 case SOLVE_MEASURE_NORM_SOLUTION: 00269 v = x; 00270 break; 00271 case SOLVE_MEASURE_NORM_INIT_RESIDUAL: 00272 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||r0||!)") 00273 break; 00274 case SOLVE_MEASURE_NORM_RHS: 00275 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||b||!)"); 00276 break; 00277 TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT(); 00278 } 00279 if (rtn >= SMT::zero()) { 00280 // We already have what we need! 00281 } 00282 else if (nonnull(v) && rtn < SMT::zero()) { 00283 if (nonnull(reductFunc)) { 00284 rtn = reductFunc->reduce(*v); 00285 } 00286 else { 00287 rtn = norm(*v); 00288 } 00289 } 00290 TEUCHOS_IF_ELSE_DEBUG_ASSERT(); 00291 return rtn; 00292 } 00293 00294 00295 template <class Scalar> 00296 void 00297 GeneralSolveCriteriaBelosStatusTest<Scalar>::printRhsStatus( 00298 const int currIter, const int j, std::ostream &out, 00299 int indent 00300 ) const 00301 { 00302 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j]; 00303 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol); 00304 for (int i = 0; i < indent; ++i) { out << " "; } 00305 out 00306 << "["<<currIter<<"] " 00307 << "gN(vN("<<j<<"))/gD(vD("<<j<<")) = " 00308 << lastNumerator_[j] << "/" << lastDenominator_[j] << " = " 00309 << convRatio << " <= " << solveCriteria_.requestedTol << " : " 00310 << (sys_converged_j ? " true" : "false") 00311 << "\n"; 00312 } 00313 00314 00315 } // namespace Thyra 00316 00317 00318 #endif // THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
1.7.6.1