ThreeB 1.1
|
00001 #ifndef INC_STORM_H 00002 #define INC_STORM_H 00003 00004 #include <TooN/TooN.h> 00005 #include <cvd/image.h> 00006 #include <utility> 00007 #include <tr1/tuple> 00008 #include "utility.h" 00009 00010 00011 /**See spot_shape() 00012 @param x \f$\Vec{x}\f$ 00013 @param phi \f$\Vec{\phi}\f$ 00014 @return \f$s(\Vec{x}, \Vec{\phi}) \f$ 00015 @ingroup gStorm 00016 */ 00017 template<class B> double spot_shape_s(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00018 { 00019 return -norm_sq(x - phi.template slice<2,2>()) / (2*phi[1]*phi[1]); 00020 } 00021 00022 /** Compute the spot shape and its derivative with respect to posision. See also spot_shape() 00023 @param x \f$\Vec{x}\f$ 00024 @param phi \f$\Vec{\phi}\f$ 00025 @ingroup gStorm 00026 */ 00027 template<class B> std::pair<double, TooN::Vector<4> > spot_shape_diff_position(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00028 { 00029 using namespace TooN; 00030 00031 double s = spot_shape_s(x, phi); 00032 double r_2_pi = sqrt(2*M_PI); 00033 00034 double prob = exp(s) * phi[0]/(phi[1]*r_2_pi); 00035 00036 Vector<4> deriv = (exp(s) / (phi[1]*r_2_pi)) * 00037 makeVector(1, 00038 -phi[0] * (1 + 2*s)/phi[1], 00039 (x[0] - phi[2])*(phi[0]/sq(phi[1])), 00040 (x[1] - phi[3])*(phi[0]/sq(phi[1]))); 00041 return std::make_pair(prob, deriv); 00042 } 00043 00044 /** Compute the spot shape and its Hessian with respect to posision. See also spot_shape() 00045 @param x \f$\Vec{x}\f$ 00046 @param phi \f$\Vec{\phi}\f$ 00047 @ingroup gStorm 00048 */ 00049 template<class B> std::tr1::tuple<double, TooN::Vector<4>, TooN::Matrix<4> > spot_shape_hess_position(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00050 { 00051 using namespace TooN; 00052 using namespace std::tr1; 00053 00054 double s = spot_shape_s(x, phi); 00055 double r_2_pi = sqrt(2*M_PI); 00056 00057 double es = exp(s); 00058 00059 double prob = es * phi[0]/(phi[1]*r_2_pi); 00060 00061 Vector<4> deriv = (es / (phi[1]*r_2_pi)) * 00062 makeVector(1, 00063 -phi[0] * (1 + 2*s)/phi[1], 00064 (x[0] - phi[2])*(phi[0]/sq(phi[1])), 00065 (x[1] - phi[3])*(phi[0]/sq(phi[1]))); 00066 00067 Matrix<4> hess; 00068 hess[0][0] = 0; 00069 00070 hess[0][1] = -es*(1+2*s) / (phi[1] * phi[1] * r_2_pi); 00071 hess[1][0] = hess[0][1]; 00072 00073 hess[0][2] = es * (x[0] - phi[2]) / (pow(phi[1], 3)*r_2_pi); 00074 hess[2][0] = es * (x[0] - phi[2]) / (pow(phi[1], 3)*r_2_pi); 00075 00076 hess[0][3] = es * (x[1] - phi[3]) / (pow(phi[1], 3)*r_2_pi); 00077 hess[3][0] = es * (x[1] - phi[3]) / (pow(phi[1], 3)*r_2_pi); 00078 00079 hess[1][1] = 2*phi[0]*es*(1 + 5*s + 2*s*s) / ( pow(phi[1], 3) * r_2_pi); 00080 00081 hess[1][2] = -phi[0] * es * (3 + 2*s) * (x[0] - phi[2]) / (pow(phi[1], 4) * r_2_pi); 00082 hess[1][3] = -phi[0] * es * (3 + 2*s) * (x[1] - phi[3]) / (pow(phi[1], 4) * r_2_pi); 00083 00084 hess[2][1] = hess[1][2]; 00085 hess[3][1] = hess[1][3]; 00086 00087 hess[2][2] = phi[0] * es * (sq(x[0] - phi[2]) - sq(phi[1])) / (r_2_pi * pow(phi[1], 5)); 00088 hess[3][3] = phi[0] * es * (sq(x[1] - phi[3]) - sq(phi[1])) / (r_2_pi * pow(phi[1], 5)); 00089 00090 hess[2][3] = phi[0] * es * (x[0] - phi[2])*(x[1] - phi[3]) / (r_2_pi * pow(phi[1], 5)); 00091 hess[3][2] = hess[2][3]; 00092 00093 00094 return make_tuple(prob, deriv, hess); 00095 } 00096 00097 /**Value of the spot, given the parameters and input location. 00098 The spot is described by the following formula: 00099 \f[ 00100 \mu(\Vec{x}, \Vec{\phi}) = \frac{\phi_1}{\phi_2\sqrt(2\pi)} e^s, 00101 \f] 00102 where 00103 \f[ 00104 s = -\frac{(x_1 - \phi_3)^2 + (x_2 - \phi_4)^2}{2\phi_2^2}. 00105 \f] 00106 This describes a generic blobby spot function of a variable size. The light output 00107 can be tuned by varying \f$\phi_1\f$, and the level of blur can be changed independently 00108 by varying \f$\phi_2\f$. The derivative is: 00109 \f{eqnarray}{ 00110 \frac{\partial \mu}{\partial \phi_1} &=& \frac{1}{\phi_2\sqrt{2\pi}}e^s\\ 00111 \frac{\partial \mu}{\partial \phi_2} &=& -\frac{\phi_1}{\phi_2^2\sqrt{2\pi}}e^s(1 + 2s) 00112 \f} 00113 And the hessian is: 00114 \f{eqnarray}{ 00115 \frac{\partial^2 \mu}{\partial \phi_1^2} &=& 0\\ 00116 \frac{\partial^2 \mu}{\partial\phi_1 \partial \phi_2} &=& -\frac{1}{\phi_2^2\sqrt{2\pi}}e^s(1 + 2s)\\ 00117 \frac{\partial^2 \mu}{\partial \phi_2^2} &=& \frac{2\phi_1}{\phi_2^3\sqrt{2\pi}}e^s(1 + 5s + 2s^2) 00118 \f} 00119 @param x \f$\Vec{x}\f$ 00120 @param phi \f$\Vec{\phi}\f$ 00121 @return \f$\mu(\Vec{x}, \Vec{\phi}) \f$ 00122 @ingroup gStorm 00123 */ 00124 template<class B> std::tr1::tuple<double, TooN::Vector<2>, TooN::Matrix<2> > spot_shape_hess(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00125 { 00126 double s = spot_shape_s(x, phi); 00127 double r_2_pi = sqrt(2*M_PI); 00128 00129 double prob = exp(s) * phi[0]/(phi[1]*r_2_pi); 00130 TooN::Vector<2> deriv = (exp(s) / (phi[1]*r_2_pi)) * TooN::makeVector(1, -phi[0] * (1 + 2*s)/phi[1]); 00131 TooN::Matrix<2> hess; 00132 00133 hess[0][0] = 0; 00134 hess[0][1] = -exp(s)*(1+2*s) / (phi[1] * phi[1] * r_2_pi); 00135 hess[1][0] = hess[0][1]; 00136 hess[1][1] = 2*phi[0]*exp(s)*(1 + 5*s + 2*s*s) / ( pow(phi[1], 3) * r_2_pi); 00137 00138 return std::tr1::make_tuple(prob, deriv, hess); 00139 } 00140 /** see spot_shape_hess() 00141 @param x \f$\Vec{x}\f$ 00142 @param phi \f$\Vec{\phi}\f$ 00143 @return \f$\mu(\Vec{x}, \Vec{\phi}) \f$ 00144 @ingroup gStorm 00145 */ 00146 template<class B> std::pair<double, TooN::Vector<2> > spot_shape_diff(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00147 { 00148 double s = spot_shape_s(x, phi); 00149 double r_2_pi = sqrt(2*M_PI); 00150 00151 double prob = exp(s) * phi[0]/(phi[1]*r_2_pi); 00152 TooN::Vector<2> deriv = (exp(s) / (phi[1]*r_2_pi)) * TooN::makeVector(1, -phi[0] * (1 + 2*s)/phi[1]); 00153 return std::make_pair(prob, deriv); 00154 } 00155 00156 /** see spot_shape_hess() 00157 @param x \f$\Vec{x}\f$ 00158 @param phi \f$\Vec{\phi}\f$ 00159 @return \f$\mu(\Vec{x}, \Vec{\phi}) \f$ 00160 @ingroup gStorm 00161 */ 00162 template<class B> double spot_shape(const TooN::Vector<2>& x, const TooN::Vector<4, double, B>& phi) 00163 { 00164 double s = spot_shape_s(x, phi); 00165 double r_2_pi = sqrt(2*M_PI); 00166 00167 // FIXME FIXME FIXME and don't forget to fix the HESSIAN AND DERIVATIVE 00168 // Should be: 1/(2 pi s^2) for two dimensions 00169 // vvvvvvvvvvvvv http://lol.i.trollyou.com/ 00170 double prob = exp(s) * phi[0]/(phi[1]*r_2_pi); 00171 00172 00173 return prob; 00174 } 00175 00176 /**Find the log probability of an image patch, 00177 assuming zero mean and the given variance, and no spot present. 00178 See also log_probability_spot() 00179 @param im Image 00180 @param variance variance 00181 @returns The log probability 00182 */ 00183 inline double log_probability_no_spot(const CVD::SubImage<float>& im, double variance) 00184 { 00185 double logprob_part=0; 00186 for(int y=0; y < im.size().y; y++) 00187 for(int x=0; x < im.size().x; x++) 00188 logprob_part -= im[y][x] * im[y][x]; 00189 return logprob_part/(2*variance) - im.size().area() * log(2*M_PI*variance)/2; 00190 00191 } 00192 00193 /**Find the log probability of an image patch, assuming zero base-line mean and 00194 the given variance. This function makes use of the spot shape. It is assumed that the 00195 centre pixel of the image is at 0,0. Since the noise is Gaussian: 00196 \f{eqnarray}{ 00197 P(\text{image}) &=& \prod_{\Vec{x} \in \text{pixels}} \frac{1}{\sqrt{2\pi\sigma^2}}e^{-\frac{(I(\Vec{x}) - \mu(\Vec{x}, \Vec{\phi}))^2}{2\sigma^2}} \\ 00198 \ln P(\text{image}) &=& \sum_{\Vec{x} \in \text{pixels}} -\frac{(I(\Vec{x}) - \mu(\Vec{x}, \Vec{\phi}))^2}{2\sigma^2} - \frac{N}{2} \ln {2 \pi \sigma^2}, 00199 \f} 00200 where \e I is the image, and \e N is the number of pixels. See also ::log_probability_no_spot and \f$\mu\f$ (::spot_shape). 00201 The derivatives are: 00202 \f{eqnarray}{ 00203 \frac{\partial \ln P(I)}{\partial \phi_0} &=& \frac{1}{\sigma^2} \sum_{\Vec{x}}(I_{\Vec{x}} - \mu(\Vec{x},\Vec{\phi})) 00204 \frac{\partial}{\partial \phi_0}\mu(\Vec{x}, \Vec{\phi})\\ 00205 \frac{\partial^2 \ln P(I)}{\partial \phi_0 \partial \phi_1} &=& 00206 \frac{1}{\sigma^2} \sum_{\Vec{x}}(I_{\Vec{x}} - \mu(\Vec{x},\Vec{\phi})) 00207 \frac{\partial^2}{\partial \phi_0 \partial \phi_1}\mu(\Vec{x}, \Vec{\phi}) - 00208 \frac{\partial}{\partial \phi_0}\mu(\Vec{x},\Vec{\phi}) 00209 \frac{\partial}{\partial \phi_1}\mu(\Vec{x},\Vec{\phi}) 00210 \f} 00211 @ingroup gStorm 00212 @param im Image 00213 @param variance \f$\sigma^2\f$ 00214 @param spot_parameters \f$\Vec{\phi}\f$ 00215 @returns The log probability 00216 */ 00217 template<class Base> std::tr1::tuple<double, TooN::Vector<2>, TooN::Matrix<2> > log_probability_spot_hess(const CVD::SubImage<float>& im, double variance, const TooN::Vector<4, double, Base>& spot_parameters) 00218 { 00219 using namespace TooN; 00220 using namespace std::tr1; 00221 00222 //-1 because if the image is 3x3, ie 0,1,2 then 1,1 is the centre. 00223 //If it is 2x2, ie 0,1 then .5,.5 is the centre 00224 Vector<2> centre = makeVector((im.size().x-1) / 2.0, (im.size().y-1) / 2.0); 00225 00226 double logprob_part=0; 00227 Vector<2> diff = Zeros; 00228 Matrix<2> hess = Zeros; 00229 for(int y=0; y < im.size().y; y++) 00230 for(int x=0; x < im.size().x; x++) 00231 { 00232 Vector<2> d = TooN::makeVector(x, y) - centre; 00233 00234 double mu; 00235 Vector<2> diff_mu; 00236 Matrix<2> hess_mu; 00237 tie(mu, diff_mu, hess_mu) = spot_shape_hess(d, spot_parameters); 00238 00239 double e = im[y][x] - mu; 00240 00241 logprob_part += -sq(e); 00242 diff += diff_mu * e; 00243 hess += e * hess_mu - diff_mu.as_col() * diff_mu.as_row(); 00244 } 00245 return make_tuple( logprob_part / (2*variance) - im.size().area() * log(2*M_PI*variance)/2, 00246 diff / variance, 00247 hess / variance); 00248 } 00249 00250 /** See log_probability_spot_hess 00251 @ingroup gStorm 00252 @param im Image 00253 @param variance \f$\sigma^2\f$ 00254 @param spot_parameters \f$\Vec{\phi}\f$ 00255 @returns The log probability 00256 */ 00257 template<class Base> std::pair<double, TooN::Vector<2> > log_probability_spot_diff(const CVD::SubImage<float>& im, double variance, const TooN::Vector<4, double, Base>& spot_parameters) 00258 { 00259 using namespace TooN; 00260 using namespace std::tr1; 00261 using namespace std; 00262 //-1 because if the image is 3x3, ie 0,1,2 then 1,1 is the centre. 00263 //If it is 2x2, ie 0,1 then .5,.5 is the centre 00264 Vector<2> centre = makeVector((im.size().x-1) / 2.0, (im.size().y-1) / 2.0); 00265 00266 double logprob_part=0; 00267 Vector<2> diff = Zeros; 00268 for(int y=0; y < im.size().y; y++) 00269 for(int x=0; x < im.size().x; x++) 00270 { 00271 Vector<2> d = makeVector(x, y) - centre; 00272 00273 double mu; 00274 Vector<2> diff_mu; 00275 tie(mu, diff_mu) = spot_shape_diff(d, spot_parameters); 00276 00277 double e = im[y][x] - mu; 00278 00279 logprob_part += -sq(e); 00280 diff += diff_mu * e; 00281 } 00282 return make_pair(logprob_part / (2*variance) - im.size().area() * log(2*M_PI*variance)/2, diff / variance); 00283 } 00284 00285 /** See log_probability_spot_hess 00286 @ingroup gStorm 00287 @param im Image 00288 @param variance \f$\sigma^2\f$ 00289 @param spot_parameters \f$\Vec{\phi}\f$ 00290 @returns The log probability 00291 */ 00292 template<class Base> double log_probability_spot(const CVD::SubImage<float>& im, double variance, const TooN::Vector<4, double, Base>& spot_parameters) 00293 { 00294 //-1 because if the image is 3x3, ie 0,1,2 then 1,1 is the centre. 00295 //If it is 2x2, ie 0,1 then .5,.5 is the centre 00296 TooN::Vector<2> centre = TooN::makeVector((im.size().x-1) / 2.0, (im.size().y-1) / 2.0); 00297 00298 double logprob_part=0; 00299 for(int y=0; y < im.size().y; y++) 00300 for(int x=0; x < im.size().x; x++) 00301 { 00302 TooN::Vector<2> d = TooN::makeVector(x, y) - centre; 00303 00304 double mu = spot_shape(d, spot_parameters); 00305 00306 double e = im[y][x] - mu; 00307 00308 logprob_part += -sq(e); 00309 } 00310 return logprob_part / (2*variance) - im.size().area() * log(2*M_PI*variance)/2; 00311 } 00312 00313 /**Compute the standard deviation of a log-normal distribution. 00314 00315 See log_normal(). 00316 \f{equation} 00317 \mathrm{Var}[P(x)] = (e^(\sigma^2)-1)e^(2*\mu+\sigma^2) 00318 \f} 00319 @param sigma \f$ \sigma\f$ 00320 @param mu \f$ \mu\f$ 00321 @returns The standard deviation 00322 @ingroup gStorm 00323 */ 00324 inline double log_normal_std(double mu, double sigma) 00325 { 00326 return sqrt((exp(sq(sigma)) - 1) * exp(2*mu + sq(sigma))); 00327 } 00328 00329 /**Compute the mode of a log-normal distribution. 00330 00331 See log_normal(). 00332 \f{equation} 00333 \mathrm{Mode}[P(x)] = e^(\mu-\sigma^2) 00334 \f} 00335 @param sigma \f$ \sigma\f$ 00336 @param mu \f$ \mu\f$ 00337 @returns The mode 00338 @ingroup gStorm 00339 */ 00340 inline double log_normal_mode(double mu, double sigma) 00341 { 00342 return exp(mu - sigma * sigma); 00343 } 00344 00345 /**Log-normal distribution. This is given by: 00346 \f{eqnarray}{ 00347 P(x) &=& \frac{1}{x\sigma\sqrt{2\pi}} e^{-\frac{(\ln x - \mu)^2}{s\sigma^2}}\\ 00348 \ln P(x) &=& -\frac{(\ln x - \mu)^2}{s\sigma^2} - \ln x - \ln\sigma\sqrt{2\pi}. 00349 \f} 00350 @param x \e x 00351 @param mu \f$\mu\f$ 00352 @param sigma \f$\sigma\f$ 00353 @ingroup gStorm 00354 */ 00355 inline double log_log_normal(double x, double mu, double sigma) 00356 { 00357 return -sq(ln(x) - mu) / (2*sq(sigma)) - ln(x) - ln(sigma * sqrt(2*M_PI)); 00358 } 00359 00360 /**Derivative of the log of the log-normal distribution: 00361 \f[ 00362 \frac{\partial \ln P(x)}{\partial x} = -\frac{1}{x}\left(1 + \frac{\ln x - \mu}{\sigma^2}\right). 00363 \f] 00364 @param x \e x 00365 @param mu \f$\mu\f$ 00366 @param sigma \f$\sigma\f$ 00367 @ingroup gStorm 00368 */ 00369 inline double diff_log_log_normal(double x, double mu, double sigma) 00370 { 00371 return -(1 + (ln(x) - mu)/sq(sigma)) / x; 00372 } 00373 00374 00375 /**Second derivative of the log of the log-normal distribution: 00376 \f[ 00377 \frac{\partial^2 \ln P(x)}{\partial x^2} = \frac{1}{x^2}\left(1 + \frac{\ln x - \mu}{\sigma^2} - \frac{1}{\sigma^2}\right). 00378 \f] 00379 @param x \e x 00380 @param mu \f$\mu\f$ 00381 @param sigma \f$\sigma\f$ 00382 @ingroup gStorm 00383 */ 00384 inline double hess_log_log_normal(double x, double mu, double sigma) 00385 { 00386 return (1 + (ln(x) - mu - 1)/sq(sigma)) / sq(x); 00387 } 00388 00389 00390 #endif