|
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_TestLagrPolyMeritFunc1D.hpp" 00046 #include "GlobiPack_BrentsLineSearch.hpp" 00047 #include "Teuchos_Tuple.hpp" 00048 00049 #include "meritFuncsHelpers.hpp" 00050 00051 #include "Teuchos_UnitTestHarness.hpp" 00052 00053 00054 namespace { 00055 00056 00057 // 00058 // Helper code and declarations 00059 // 00060 00061 00062 using GlobiPack::BrentsLineSearch; 00063 using GlobiPack::brentsLineSearch; 00064 using GlobiPack::computeValue; 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 double g_tol_scale = 100.0; 00078 00079 00080 TEUCHOS_STATIC_SETUP() 00081 { 00082 Teuchos::UnitTestRepository::getCLP().setOption( 00083 "tol", &g_tol_scale, "Floating point tolerance scaling of eps." ); 00084 } 00085 00086 00087 // 00088 // Unit tests for BrentsLineSearch 00089 // 00090 00091 00092 // 00093 // Check that object can exactly interplate a quadratic merit function. This 00094 // takes more than one iteration because of the golden search stuff. 00095 // 00096 00097 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( BrentsLineSearch, quadExact, Scalar ) 00098 { 00099 00100 typedef Teuchos::ScalarTraits<Scalar> ST; 00101 typedef typename ST::magnitudeType ScalarMag; 00102 00103 const RCP<TestLagrPolyMeritFunc1D<Scalar> > phi = quadPhi<Scalar>(); 00104 00105 RCP<BrentsLineSearch<Scalar> > linesearch = brentsLineSearch<Scalar>(); 00106 00107 linesearch->setOStream(rcpFromRef(out)); 00108 00109 const PointEval1D<Scalar> point_k = computePoint(*phi, ST::zero()); 00110 PointEval1D<Scalar> point_kp1 = computePoint(*phi, as<Scalar>(8.0)); 00111 00112 int numIters = -1; 00113 const bool linesearchResult = linesearch->doLineSearch( 00114 *phi, point_k, inOutArg(point_kp1), outArg(numIters) ); 00115 00116 TEST_ASSERT(linesearchResult); 00117 TEST_FLOATING_EQUALITY(point_kp1.alpha, as<Scalar>(2.0), 00118 as<Scalar>(g_tol_scale*ST::squareroot(ST::eps()))); 00119 TEST_FLOATING_EQUALITY(point_kp1.phi, as<ScalarMag>(3.0), 00120 as<Scalar>(g_tol_scale)*ST::eps()); 00121 00122 } 00123 00124 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES( BrentsLineSearch, quadExact ) 00125 00126 00127 // 00128 // Check that object can approximately mimimize a cubic function. 00129 // 00130 00131 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( BrentsLineSearch, cubicApprox, Scalar ) 00132 { 00133 00134 typedef Teuchos::ScalarTraits<Scalar> ST; 00135 typedef typename ST::magnitudeType ScalarMag; 00136 00137 const RCP<TestLagrPolyMeritFunc1D<Scalar> > phi = quadPhi<Scalar>(); 00138 00139 RCP<BrentsLineSearch<Scalar> > linesearch = brentsLineSearch<Scalar>(); 00140 00141 const RCP<ParameterList> pl = parameterList(); 00142 pl->sublist("Minimize").set("Relative Tol", as<double>(g_tol_scale*ST::eps())); 00143 pl->sublist("Minimize").set("Bracket Tol", as<double>(ST::eps())); 00144 linesearch->setParameterList(pl); 00145 00146 linesearch->setOStream(rcpFromRef(out)); 00147 00148 const PointEval1D<Scalar> point_k = computePoint(*phi, ST::zero()); 00149 PointEval1D<Scalar> point_kp1 = computePoint(*phi, as<Scalar>(8.0)); 00150 00151 int numIters = -1; 00152 const bool linesearchResult = linesearch->doLineSearch( 00153 *phi, point_k, inOutArg(point_kp1), outArg(numIters) ); 00154 00155 TEST_ASSERT(linesearchResult); 00156 TEST_FLOATING_EQUALITY(point_kp1.alpha, as<Scalar>(2.0), 00157 as<Scalar>(g_tol_scale*ST::squareroot(ST::eps()))); 00158 TEST_FLOATING_EQUALITY(point_kp1.phi, as<ScalarMag>(3.0), 00159 as<Scalar>(g_tol_scale)*ST::eps()); 00160 00161 } 00162 00163 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES( BrentsLineSearch, cubicApprox ) 00164 00165 00166 } // namespace
1.7.6.1