|
GlobiPack Package Browser (Single Doxygen Collection)
Version of the Day
|
00001 /* 00002 // @HEADER 00003 // *********************************************************************** 00004 // 00005 // GlobiPack: Collection of Scalar 1D globalizaton utilities 00006 // Copyright (2009) 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 #include "GlobiPack_GoldenQuadInterpBracket.hpp" 00046 #include "GlobiPack_TestLagrPolyMeritFunc1D.hpp" 00047 #include "Teuchos_Tuple.hpp" 00048 #include "Teuchos_UnitTestHarness.hpp" 00049 00050 00051 namespace { 00052 00053 00054 // 00055 // Helper code and declarations 00056 // 00057 00058 00059 using GlobiPack::TestLagrPolyMeritFunc1D; 00060 using GlobiPack::testLagrPolyMeritFunc1D; 00061 using GlobiPack::GoldenQuadInterpBracket; 00062 using GlobiPack::goldenQuadInterpBracket; 00063 using GlobiPack::PointEval1D; 00064 using GlobiPack::computePoint; 00065 using Teuchos::as; 00066 using Teuchos::inOutArg; 00067 using Teuchos::outArg; 00068 using Teuchos::null; 00069 using Teuchos::RCP; 00070 using Teuchos::rcpFromRef; 00071 using Teuchos::Array; 00072 using Teuchos::tuple; 00073 using Teuchos::ParameterList; 00074 using Teuchos::parameterList; 00075 00076 00077 template<class Scalar> 00078 inline Scalar sqr(const Scalar &x) { return x*x; } 00079 00080 00081 // Set up a quadratic merit function with minimizer at alpha=2.0, phi=3.0; 00082 template<class Scalar> 00083 const RCP<TestLagrPolyMeritFunc1D<Scalar> > quadPhi() 00084 { 00085 typedef Teuchos::ScalarTraits<Scalar> ST; 00086 typedef typename ST::magnitudeType ScalarMag; 00087 Array<Scalar> alphaPoints = tuple<Scalar>(0.0, 2.0, 4.0); 00088 Array<ScalarMag> phiPoints = tuple<ScalarMag>(6.0, 3.0, 6.0); 00089 return testLagrPolyMeritFunc1D<Scalar>(alphaPoints, phiPoints); 00090 } 00091 00092 00093 double g_tol = Teuchos::ScalarTraits<double>::eps()*100.0; 00094 00095 00096 TEUCHOS_STATIC_SETUP() 00097 { 00098 Teuchos::UnitTestRepository::getCLP().setOption( 00099 "tol", &g_tol, "Floating point tolerance" ); 00100 } 00101 00102 00103 // 00104 // Unit tests for GoldenQuadInterpBracket 00105 // 00106 00107 00108 // 00109 // Check that object can exactly interplate a quadratic merit function at the 00110 // very first iteration 00111 // 00112 00113 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( GoldenQuadInterpBracket, bracket, Scalar ) 00114 { 00115 00116 typedef Teuchos::ScalarTraits<Scalar> ST; 00117 using std::min; 00118 00119 const RCP<TestLagrPolyMeritFunc1D<Scalar> > phi = quadPhi<Scalar>(); 00120 00121 RCP<GoldenQuadInterpBracket<Scalar> > bracket = goldenQuadInterpBracket<Scalar>(); 00122 00123 bracket->setOStream(rcpFromRef(out)); 00124 00125 const Array<Scalar> alpha = 00126 tuple<Scalar>( 00127 1e-14, 1e-10, 1e-7, 1e-4, 0.1, 1.0, 1.1, 1.5, 00128 1.9, 2.0, 2.1, 4.0, 8.0, 30.0); 00129 00130 for (int i = 0; i < as<int>(alpha.size()); ++i ) { 00131 00132 PointEval1D<Scalar> p_l = computePoint(*phi, ST::zero()); 00133 PointEval1D<Scalar> p_m = computePoint(*phi, alpha[i]); 00134 PointEval1D<Scalar> p_u; 00135 int numIters = -1; 00136 00137 if (alpha[i] > ST::eps()) { 00138 00139 const bool bracketResult = bracket->bracketMinimum( 00140 *phi, inOutArg(p_l), inOutArg(p_m), outArg(p_u), outArg(numIters) ); 00141 00142 TEST_ASSERT(bracketResult); 00143 TEST_COMPARE(p_l.alpha, <, p_m.alpha); 00144 TEST_COMPARE(p_m.alpha, <, p_u.alpha); 00145 TEST_COMPARE(p_l.phi, >, p_m.phi); 00146 TEST_COMPARE(p_m.phi, <, p_u.phi); 00147 00148 } 00149 else { 00150 00151 // If alpha[i] < eps, then there will not be enough accuracy in the 00152 // merit function to allow for a successfull line search. 00153 00154 } 00155 00156 } 00157 00158 } 00159 00160 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES( GoldenQuadInterpBracket, bracket ) 00161 00162 00163 } // namespace
1.7.6.1