Sacado_rad.hpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 //
00004 //                           Sacado Package
00005 //                 Copyright (2006) Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact David M. Gay (dmgay@sandia.gov) or Eric T. Phipps
00025 // (etphipp@sandia.gov).
00026 //
00027 // ***********************************************************************
00028 // @HEADER
00029 
00030 // RAD package (Reverse Automatic Differentiation) --
00031 // a package specialized for function and gradient evaluations.
00032 // Written in 2004 by David M. Gay at Sandia National Labs, Albuquerque, NM.
00033 
00034 #ifndef SACADO_RAD_H
00035 #define SACADO_RAD_H
00036 
00037 #include <stddef.h>
00038 #include <cmath>
00039 #include <math.h>
00040 
00041 namespace Sacado {
00042 namespace Radnt { // nontemplated RAD
00043 
00044 // -DRAD_NO_USING_STDCC is needed, e.g., with Sun CC 5.7
00045 #ifndef RAD_NO_USING_STDCC
00046   // Bring math functions into scope
00047   using std::exp;
00048   using std::log;
00049   using std::log10;
00050   using std::sqrt;
00051   using std::cos;
00052   using std::sin;
00053   using std::tan;
00054   using std::acos;
00055   using std::asin;
00056   using std::atan;
00057   using std::cosh;
00058   using std::sinh;
00059   using std::tanh;
00060   using std::abs;
00061   using std::fabs;
00062   using std::atan2;
00063   using std::pow;
00064 #endif 
00065 
00066  class ADvar;
00067  class ADvari;
00068  class ADvar1;
00069  class ADvar2;
00070  class ConstADvar;
00071  class Derp;
00072  class IndepADvar;
00073 
00074  extern ADvari& ADf1(double f, double g, const ADvari &x);
00075  extern ADvari& ADf2(double f, double gx, double gy, const ADvari &x, const ADvari &y);
00076  extern ADvari& ADfn(double f, int n, const ADvar *x, const double *g);
00077 
00078  struct
00079 ADmemblock {  // We get memory in ADmemblock chunks and never give it back,
00080     // but reuse it once computations start anew after call(s) on
00081     // ADcontext::Gradcomp() or ADcontext::Weighted_Gradcomp().
00082   ADmemblock *next;
00083   double memblk[1000];
00084   };
00085 
00086  class
00087 ADcontext { // A singleton class: one instance in radops.c
00088   ADmemblock *Busy, *Free;
00089   char *Mbase;
00090   size_t Mleft;
00091   ADmemblock First;
00092   void *new_ADmemblock(size_t);
00093  public:
00094   ADcontext();
00095   void *Memalloc(size_t len);
00096   static void Gradcomp();
00097   static void Weighted_Gradcomp(int, ADvar**, double*);
00098   };
00099 
00100 inline void *ADcontext::Memalloc(size_t len) {
00101     if (Mleft >= len)
00102       return Mbase + (Mleft -= len);
00103     return new_ADmemblock(len);
00104     }
00105 
00106  class
00107 CADcontext: public ADcontext {
00108  protected:
00109   bool fpval_implies_const;
00110  public:
00111   friend class ADvar;
00112   CADcontext(): ADcontext() { fpval_implies_const = false; }
00113   static const double One, negOne;
00114   };
00115 
00116  class
00117 Derp {    // one derivative-propagation operation
00118  public:
00119   friend class ADvarn;
00120   static Derp *LastDerp;
00121   Derp *next;
00122   const double *a;
00123   const ADvari *b;
00124   const ADvari *c;
00125   void *operator new(size_t);
00126   void operator delete(void*) {} /*Should never be called.*/
00127   inline Derp(){};
00128   Derp(const ADvari *);
00129   Derp(const double *, const ADvari *);
00130   Derp(const double *, const ADvari *, const ADvari *);
00131   /* c->aval += a * b->aval; */
00132   };
00133 
00134 inline Derp::Derp(const ADvari *c1): c(c1) {
00135     next = LastDerp;
00136     LastDerp = this;
00137     }
00138 
00139 inline Derp::Derp(const double *a1, const ADvari *c1): a(a1), c(c1) {
00140     next = LastDerp;
00141     LastDerp = this;
00142     }
00143 
00144 inline Derp::Derp(const double *a1, const ADvari *b1, const ADvari *c1): a(a1), b(b1), c(c1) {
00145     next = LastDerp;
00146     LastDerp = this;
00147     }
00148 
00149  class
00150 ADvari {  // implementation of an ADvar
00151  public:
00152   double Val;   // result of this operation
00153   mutable double aval;  // adjoint -- partial of final result w.r.t. this Val
00154   void *operator new(size_t len) { return ADvari::adc.Memalloc(len); }
00155   void operator delete(void*) {} /*Should never be called.*/
00156   inline ADvari(double t): Val(t), aval(0.) {}
00157   inline ADvari(double t, double ta): Val(t), aval(ta) {}
00158   inline ADvari(): Val(0.), aval(0.) {}
00159   static ADcontext adc;
00160 #ifdef RAD_AUTO_AD_Const
00161   friend class ADcontext;
00162   friend class ADvar1;
00163   friend class ADvar;
00164   friend class ConstADvar;
00165   friend class IndepADvar;
00166  private:
00167   ADvari *Next;
00168   static ADvari *First_ADvari, **Last_ADvari;
00169  protected:
00170   IndepADvar *padv;
00171  public:
00172   ADvari(const IndepADvar *, double);
00173 #endif
00174 #define F friend
00175 #define R ADvari&
00176 #define Ai const ADvari&
00177 #define T1(r,f) F r f(Ai);
00178 #define T2(r,f) \
00179 F r f(Ai,Ai); \
00180 F r f(double,Ai); \
00181 F r f(Ai,double);
00182   T1(R,operator+)
00183   T2(R,operator+)
00184   T1(R,operator-)
00185   T2(R,operator-)
00186   T2(R,operator*)
00187   T2(R,operator/)
00188   T1(R,abs)
00189   T1(R,acos)
00190   T1(R,acosh)
00191   T1(R,asin)
00192   T1(R,asinh)
00193   T1(R,atan)
00194   T1(R,atanh)
00195   T2(R,atan2)
00196   T2(R,max)
00197   T2(R,min)
00198   T1(R,cos)
00199   T1(R,cosh)
00200   T1(R,exp)
00201   T1(R,log)
00202   T1(R,log10)
00203   T2(R,pow)
00204   T1(R,sin)
00205   T1(R,sinh)
00206   T1(R,sqrt)
00207   T1(R,tan)
00208   T1(R,tanh)
00209   T1(R,fabs)
00210   T1(R,copy)
00211   T2(int,operator<)
00212   T2(int,operator<=)
00213   T2(int,operator==)
00214   T2(int,operator!=)
00215   T2(int,operator>=)
00216   T2(int,operator>)
00217 #undef T2
00218 #undef T1
00219 #undef Ai
00220 #undef R
00221 #undef F
00222 
00223   friend ADvari& ADf1(double f, double g, const ADvari &x);
00224   friend ADvari& ADf2(double f, double gx, double gy, const ADvari &x, const ADvari &y);
00225   friend ADvari& ADfn(double f, int n, const ADvar *x, const double *g);
00226   };
00227 
00228  inline void* Derp::operator new(size_t len) { return ADvari::adc.Memalloc(len); }
00229 
00230 
00231  class
00232 ADvar1: public ADvari { // simplest unary ops
00233  public:
00234   Derp d;
00235   ADvar1(double val1): ADvari(val1) {}
00236   ADvar1(double val1, const ADvari *c1): d(c1) { Val = val1; }
00237   ADvar1(double val1, const double *a1, const ADvari *c1): d(a1,this,c1) { Val = val1; }
00238 #ifdef RAD_AUTO_AD_Const
00239   ADvar1(const IndepADvar *, const IndepADvar &);
00240   ADvar1(const IndepADvar *, const ADvari &);
00241 #endif
00242   };
00243 
00244  class
00245 ConstADvari: public ADvari {
00246  private:
00247   ConstADvari *prevcad;
00248   ConstADvari() {}; // prevent construction without value (?)
00249   static ConstADvari *lastcad;
00250  public:
00251   static CADcontext cadc;
00252   inline void *operator new(size_t len) { return ConstADvari::cadc.Memalloc(len); }
00253   inline ConstADvari(double t): ADvari(t) { prevcad = lastcad; lastcad = this; }
00254   static void aval_reset(void);
00255   };
00256 
00257  class
00258 IndepADvar
00259 {
00260  private:
00261   inline IndepADvar& operator=(const IndepADvar &x) {
00262     /* private to prevent assignment */
00263 #ifdef RAD_AUTO_AD_Const
00264     if (cv)
00265       cv->padv = 0;
00266     cv = new ADvar1(this,x);
00267     return *this;
00268 #else
00269 #ifdef RAD_EQ_ALIAS
00270     cv = x.cv;
00271     return *this;
00272 #else
00273     return ADvar_operatoreq(this,*x.cv);
00274 #endif
00275 #endif /* RAD_AUTO_AD_Const */
00276     }
00277  protected:
00278   static void AD_Const(const IndepADvar&);
00279   ADvari *cv;
00280  public:
00281   typedef double value_type;
00282   friend class ADvar;
00283   friend class ADvar1;
00284   friend class ADvarn;
00285   friend class ADcontext;
00286   IndepADvar(double);
00287   IndepADvar(int);
00288   IndepADvar(long);
00289   IndepADvar& operator=(double);
00290 #ifdef RAD_AUTO_AD_Const
00291   inline IndepADvar(const IndepADvar &x) { cv = x.cv ? new ADvar1(this, x) : 0; };
00292   inline IndepADvar() { cv = 0; }
00293   inline ~IndepADvar() {
00294       if (cv)
00295         cv->padv = 0;
00296     }
00297 #else
00298   inline IndepADvar() {
00299 #ifndef RAD_EQ_ALIAS
00300     cv = 0;
00301 #endif
00302     }
00303   inline ~IndepADvar() {}
00304   friend IndepADvar& ADvar_operatoreq(IndepADvar*, const ADvari&);
00305 #endif
00306 
00307   friend void AD_Const(const IndepADvar&);
00308 
00309   inline operator ADvari&() { return *cv; }
00310   inline operator ADvari&() const { return *(ADvari*)cv; }
00311 
00312   inline double val() const { return cv->Val; }
00313   inline double adj() const { return cv->aval; }
00314   static inline void Gradcomp() { ADcontext::Gradcomp(); }
00315   static inline void aval_reset() { ConstADvari::aval_reset(); }
00316   static inline void Weighted_Gradcomp(int n, ADvar **v, double *w)
00317         { ADcontext::Weighted_Gradcomp(n, v, w); }
00318 
00319 
00320 #define Ai const ADvari&
00321 #define AI const IndepADvar&
00322 #define T2(r,f) \
00323  r f(AI,AI);\
00324  r f(Ai,AI);\
00325  r f(AI,Ai);\
00326  r f(double,AI);\
00327  r f(AI,double);
00328 #define T1(f) friend ADvari& f(AI);
00329 
00330 #define F friend ADvari&
00331 T2(F, operator+)
00332 T2(F, operator-)
00333 T2(F, operator*)
00334 T2(F, operator/)
00335 T2(F, atan2)
00336 T2(F, max)
00337 T2(F, min)
00338 T2(F, pow)
00339 #undef F
00340 #define F friend int
00341 T2(F, operator<)
00342 T2(F, operator<=)
00343 T2(F, operator==)
00344 T2(F, operator!=)
00345 T2(F, operator>=)
00346 T2(F, operator>)
00347 
00348 T1(operator+)
00349 T1(operator-)
00350 T1(abs)
00351 T1(acos)
00352 T1(acosh)
00353 T1(asin)
00354 T1(asinh)
00355 T1(atan)
00356 T1(atanh)
00357 T1(cos)
00358 T1(cosh)
00359 T1(exp)
00360 T1(log)
00361 T1(log10)
00362 T1(sin)
00363 T1(sinh)
00364 T1(sqrt)
00365 T1(tan)
00366 T1(tanh)
00367 T1(fabs)
00368 T1(copy)
00369 
00370 #undef T1
00371 #undef T2
00372 #undef F
00373 #undef Ai
00374 #undef AI
00375 
00376   };
00377 
00378  class
00379 ADvar: public IndepADvar {    // an "active" variable
00380   void ADvar_ctr(double d);
00381  public:
00382   inline ADvar() { /* cv = 0; */ }
00383   inline ADvar(double d) { ADvar_ctr(d); }
00384   inline ADvar(int i) { ADvar_ctr((double)i); }
00385   inline ADvar(long i)  { ADvar_ctr((double)i); }
00386   inline ~ADvar() {}
00387 #ifdef RAD_AUTO_AD_Const
00388   friend class ADvar1;
00389   inline ADvar(const IndepADvar &x) { cv = x.cv ? new ADvar1(this, x) : 0; }
00390   inline ADvar(ADvari &x) { cv = &x; x.padv = this; }
00391   inline ADvar& operator=(const IndepADvar &x) {
00392     if (cv)
00393       cv->padv = 0;
00394     cv = new ADvar1(this,x);
00395     return *this;
00396     }
00397   inline ADvar& operator=(const ADvari &x) {
00398     if (cv)
00399       cv->padv = 0;
00400     cv = new ADvar1(this, x);
00401     return *this;
00402     }
00403 #else
00404   friend ADvar& ADvar_operatoreq(ADvar*, const ADvari&);
00405 #ifdef RAD_EQ_ALIAS
00406   /* allow aliasing v and w after "v = w;" */
00407   inline ADvar(const IndepADvar &x) { cv = x.cv; }
00408   inline ADvar(const ADvari &x) { cv = (ADvari*)&x; }
00409   inline ADvar& operator=(const ADvari &x) { cv = (ADvari*)&x;   return *this; }
00410   inline ADvar& operator=(const IndepADvar &x)  { cv = (ADvari*)x.cv; return *this; }
00411 #else
00412   ADvar(const IndepADvar &x) { cv = x.cv ? new ADvar1(x.cv->Val, &CADcontext::One, x.cv) : 0; }
00413   ADvar(const ADvari &x) { cv = new ADvar1(x.Val, &CADcontext::One, &x); }
00414   inline ADvar& operator=(const ADvari &x) { return ADvar_operatoreq(this,x); }
00415   inline ADvar& operator=(const IndepADvar &x)    { return ADvar_operatoreq(this,*x.cv); }
00416 #endif /* RAD_EQ_ALIAS */
00417 #endif /* RAD_AUTO_AD_Const */
00418   ADvar& operator=(double);
00419   ADvar& operator+=(const ADvari&);
00420   ADvar& operator+=(double);
00421   ADvar& operator-=(const ADvari&);
00422   ADvar& operator-=(double);
00423   ADvar& operator*=(const ADvari&);
00424   ADvar& operator*=(double);
00425   ADvar& operator/=(const ADvari&);
00426   ADvar& operator/=(double);
00427   inline static bool get_fpval_implies_const(void)
00428     { return ConstADvari::cadc.fpval_implies_const; }
00429   inline static void set_fpval_implies_const(bool newval)
00430     { ConstADvari::cadc.fpval_implies_const = newval; }
00431   inline static bool setget_fpval_implies_const(bool newval) {
00432     bool oldval = ConstADvari::cadc.fpval_implies_const;
00433     ConstADvari::cadc.fpval_implies_const = newval;
00434     return oldval;
00435     }
00436   static inline void Gradcomp() { ADcontext::Gradcomp(); }
00437   static inline void aval_reset() { ConstADvari::aval_reset(); }
00438   static inline void Weighted_Gradcomp(int n, ADvar **v, double *w)
00439         { ADcontext::Weighted_Gradcomp(n, v, w); }
00440   };
00441 
00442  inline void AD_Const(const IndepADvar&v) { IndepADvar::AD_Const(v); }
00443 
00444  class
00445 ConstADvar: public ADvar {
00446  private: // disable op=
00447   ConstADvar& operator+=(const ADvari&);
00448   ConstADvar& operator+=(double);
00449   ConstADvar& operator-=(const ADvari&);
00450   ConstADvar& operator-=(double);
00451   ConstADvar& operator*=(const ADvari&);
00452   ConstADvar& operator*=(double);
00453   ConstADvar& operator/=(const ADvari&);
00454   ConstADvar& operator/=(double);
00455   void ConstADvar_ctr(double);
00456  public:
00457   inline ConstADvar(double d) { ConstADvar_ctr(d); }
00458   inline ConstADvar(int i)  { ConstADvar_ctr((double)i); }
00459   inline ConstADvar(long i) { ConstADvar_ctr((double)i); }
00460   ConstADvar(const ADvari &x);
00461 #ifdef RAD_AUTO_AD_Const
00462   ConstADvar(const IndepADvar &x) { cv = new ADvar1(this,x); }
00463 #endif
00464   inline ~ConstADvar(){}
00465 #ifdef RAD_NO_CONST_UPDATE
00466  private:
00467 #endif
00468   ConstADvar();
00469   inline ConstADvar& operator=(double d) { cv->Val = d; return *this; }
00470   inline ConstADvar& operator=(const IndepADvar& d) { cv->Val = d.val(); return *this; }
00471  };
00472 
00473  class
00474 ADvar1s: public ADvar1 { // unary ops with partial "a"
00475  public:
00476   double a;
00477   ADvar1s(double val1, double a1, const ADvari *c1): ADvar1(val1,&a,c1), a(a1) {}
00478   };
00479 
00480  class
00481 ADvar2: public ADvari { // basic binary ops
00482  public:
00483   Derp dL, dR;
00484   ADvar2(double val1): ADvari(val1) {}
00485   ADvar2(double val1, const ADvari *Lcv, const double *Lc, const ADvari *Rcv,
00486       const double *Rc): ADvari(val1) {
00487     dR.next = Derp::LastDerp;
00488     dL.next = &dR;
00489     Derp::LastDerp = &dL;
00490     dL.a = Lc;
00491     dL.c = Lcv;
00492     dR.a = Rc;
00493     dR.c = Rcv;
00494     dL.b = dR.b = this;
00495     }
00496   };
00497 
00498  class
00499 ADvar2q: public ADvar2 { // binary ops with partials "a", "b"
00500  public:
00501   double a, b;
00502   ADvar2q(double val1, double Lp, double Rp, const ADvari *Lcv, const ADvari *Rcv):
00503       ADvar2(val1), a(Lp), b(Rp) {
00504     dR.next = Derp::LastDerp;
00505     dL.next = &dR;
00506     Derp::LastDerp = &dL;
00507     dL.a = &a;
00508     dL.c = Lcv;
00509     dR.a = &b;
00510     dR.c = Rcv;
00511     dL.b = dR.b = this;
00512     }
00513   };
00514 
00515  class
00516 ADvarn: public ADvari { // n-ary ops with partials "a"
00517  public:
00518   int n;
00519   double *a;
00520   Derp *Da;
00521   ADvarn(double val1, int n1, const ADvar *x, const double *g);
00522   };
00523 
00524 inline ADvari &operator+(ADvari &T) { return T; }
00525 inline ADvari &operator+(const ADvari &T) { return (ADvari&) T; }
00526 
00527 inline int operator<(const ADvari &L, const ADvari &R) { return L.Val < R.Val; }
00528 inline int operator<(const ADvari &L, double R) { return L.Val < R; }
00529 inline int operator<(double L, const ADvari &R) { return L < R.Val; }
00530 
00531 inline int operator<=(const ADvari &L, const ADvari &R) { return L.Val <= R.Val; }
00532 inline int operator<=(const ADvari &L, double R) { return L.Val <= R; }
00533 inline int operator<=(double L, const ADvari &R) { return L <= R.Val; }
00534 
00535 inline int operator==(const ADvari &L, const ADvari &R) { return L.Val == R.Val; }
00536 inline int operator==(const ADvari &L, double R) { return L.Val == R; }
00537 inline int operator==(double L, const ADvari &R) { return L == R.Val; }
00538 
00539 inline int operator!=(const ADvari &L, const ADvari &R) { return L.Val != R.Val; }
00540 inline int operator!=(const ADvari &L, double R) { return L.Val != R; }
00541 inline int operator!=(double L, const ADvari &R) { return L != R.Val; }
00542 
00543 inline int operator>=(const ADvari &L, const ADvari &R) { return L.Val >= R.Val; }
00544 inline int operator>=(const ADvari &L, double R) { return L.Val >= R; }
00545 inline int operator>=(double L, const ADvari &R) { return L >= R.Val; }
00546 
00547 inline int operator>(const ADvari &L, const ADvari &R) { return L.Val > R.Val; }
00548 inline int operator>(const ADvari &L, double R) { return L.Val > R; }
00549 inline int operator>(double L, const ADvari &R) { return L > R.Val; }
00550 
00551 inline ADvari& copy(const IndepADvar &x)
00552 { return *(new ADvar1(x.cv->Val, &CADcontext::One, x.cv)); }
00553 
00554 inline ADvari& copy(const ADvari &x)
00555 { return *(new ADvar1(x.Val, &CADcontext::One, &x)); }
00556 
00557 inline ADvari& abs(const ADvari &x)
00558 { return fabs(x); }
00559 
00560 #define A (ADvari*)
00561 #define T inline
00562 #define F ADvari&
00563 #define Ai const ADvari&
00564 #define AI const IndepADvar&
00565 #define D double
00566 #define T2(r,f) \
00567  T r f(Ai L, AI R) { return f(L, *A R.cv); }\
00568  T r f(AI L, Ai R) { return f(*A L.cv, R); }\
00569  T r f(AI L, AI R) { return f(*A L.cv, *A R.cv); }\
00570  T r f(AI L, D R) { return f(*A L.cv, R); }\
00571  T r f(D L, AI R) { return f(L, *A R.cv); }
00572 
00573 T2(F, operator+)
00574 T2(F, operator-)
00575 T2(F, operator*)
00576 T2(F, operator/)
00577 T2(F, atan2)
00578 T2(F, pow)
00579 T2(F, max)
00580 T2(F, min)
00581 T2(int, operator<)
00582 T2(int, operator<=)
00583 T2(int, operator==)
00584 T2(int, operator!=)
00585 T2(int, operator>=)
00586 T2(int, operator>)
00587 
00588 #undef T2
00589 #undef D
00590 
00591 #define T1(f)\
00592  T F f(AI x) { return f(*A x.cv); }
00593 
00594 T1(operator+)
00595 T1(operator-)
00596 T1(abs)
00597 T1(acos)
00598 T1(acosh)
00599 T1(asin)
00600 T1(asinh)
00601 T1(atan)
00602 T1(atanh)
00603 T1(cos)
00604 T1(cosh)
00605 T1(exp)
00606 T1(log)
00607 T1(log10)
00608 T1(sin)
00609 T1(sinh)
00610 T1(sqrt)
00611 T1(tan)
00612 T1(tanh)
00613 T1(fabs)
00614 
00615 #undef T1
00616 #undef AI
00617 #undef Ai
00618 #undef F
00619 #undef T
00620 #undef A
00621 
00622 } // namespace Radnt
00623 } // namespace Sacado
00624 #endif /* SACADO_RAD_H */

Generated on Wed May 12 21:59:05 2010 for Sacado Package Browser (Single Doxygen Collection) by  doxygen 1.4.7