lnsrch.hh

Go to the documentation of this file.
00001 // lnsrch.H
00002 // created by ADP
00003 // last modified 06/17/04
00004 // considerable modification 2007-2011 WWS
00005 
00006 /*************************************************************************
00007 
00008 Copyright Rice University, 2004.
00009 All rights reserved.
00010 
00011 Permission is hereby granted, free of charge, to any person obtaining a
00012 copy of this software and associated documentation files (the "Software"),
00013 to deal in the Software without restriction, including without limitation
00014 the rights to use, copy, modify, merge, publish, distribute, and/or sell
00015 copies of the Software, and to permit persons to whom the Software is
00016 furnished to do so, provided that the above copyright notice(s) and this
00017 permission notice appear in all copies of the Software and that both the
00018 above copyright notice(s) and this permission notice appear in supporting
00019 documentation.
00020 
00021 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
00022 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
00023 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT OF THIRD PARTY
00024 RIGHTS. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR HOLDERS INCLUDED IN THIS
00025 NOTICE BE LIABLE FOR ANY CLAIM, OR ANY SPECIAL INDIRECT OR CONSEQUENTIAL
00026 DAMAGES, OR ANY DAMAGES WHATSOEVER RESULTING FROM LOSS OF USE, DATA OR
00027 PROFITS, WHETHER IN AN ACTION OF CONTRACT, NEGLIGENCE OR OTHER TORTIOUS
00028 ACTION, ARISING OUT OF OR IN CONNECTION WITH THE USE OR PERFORMANCE OF
00029 THIS SOFTWARE.
00030 
00031 Except as contained in this notice, the name of a copyright holder shall
00032 not be used in advertising or otherwise to promote the sale, use or other
00033 dealings in this Software without prior written authorization of the
00034 copyright holder.
00035 
00036 **************************************************************************/
00037 
00038 
00039 #ifndef __RVL_ALG_LINE_SEARCH_
00040 #define __RVL_ALG_LINE_SEARCH_
00041 
00042 #include "alg.hh"
00043 #include "functional.hh"
00044 #include "terminator.hh"
00045 #include "table.hh"
00046 
00047 namespace RVLUmin {
00048 
00049   using namespace RVL;
00050   using namespace RVLAlg;
00051 
00052   // forward declaration;
00053   template<typename Scalar>
00054   class LineSearchAlg;
00055 
00088   template <class Scalar>
00089   class LineSearchAlgBase: public Algorithm, public Terminator {
00090 
00091   private:
00092 
00093     Scalar f0;
00094     Vector<Scalar> x0;
00095     Vector<Scalar> g0;
00096     LineSearchAlg<Scalar> const & lsalg;
00097 
00098     LineSearchAlgBase();
00099 
00100   protected:
00101 
00102     FunctionalEvaluation<Scalar> & fx;
00103     Vector<Scalar> const & dx;
00104 
00105     // virtual copy constructor, in order that LSAlg can have real one
00106     virtual LineSearchAlgBase<Scalar> * clone() const = 0;
00107 
00108   public:
00109 
00110     LineSearchAlgBase(LineSearchAlg<Scalar> const & _lsalg,
00111               FunctionalEvaluation<Scalar> & _fx,
00112               Vector<Scalar> const & _dx)
00113     : 
00114     f0(_fx.getValue()),
00115     x0(_fx.getPoint()),
00116     g0(_fx.getGradient()),
00117     lsalg(_lsalg),
00118     fx(_fx), 
00119     dx(_dx) {} //cerr<<"from LineSearchAlgBase constructor\n";}
00120     LineSearchAlgBase(const LineSearchAlgBase<Scalar> & ls)
00121     : f0(ls.f0), x0(ls.x0), g0(ls.g0), lsalg(ls.lsalg), fx(ls.fx), dx(ls.dx) {}
00122     virtual ~LineSearchAlgBase() {}
00123 
00124     virtual ostream & write(ostream & str) const {
00125       str<<"lnsrch base class"<<endl;
00126       return str;
00127     }
00128 
00129     // note that Terminator::query and Algorithm::run are abstract
00130 
00131     virtual Scalar getStep() const = 0;
00132     Scalar getBaseValue() const { return f0; }
00133     Vector<Scalar> const & getBasePoint() const { 
00134       return x0; }
00135     Vector<Scalar> const & getBaseGradient() const { 
00136       return g0; }
00137     Vector<Scalar> & getTrialPoint() const { 
00138       return fx.getPoint(); }
00139     Vector<Scalar> const & getSearchDirection() const { 
00140       return dx; }
00141     FunctionalEvaluation<Scalar> & getFunctionalEvaluation() const { 
00142       return fx; }
00143     //    bool checkMinStep() const { return lsalg.checkMinStep(); }
00144     Scalar getMinStep() const { return lsalg.getMinStep(); }
00145   };
00146 
00159   template<typename Scalar>
00160   class LineSearchAlg: public Algorithm, public Terminator {
00161   private:
00162     // the actual algorithm
00163     mutable LineSearchAlgBase<Scalar> * ls;
00164     // step used in last successful line search, initialized somehow
00165     mutable Scalar firststep;
00166     // minimum permitted step
00167     Scalar minsteptol;
00168 
00169   protected:
00175     virtual LineSearchAlgBase<Scalar> * 
00176     build(FunctionalEvaluation<Scalar> & fx,
00177       Vector<Scalar> const & dx,
00178       Scalar firststep) = 0;
00179 
00180   public:
00185     LineSearchAlg(Scalar _firststep = numeric_limits<Scalar>::max(),
00186           Scalar _minsteptol = numeric_limits<Scalar>::epsilon())
00187       : ls(NULL), firststep(_firststep), minsteptol(_minsteptol) {}
00188     LineSearchAlg(LineSearchAlg<Scalar> const & lsh)
00189       : ls(NULL), firststep(lsh.firststep), minsteptol(lsh.minsteptol) {
00190       if (lsh.ls) ls=(lsh.ls)->clone();
00191     }
00192     ~LineSearchAlg() {
00193       if (ls) delete ls;
00194     }
00195 
00197     void initialize(FunctionalEvaluation<Scalar> & fx,
00198             Vector<Scalar> const & dx) {
00199       if (ls) delete ls;
00200       ls = build(fx,dx,firststep);
00201     }
00202 
00204     void run() {
00205       if (ls) { 
00206     ls->run();
00207     firststep=ls->getStep();
00208     return;
00209       }
00210       RVLException e;
00211       e<<"Error: LineSearchAlg::run\n";
00212       e<<"this alg not initialized\n";
00213       throw e;
00214     }
00215 
00217     bool query() { 
00218       if (ls) {
00219     return ls->query();
00220       }
00221       return false;
00222     }
00223 
00224     Scalar getStep() const { return ls->getStep(); }
00225 
00226     Scalar getBaseValue() const { 
00227       if (ls) return ls->getBaseValue();
00228       RVLException e;
00229       e<<"Error: LineSearchAlg::getBaseValue()\n";
00230       e<<"this alg not initialized\n";
00231       throw e;
00232     }
00233 
00234     Vector<Scalar> const & getBasePoint() const { 
00235       if (ls) return ls->getBasePoint(); 
00236       RVLException e;
00237       e<<"Error: LineSearchAlg::getBasePoint()\n";
00238       e<<"this alg not initialized\n";
00239       throw e;
00240     }
00241 
00242     Vector<Scalar> const & getBaseGradient() const { 
00243       if (ls) return ls->getBaseGradient(); 
00244       RVLException e;
00245       e<<"Error: LineSearchAlg::getBaseGradient()\n";
00246       e<<"this alg not initialized\n";
00247       throw e;
00248     }
00249 
00250     Vector<Scalar> & getTrialPoint() const { 
00251       if (ls) return ls->getTrialPoint(); 
00252       RVLException e;
00253       e<<"Error: LineSearchAlg::getTrialPoint()\n";
00254       e<<"this alg not initialized\n";
00255       throw e;
00256     }
00257 
00258     Vector<Scalar> const & getSearchDirection() const { 
00259       if (ls) return ls->getSearchDirection(); 
00260       RVLException e;
00261       e<<"Error: LineSearchAlg::getSearchDirection()\n";
00262       e<<"this alg not initialized\n";
00263       throw e;
00264     }
00265 
00266     FunctionalEvaluation<Scalar> & getFunctionalEvaluation() const {
00267       if (ls) return ls->getFunctionalEvaluation(); 
00268       RVLException e;
00269       e<<"Error: LineSearchAlg::getFunctionalEvaluation()\n";
00270       e<<"this alg not initialized\n";
00271       throw e;
00272     }
00273 
00280     /*
00281     bool checkMinStep() const {
00282       try {
00283     Scalar foo = this->getStep();
00284     return (foo > minsteptol);
00285       }
00286       catch (RVLException & e) {
00287     e<<"\ncalled from LineSearchAlg::checkMinStep\n";
00288     throw e;
00289       }
00290     }
00291     */
00292     Scalar getMinStep() const { return minsteptol; }
00293 
00294   };
00295 
00296 }
00297 
00298 #endif

Generated on 5 Jan 2017 for RvlUmin by  doxygen 1.4.7