Jpp  master_rocky
the software that should make you happy
JShowerPointSimplex.hh
Go to the documentation of this file.
1 #ifndef JSHOWERPOINTSIMPLEX_INCLUDE
2 #define JSHOWERPOINTSIMPLEX_INCLUDE
3 
4 #include <string>
5 #include <iostream>
6 #include <set>
7 #include <vector>
8 #include <algorithm>
9 #include <memory>
10 #include <math.h>
11 
14 
15 #include "JTrigger/JHit.hh"
16 #include "JTrigger/JTimeslice.hh"
17 #include "JTrigger/JHitR1.hh"
18 #include "JTrigger/JBuildL0.hh"
19 #include "JTrigger/JBuildL2.hh"
20 #include "JTrigger/JAlgorithm.hh"
21 #include "JTrigger/JMatch3G.hh"
22 #include "JTrigger/JBind2nd.hh"
23 
24 #include "JFit/JFitToolkit.hh"
25 #include "JFit/JEstimator.hh"
26 #include "JFit/JPoint4D.hh"
27 #include "JFit/JModel.hh"
28 #include "JFit/JSimplex.hh"
29 #include "JFit/JGandalf.hh"
31 
32 #include "JReconstruction/JEvt.hh"
35 
36 #include "JLang/JPredicate.hh"
37 
38 #include "JGeometry3D/JVector3D.hh"
39 #include "JGeometry3D/JShower3D.hh"
40 #include "JGeometry3D/JShower3E.hh"
41 #include "JGeometry3D/JTrack3D.hh"
42 #include "JGeometry3D/JTrack3E.hh"
44 
46 #include "JDetector/JDetector.hh"
48 
49 #include "JTools/JPermutation.hh"
50 #include "JTools/JRange.hh"
51 
52 /**
53  * \author adomi, vcarretero
54  */
55 namespace JRECONSTRUCTION {}
56 namespace JPP { using namespace JRECONSTRUCTION; }
57 
58 namespace JRECONSTRUCTION {
59 
61  using JFIT::JPoint4D;
62  using JFIT::JGandalf;
63  using JFIT::JRegressor;
64 
65  /**
66  * class to handle the second position fit of the shower reconstruction, mainly dedicated for ORCA
67  */
69 
71  public JRegressor<JPoint4D, JGandalf>
72 
73  {
74 
76  using JRegressor_t::operator();
79 
80  public:
81 
82  /**
83  * Parameterized constructor
84  *
85  * \param parameters struct that holds default-optimized parameters for the reconstruction.
86  * \param router module router, this is built via detector file.
87  * \param debug debug
88  */
90  const JModuleRouter& router,
91  const int debug = 0):
93  JRegressor_t(parameters.sigma_ns),
94  router(router)
95  {
96  using namespace JPP;
97 
99  JRegressor_t::MAXIMUM_ITERATIONS = NMax;
100  JRegressor_t::EPSILON = 1e-4;
101  //this->estimator.reset(new JMEstimatorTukey(tukey_std_dev));
102  this->estimator.reset(getMEstimator(parameters.mestimator));
103  this->parameters.resize(4);
104  this->parameters[0] = JPoint4D::pX();
105  this->parameters[1] = JPoint4D::pY();
106  this->parameters[2] = JPoint4D::pZ();
107  this->parameters[3] = JPoint4D::pT();
108  }
109 
110  /**
111  * Declaration of the member function that actually performs the reconstruction
112  *
113  * \param event which is a JDAQEvent
114  * \param in input fits
115  */
117  {
118  using namespace std;
119  using namespace JPP;
120 
121  const JBuildL0<hit_type> buildL0;
123 
124  buffer_type dataL0;
125  buffer_type dataL1;
126 
127  JEvt out;
128 
129  buildL0(JDAQTimeslice(event, true), router, back_inserter(dataL0));
130  buildL2(JDAQTimeslice(event, false), router, back_inserter(dataL1)); // false = triggered
131 
132  for (JEvt::const_iterator shower = in.begin(); shower != in.end(); ++shower) {
133 
134  JPoint4D vx(getPosition(*shower), shower->getT());
135 
136 
139 
140  data.reserve(dataL0.size() + dataL1.size());
141  for (buffer_type::const_iterator i = dataL1.begin(); i != dataL1.end(); ++i) {
142 
143  if (match(*i)) {
144  data.push_back(*i);
145  }
146  }
147 
148  buffer_type::iterator __end = data.end();
149 
150  for (buffer_type::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
151  if (find_if(data.begin(), __end, make_predicate(&hit_type::getModuleID, i->getModuleID())) == __end) {
152  if (match(*i)) {
153  data.push_back(*i);
154  }
155  }
156  }
157 
158  const int NDF = getCount(data.begin(), data.end()) - this->parameters.size();
159 
160  if (NDF > 0) {
161 
162  double chi2 = (*this)(vx, data.begin(), data.end());
163 
164  JShower3E sh_fit(this->value.getPosition(),
165  JDirection3D(),
166  this->value.getT(),
167  shower->getE());
168  out.push_back(getFit(JHistory(shower->getHistory()).add(JSHOWERPOINTSIMPLEX), sh_fit, getQuality(chi2, NDF), NDF, sh_fit.getE()));
169  }
170  }
171 
172 
173  return out;
174  }
176  };
177 }
178 #endif
179 
Algorithms for hit clustering and sorting.
Auxiliary class to extract a subset of optical modules from a detector.
Data structure for detector geometry and calibration.
Linear fit methods.
Auxiliary methods to evaluate Poisson probabilities and chi2.
Reduced data structure for L1 hit.
Match operator for Cherenkov light from shower in any direction.
int debug
debug level
Definition: JSirene.cc:69
Direct access to module in detector data structure.
Auxiliary class to define a range between two values.
Basic data structure for time and time over threshold information of hit.
Router for direct addressing of module data in detector data structure.
Data structure for set of track fit results.
Fit method based on the Levenberg-Marquardt method.
Definition: JGandalf.hh:86
Data structure for vertex fit.
Definition: JPoint4D.hh:24
static parameter_type pT()
Definition: JPoint4D.hh:61
static parameter_type pZ()
Definition: JPoint4D.hh:60
static parameter_type pX()
Definition: JPoint4D.hh:58
static parameter_type pY()
Definition: JPoint4D.hh:59
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:35
3D track with energy.
Definition: JTrack3E.hh:32
class to handle the second position fit of the shower reconstruction, mainly dedicated for ORCA
JRegressor< JPoint4D, JGandalf > JRegressor_t
JShowerPointSimplex(const JShowerPointSimplexParameters_t &parameters, const JModuleRouter &router, const int debug=0)
Parameterized constructor.
JEvt operator()(const KM3NETDAQ::JDAQEvent &event, const JFIT::JEvt &in)
Declaration of the member function that actually performs the reconstruction.
Template L0 hit builder.
Definition: JBuildL0.hh:38
Template L2 builder.
Definition: JBuildL2.hh:49
Reduced data structure for L1 hit.
Definition: JHitR1.hh:35
int getModuleID() const
Get module identifier.
static const int JSHOWERPOINTSIMPLEX
JPosition3D getPosition(const Vec &pos)
Get position.
double getQuality(const double chi2, const int NDF)
Get quality of fit.
JFit getFit(const int id, const JMODEL::JString &string)
Get fit parameters of string.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
JMEstimator * getMEstimator(const int type)
Get M-Estimator.
Definition: JMEstimator.hh:203
size_t getCount(const array_type< T > &buffer, const JCompare_t &compare)
Count number of unique values.
Definition: JVectorize.hh:261
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Model fits to data.
JFIT::JHistory JHistory
Definition: JHistory.hh:354
Definition: JSTDTypes.hh:14
JHistory & add(const int type)
Add event to history.
Definition: JHistory.hh:295
Template specialisation of class JModel to match hit with bright point.
Definition: JFit/JModel.hh:123
Template definition of a data regressor of given model.
Definition: JRegressor.hh:70
double TMin_ns
minimum time for local coincidences [ns]
double TMaxLocal_ns
time window for local coincidences [ns]
double TMax_ns
maximum time for local coincidences [ns]
double ctMin
minimal cosine space angle between PMT axes
double DMax_m
maximal distance to optical module [m]
Data structure for L2 parameters.