Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042 #ifndef SUNDANCE_ADREAL_H
00043 #define SUNDANCE_ADREAL_H
00044
00045
00046 #include "SundanceDefs.hpp"
00047 #include "SundancePoint.hpp"
00048 #include "SundanceStdMathFunctors.hpp"
00049 #include "Teuchos_RefCountPtr.hpp"
00050
00051 namespace Sundance
00052 {
00053
00054
00055
00056
00057 class ADReal
00058 {
00059 public:
00060
00061 ADReal() : value_(0), gradient_(0.0, 0.0, 0.0){;}
00062
00063
00064 ADReal(double value, int direction, int spatialDimension)
00065 : value_(value), gradient_()
00066 {
00067 if (spatialDimension==1) gradient_ = Point(0.0);
00068 if (spatialDimension==2) gradient_ = Point(0.0, 0.0);
00069 if (spatialDimension==3) gradient_ = Point(0.0, 0.0, 0.0);
00070 gradient_[direction] = 1.0;
00071 }
00072
00073 ADReal(double value, int spatialDimension)
00074 : value_(value), gradient_()
00075 {
00076 if (spatialDimension==1) gradient_ = Point(0.0);
00077 if (spatialDimension==2) gradient_ = Point(0.0, 0.0);
00078 if (spatialDimension==3) gradient_ = Point(0.0, 0.0, 0.0);
00079 }
00080
00081
00082 ADReal(double value, const Point& gradient)
00083 : value_(value), gradient_(gradient) {;}
00084
00085
00086 ADReal operator-() const ;
00087
00088 ADReal& operator+=(const ADReal& other) ;
00089
00090 ADReal& operator-=(const ADReal& other) ;
00091
00092 ADReal& operator*=(const ADReal& other) ;
00093
00094 ADReal& operator/=(const ADReal& other) ;
00095
00096
00097 ADReal& operator+=(const double& scalar) ;
00098
00099 ADReal& operator-=(const double& scalar) ;
00100
00101 ADReal& operator*=(const double& scalar) ;
00102
00103 ADReal& operator/=(const double& scalar) ;
00104
00105
00106 ADReal operator+(const ADReal& other) const ;
00107
00108 ADReal operator-(const ADReal& other) const ;
00109
00110 ADReal operator*(const ADReal& other) const ;
00111
00112 ADReal operator/(const ADReal& other) const ;
00113
00114
00115 ADReal operator+(const double& scalar) const ;
00116
00117 ADReal operator-(const double& scalar) const ;
00118
00119 ADReal operator*(const double& scalar) const ;
00120
00121 ADReal operator/(const double& scalar) const ;
00122
00123
00124 const double& value() const {return value_;}
00125
00126 const Point& gradient() const {return gradient_;}
00127
00128 void reciprocate() ;
00129
00130 static double& totalFlops() {static double rtn = 0; return rtn;}
00131
00132
00133 static void addFlops(const double& flops) {totalFlops() += flops;}
00134 private:
00135 double value_;
00136 Point gradient_;
00137 };
00138
00139
00140
00141
00142 ADReal operator+(const double& scalar,
00143 const ADReal& a);
00144
00145 ADReal operator-(const double& scalar,
00146 const ADReal& a);
00147
00148 ADReal operator*(const double& scalar,
00149 const ADReal& a);
00150
00151 ADReal operator/(const double& scalar,
00152 const ADReal& a);
00153
00154
00155 inline ADReal pow(const ADReal& x, const double& y)
00156 {
00157 Teuchos::RCP<UnaryFunctor> func = Teuchos::rcp(new PowerFunctor(y));
00158 double f;
00159 double df;
00160 double xVal = x.value();
00161 func->eval1(&xVal, 1, &f, &df);
00162 return ADReal(f, df*x.gradient());
00163 }
00164
00165 #define SUNDANCE_AD_FUNCTOR(opName, functorName) \
00166 inline ADReal opName(const ADReal& x)\
00167 {\
00168 Teuchos::RCP<UnaryFunctor> func = Teuchos::rcp(new functorName());\
00169 double f;\
00170 double df;\
00171 double xVal = x.value();\
00172 func->eval1(&xVal, 1, &f, &df);\
00173 return ADReal(f, df*x.gradient());}
00174
00175
00176 SUNDANCE_AD_FUNCTOR(exp, StdExp)
00177
00178 SUNDANCE_AD_FUNCTOR(log, StdLog)
00179
00180 SUNDANCE_AD_FUNCTOR(sqrt, StdSqrt)
00181
00182 SUNDANCE_AD_FUNCTOR(sin, StdSin)
00183
00184 SUNDANCE_AD_FUNCTOR(cos, StdCos)
00185
00186
00187
00188 }
00189
00190
00191 namespace std
00192 {
00193 inline ostream& operator<<(std::ostream& os, const Sundance::ADReal& x)
00194 {
00195 std::cerr << "AD[" << x.value() << ", grad="<< x.gradient() << "]";
00196 return os;
00197 }
00198 }
00199
00200
00201
00202 #endif
00203
00204
00205