Jpp master_rocky-44-g75b7c4f75
the software that should make you happy
Loading...
Searching...
No Matches
JLine3ZRegressor.hh
Go to the documentation of this file.
1#ifndef __JFIT__JLINE3ZREGRESSOR__
2#define __JFIT__JLINE3ZREGRESSOR__
3
4#include <string>
5
10#include "JPhysics/JGeane.hh"
14#include "JTools/JRange.hh"
17#include "JMath/JZero.hh"
19#include "JFit/JTimeRange.hh"
20#include "JFit/JLine3Z.hh"
21#include "JFit/JPMTW0.hh"
22#include "JFit/JSimplex.hh"
23#include "JFit/JGandalf.hh"
24#include "JFit/JMEstimator.hh"
25#include "JFit/JRegressor.hh"
27#include "Jeep/JMessage.hh"
28
29
30/**
31 * \file
32 * Data regression method for JFIT::JLine3Z.
33 * \author mdejong
34 */
35namespace JFIT {}
36namespace JPP { using namespace JFIT; }
37
38namespace JFIT {
39
40 /**
41 * Regressor function object for JLine3Z fit using JSimplex minimiser.
42 */
43 template<>
44 struct JRegressor<JLine3Z, JSimplex> :
45 public JAbstractRegressor<JLine3Z, JSimplex>
46 {
47 using JAbstractRegressor<JLine3Z, JSimplex>::operator();
48
49 /**
50 * Constructor.
51 *
52 * \param sigma time resolution of hit [ns]
53 */
54 JRegressor(double sigma) :
55 estimator(new JMEstimatorNormal())
56 {
57 this->sigma = sigma;
58 }
59
60
61 /**
62 * Fit function.
63 * This method is used to determine the chi2 of given hit with respect to trajectory of muon.
64 *
65 * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
66 * - double getX(); // [m]
67 * - double getY(); // [m]
68 * - double getZ(); // [m]
69 * - double getT(); // [ns]
70 *
71 * \param track track
72 * \param hit hit
73 * \return chi2
74 */
75 template<class JHit_t>
76 double operator()(const JLine3Z& track, const JHit_t& hit) const
77 {
78 using namespace JPP;
79
80 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
81
82 D.sub(track.getPosition());
83
84 const double z = D.getDot(track.getDirection());
85 const double R = sqrt(D.getLengthSquared() - z*z);
86
87 const double t1 = track.getT() + (z + R * getKappaC()) * getInverseSpeedOfLight();
88
89 const double u = (t1 - hit.getT()) / sigma;
90
91 return estimator->getRho(u) * hit.getW();
92 }
93
94 JLANG::JSharedPointer<JMEstimator> estimator; //!< M-Estimator function
95 double sigma; //!< Time resolution [ns]
96 };
97
98
99 /**
100 * Regressor function object for JLine3Z fit using JGandalf minimiser.
101 */
102 template<>
103 struct JRegressor<JLine3Z, JGandalf> :
104 public JAbstractRegressor<JLine3Z, JGandalf>,
105 public JRegressorStorage <JLine3Z, JGandalf>
106 {
107 using JAbstractRegressor<JLine3Z, JGandalf>::operator();
108
110
111 /**
112 * Default constructor
113 */
114 JRegressor() :
116 pdf(getPDF()),
117 npe(getNPE())
118 {}
119
120
121 /**
122 * Constructor.
123 *
124 * The PDF file descriptor should contain the wild card character JPHYSICS::WILDCARD which
125 * will be replaced by the corresponding PDF types listed in JRegressor<JLine3Z, JGandalf>::pdf_t.
126 *
127 * The <tt>TTS</tt> corresponds to the additional time smearing applied to the PDFs.
128 *
129 * \param fileDescriptor PDF file descriptor
130 * \param TTS TTS [ns]
131 * \param numberOfPoints number of points for Gauss-Hermite integration of TTS
132 * \param epsilon precision for Gauss-Hermite integration of TTS
133 */
134 JRegressor(const std::string& fileDescriptor,
135 const double TTS,
136 const int numberOfPoints = 25,
137 const double epsilon = 1.0e-10) :
138 JRegressorStorage_t(fileDescriptor, TTS, numberOfPoints, epsilon),
139 pdf(getPDF()),
140 npe(getNPE())
141 {}
142
143
144 /**
145 * Constructor.
146 *
147 * \param storage PDF storage
148 */
150 pdf(storage.getPDF()),
151 npe(storage.getNPE())
152 {}
153
154
155 /**
156 * Fit function.
157 * This method is used to determine the chi2 and gradient of given hit with respect to trajectory of muon.
158 *
159 * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
160 * - double getX(); // [m]
161 * - double getY(); // [m]
162 * - double getZ(); // [m]
163 * - double getDX(); // [u]
164 * - double getDY(); // [u]
165 * - double getDZ(); // [u]
166 * - double getT(); // [ns]
167 *
168 * \param track track
169 * \param hit hit
170 * \return chi2 and gradient
171 */
172 template<class JHit_t>
173 result_type operator()(const JLine3Z& track, const JHit_t& hit) const
174 {
175 using namespace JPP;
176
177 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
178 JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
179
180 D.sub(track.getPosition());
181
182 const double z = D.getDot(track.getDirection());
183 const double x = D.getX() - z * track.getDX();
184 const double y = D.getY() - z * track.getDY();
185 const double R2 = D.getLengthSquared() - z*z;
186 const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
187
188 const double t1 = track.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
189
190 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
191
192 const double theta = U.getTheta();
193 const double phi = fabs(U.getPhi());
194
195 const double E = gWater.getE(E_GeV, z);
196 const double dt = T_ns.constrain(hit.getT() - t1);
197
198 JPDF_t::result_type H0 = getH0(hit.getR(), dt);
199 JPDF_t::result_type H1 = getH1(E, R, theta, phi, dt);
200
201 if (H1.V >= Vmax_npe) {
202 H1 *= Vmax_npe / H1.V;
203 }
204
205 H1 += H0; // signal + background
206
207 result_type result;
208
209 result.chi2 = H1.getChi2() - H0.getChi2(); // Likelihood ratio
210
211 const double wc = 1.0 - getTanThetaC() * z / R; // d(ct1)/d(z)
212
213 result.gradient = JLine3Z(JLine1Z(JPosition3D(-getTanThetaC() * D.getX() / R, // d(ct1)/d(x)
214 -getTanThetaC() * D.getY() / R, // d(ct1)/d(y)
215 0.0),
216 getSpeedOfLight()), // d(ct1)/d(t)
217 JVersor3Z(wc * (D.getX() - D.getZ()*track.getDX()/track.getDZ()), // d(ct1)/d(dx)
218 wc * (D.getY() - D.getZ()*track.getDY()/track.getDZ()))); // d(ct1)/d(dy)
219
220 result.gradient.mul(getInverseSpeedOfLight() * (H1.getDerivativeOfChi2() -
221 H0.getDerivativeOfChi2())); // x d(chi2)/d(ct1)
222
223 return result;
224 }
225
226
227 /**
228 * Fit function.
229 * This method is used to determine the chi2 and gradient of given PMT with respect to trajectory of muon.
230 *
231 * \param track track
232 * \param pmt pmt
233 * \return chi2 and gradient
234 */
235 result_type operator()(const JLine3Z& track, const JPMTW0& pmt) const
236 {
237 using namespace JGEOMETRY3D;
238 using namespace JPHYSICS;
239
240 JPosition3D D(pmt.getPosition());
241 JDirection3D U(pmt.getDirection());
242
243 D.sub(track.getPosition());
244
245 const double z = D.getDot(track.getDirection());
246 const double x = D.getX() - z * track.getDX();
247 const double y = D.getY() - z * track.getDY();
248 const double R2 = D.getLengthSquared() - z*z;
249 const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
250
251 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
252
253 const double theta = U.getTheta();
254 const double phi = fabs(U.getPhi());
255
256 const double E = gWater.getE(E_GeV, z);
257
258 JNPE_t::result_type H0 = getH0(pmt.getR());
259 JNPE_t::result_type H1 = getH1(E, R, theta, phi);
260
261 if (H1.f >= Vmax_npe) {
262 H1 *= Vmax_npe / H1.f;
263 }
264
265 H1 += H0;
266
267 const bool hit = pmt.getN() != 0;
268 const double u = H1.getChi2(hit);
269
270 result_type result;
271
272 result.chi2 = estimator->getRho(u);
273
274 result.gradient = JLine3Z(JLine1Z(JPosition3D(-D.getX() / R, // d(R)/d(x)
275 -D.getY() / R, // d(R)/d(y)
276 0.0),
277 0.0),
278 JVersor3Z(-z * D.getX() / R, // d(R)/d(dx)
279 -z * D.getY() / R)); // d(R)/d(dy)
280
281 result.gradient.mul(estimator->getPsi(u));
282 result.gradient.mul(H1.getDerivativeOfChi2(hit)); // x d(chi2)/d(R)
283
284 return result;
285 }
286
287
288 /**
289 * Get background hypothesis value for time differentiated PDF.
290 *
291 * \param R_Hz rate [Hz]
292 * \param t1 time [ns]
293 * \return hypothesis value
294 */
295 JPDF_t::result_type getH0(const double R_Hz,
296 const double t1) const
297 {
298 return JPDF_t::result_type(R_Hz * 1e-9, t1, T_ns);
299 }
300
301
302 /**
303 * Get signal hypothesis value for time differentiated PDF.
304 *
305 * \param E muon energy at minimum distance of approach [GeV]
306 * \param R minimum distance of approach [m]
307 * \param theta PMT zenith angle [rad]
308 * \param phi PMT azimuth angle [rad]
309 * \param t1 arrival time relative to Cherenkov hypothesis [ns]
310 * \return hypothesis value
311 */
312 JPDF_t::result_type getH1(const double E,
313 const double R,
314 const double theta,
315 const double phi,
316 const double t1) const
317 {
318 using namespace std;
319 using namespace JPP;
320
321 JPDF_t::result_type h1 = zero;
322
323 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
324
325 if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
326
327 JPDF_t::result_type y1 = pdf[i](max(R, pdf[i].getXmin()), theta, phi, t1);
328
329 // safety measures
330
331 if (y1.f <= 0.0) {
332 y1.f = 0.0;
333 y1.fp = 0.0;
334 }
335
336 if (y1.v <= 0.0) {
337 y1.v = 0.0;
338 }
339
340 // energy dependence
341
342 if (is_deltarays(pdf_t[i])) {
343 y1 *= getDeltaRaysFromMuon(E);
344 } else if (is_bremsstrahlung(pdf_t[i])) {
345 y1 *= E;
346 }
347
348 h1 += y1;
349 }
350 }
351
352 return h1;
353 }
354
355
356 /**
357 * Get background hypothesis value for time integrated PDF.
358 *
359 * \param R_Hz rate [Hz]
360 * \return hypothesis value
361 */
362 JNPE_t::result_type getH0(const double R_Hz) const
363 {
364 return JNPE_t::result_type(R_Hz * 1e-9 * T_ns.getLength(), 0.0);
365 }
366
367
368 /**
369 * Get signal hypothesis value for time integrated PDF.
370 *
371 * \param E muon energy at minimum distance of approach [GeV]
372 * \param R minimum distance of approach [m]
373 * \param theta PMT zenith angle [rad]
374 * \param phi PMT azimuth angle [rad]
375 * \return hypothesis value
376 */
377 JNPE_t::result_type getH1(const double E,
378 const double R,
379 const double theta,
380 const double phi) const
381 {
382 using namespace std;
383 using namespace JPP;
384
385 JNPE_t::result_type h1 = JMATH::zero;
386
387 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
388
389 if (!npe[i].empty() && R <= npe[i].getXmax()) {
390
391 try {
392
393 JNPE_t::result_type y1 = npe[i](max(R, npe[i].getXmin()), theta, phi);
394
395 // safety measures
396
397 if (y1.f > 0.0) {
398
399 // energy dependence
400
401 if (is_deltarays(pdf_t[i])) {
402 y1 *= getDeltaRaysFromMuon(E);
403 } else if (is_bremsstrahlung(pdf_t[i])) {
404 y1 *= E;
405 }
406
407 h1 += y1;
408 }
409 }
410 catch(JException& error) {
411 ERROR(error << endl);
412 }
413 }
414 }
415
416 return h1;
417 }
418
419
420 /**
421 * Get maximal road width of PDF.
422 *
423 * \return road width [m]
424 */
425 inline double getRmax() const
426 {
427 double xmax = 0.0;
428
429 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
430 if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
431 xmax = pdf[i].getXmax();
432 }
433 }
434
435 return xmax;
436 }
437
438
439 static JTimeRange T_ns; //!< Time window with respect to Cherenkov hypothesis [ns]
440 static double Vmax_npe; //!< Maximal integral of PDF [npe]
441 static double Rmin_m; //!< Minimal distance of [m]
442
443 const JPDFs_t& pdf; //!< PDF
444 const JNPEs_t& npe; //!< PDF
445
446 double E_GeV = 0.0; //!< Energy of muon at vertex [GeV]
447
448 JLANG::JSharedPointer<JMEstimator> estimator = new JMEstimatorNormal(); //!< M-Estimator function
449 };
450
451
452 /**
453 * Default values.
454 */
456 double JRegressor<JLine3Z, JGandalf>::Vmax_npe = std::numeric_limits<double>::max();
458}
459
460#endif
Various implementations of functional maps.
Energy loss of muon.
Maximum likelihood estimator (M-estimators).
General purpose messaging.
#define ERROR(A)
Definition JMessage.hh:66
Auxiliary methods for PDF calculations.
Numbering scheme for PDF types.
Physics constants.
Auxiliary class to define a range between two values.
General purpose data regression method.
int numberOfPoints
Definition JResultPDF.cc:22
Definition of zero value for any class.
Fit method based on the Levenberg-Marquardt method.
Data structure for fit of straight line paralel to z-axis.
Definition JLine1Z.hh:29
Data structure for fit of straight line in positive z-direction.
Definition JLine3Z.hh:40
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
Definition JLine3Z.hh:255
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition JLine3Z.hh:234
Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++,...
Definition JSimplex.hh:44
Data structure for direction in three dimensions.
const JDirection3D & getDirection() const
Get direction.
Data structure for position in three dimensions.
const JPosition3D & getPosition() const
Get position.
Rotation around Z-axis.
Data structure for normalised vector in positive z-direction.
Definition JVersor3Z.hh:41
double getDZ() const
Get z direction.
Definition JVersor3Z.hh:169
double getDY() const
Get y direction.
Definition JVersor3Z.hh:158
double getDX() const
Get x direction.
Definition JVersor3Z.hh:147
General exception.
Definition JException.hh:24
The template JSharedPointer class can be used to share a pointer to an object.
const double sigma[]
const double xmax
const double epsilon
double getNPE(const Hit &hit)
Get true charge of hit.
Auxiliary classes and methods for linear and iterative data regression.
Definition JEnergy.hh:15
Auxiliary classes and methods for 3D geometrical objects and operations.
Definition JAngle3D.hh:19
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
Auxiliary methods for light properties of deep-sea water.
double getDeltaRaysFromMuon(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit muon track length.
double getKappaC()
Get average R-dependence of arrival time of Cherenkov light (a.k.a.
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Definition JPDFTypes.hh:151
bool is_bremsstrahlung(const int pdf)
Test if given PDF type corresponds to Cherenkov light from Bremsstrahlung.
Definition JPDFTypes.hh:137
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
double u[N+1]
Definition JPolint.hh:865
Abstract class for global fit method.
Definition JRegressor.hh:79
Normal M-estimator.
Auxiliary class for handling PMT geometry, rate and response.
Definition JPMTW0.hh:24
int getN() const
Get number of hits.
Definition JPMTW0.hh:67
double getR() const
Get rate.
Definition JPMTW0.hh:56
Template data structure for storage for PDF tables.
Template definition of a data regressor of given model.
Definition JRegressor.hh:70
Auxiliary class to set-up Hit.
Definition JSirene.hh:58