Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages   Examples  

ADReal.h

Go to the documentation of this file.
00001 #ifndef ADREAL_H
00002 #define ADREAL_H
00003 
00004 #include "DenseSerialVector.h"
00005 
00006 
00007 namespace Sundance
00008 {
00009   using namespace TSF;
00010 
00011   /**\ingroup Numerical
00012    * First-order automatic differentiation of a multivariable function.
00013    * @author Kevin Long
00014    */
00015 
00016   class ADReal
00017     {
00018     public:
00019       /** Create an ADReal with value and gradient = 0 */
00020       ADReal() : value_(0), gradient_(0){;}
00021       /** Create an ADReal object that varies linearly with
00022        * one coordinate direction  */
00023       ADReal(double value, int direction, int spatialDimension)
00024         : value_(value), gradient_(spatialDimension, 0.0)
00025         {
00026           gradient_[direction] = 1.0;
00027         }
00028       /** Create a constant-valued ADReal object in a multidimensional space */
00029       ADReal(double value, int spatialDimension)
00030         : value_(value), gradient_(spatialDimension, 0.0)
00031         {;}
00032 
00033       /** unary minus */
00034       ADReal operator-() const ;
00035       /** reflexive addition */
00036       ADReal& operator+=(const ADReal& other) ;
00037       /** reflexive subtraction */
00038       ADReal& operator-=(const ADReal& other) ;
00039       /** reflexive multiplication */
00040       ADReal& operator*=(const ADReal& other) ;
00041       /** reflexive division */
00042       ADReal& operator/=(const ADReal& other) ;
00043 
00044       /** reflexive scalar addition */
00045       ADReal& operator+=(const double& scalar) ;
00046       /** reflexive scalar subtraction */
00047       ADReal& operator-=(const double& scalar) ;
00048       /** reflexive scalar multiplication */
00049       ADReal& operator*=(const double& scalar) ;
00050       /** reflexive scalar division */
00051       ADReal& operator/=(const double& scalar) ;
00052 
00053       /** addition */
00054       ADReal operator+(const ADReal& other) const ;
00055       /** subtraction */
00056       ADReal operator-(const ADReal& other) const ;
00057       /** multiplication */
00058       ADReal operator*(const ADReal& other) const ;
00059       /** division */
00060       ADReal operator/(const ADReal& other) const ;
00061 
00062       /** scalar addition */
00063       ADReal operator+(const double& scalar) const ;
00064       /** scalar subtraction */
00065       ADReal operator-(const double& scalar) const ;
00066       /** scalar multiplication */
00067       ADReal operator*(const double& scalar) const ;
00068       /** scalar division */
00069       ADReal operator/(const double& scalar) const ;
00070 
00071       /** get the value */
00072       const double& value() const {return value_;}
00073       /** get the gradient */
00074       const DenseSerialVector& gradient() const {return gradient_;}
00075 
00076       void reciprocate() ;
00077 
00078       static bool unitTest() ;
00079     private:
00080       double value_;
00081       DenseSerialVector gradient_;
00082     };
00083 
00084   ADReal operator+(const double& scalar,
00085                    const ADReal& a);
00086   ADReal operator-(const double& scalar,
00087                    const ADReal& a);
00088   ADReal operator*(const double& scalar,
00089                    const ADReal& a);
00090   ADReal operator/(const double& scalar,
00091                    const ADReal& a);
00092 
00093 
00094 }
00095 
00096 
00097 
00098 #endif
00099 
00100 
00101 

Contact:
Kevin Long (krlong@ca.sandia.gov)


Documentation generated by