Jpp  17.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JMuonPDF.cc File Reference

Program to evaluate hit probabilities. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JParallelFileScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JPhysics/JPDF_t.hh"
#include "JPhysics/JNPE_t.hh"
#include "JFit/JLine1Z.hh"
#include "JFit/JModel.hh"
#include "JReconstruction/JHitW0.hh"
#include "JReconstruction/JEvt.hh"
#include "JReconstruction/JEvtToolkit.hh"
#include "JReconstruction/JMuonGandalfParameters_t.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Program to evaluate hit probabilities.

Author
mdejong

Definition in file JMuonPDF.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 69 of file JMuonPDF.cc.

70 {
71  using namespace std;
72  using namespace JPP;
73  using namespace KM3NETDAQ;
74 
75  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
76  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
77 
78  JParallelFileScanner_t inputFile;
79  JLimit_t& numberOfEvents = inputFile.getLimit();
80  string detectorFile;
81  string pdfFile;
83  int debug;
84 
85  try {
86 
87  parameters.numberOfPrefits = 1;
88 
89  JParser<> zap("Program to evaluate hit probabilities.");
90 
91  zap['f'] = make_field(inputFile);
92  zap['a'] = make_field(detectorFile);
93  zap['n'] = make_field(numberOfEvents) = JLimit::max();
94  zap['P'] = make_field(pdfFile);
96  zap['d'] = make_field(debug) = 1;
97 
98  zap(argc, argv);
99  }
100  catch(const exception& error) {
101  FATAL(error.what() << endl);
102  }
103 
104  //setDAQLongprint(debug >= JEEP::debug_t);
105 
106  cout.tie(&cerr);
107 
108  if (parameters.numberOfPrefits != 1) {
109  WARNING("Number of prefits " << parameters.numberOfPrefits << " != " << 1 << endl);
110  }
111 
113 
114  try {
115  load(detectorFile, detector);
116  }
117  catch(const JException& error) {
118  FATAL(error);
119  }
120 
121  const JModuleRouter router(detector);
122 
123  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
124 
125  const JMuonPDF_t pdf(pdfFile, parameters.TTS_ns);
126  const JMuonNPE_t npe(pdfFile);
127 
128  const JTimeRange T_ns(parameters.TMin_ns, parameters.TMax_ns);
129 
130  typedef vector<JHitL0> JDataL0_t;
131  typedef vector<JHitW0> JDataW0_t;
132 
133  const JBuildL0<JHitL0> buildL0;
134 
135  while (inputFile.hasNext()) {
136 
137  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
138 
139  multi_pointer_type ps = inputFile.next();
140 
141  JDAQEvent* tev = ps;
142  JEvt* in = ps;
143 
144  DEBUG("event: " << *tev << endl);
145 
146  summary.update(*tev);
147 
148  in->select(parameters.numberOfPrefits, qualitySorter);
149 
150  JDataL0_t dataL0;
151 
152  buildL0(*tev, router, true, back_inserter(dataL0));
153 
154  for (JEvt::const_iterator track = in->begin(); track != in->end(); ++track) {
155 
156  DEBUG("track: " << *track << endl);
157 
158  JRotation3D R (getDirection(*track));
159  JLine1Z tz(getPosition (*track).rotate(R), track->getT());
160  JRange<double> Z_m;
161 
162  if (track->hasW(JSTART_LENGTH_METRES)) {
163  Z_m = JRange<double>(0.0, track->getW(JSTART_LENGTH_METRES)) + JRange<double>(parameters.ZMin_m, parameters.ZMax_m);
164  }
165 
166  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns, Z_m);
167 
168  // hit selection based on start value
169 
170  JDataW0_t data;
171 
172  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
173 
174  JHitW0 hit(*i, summary.getRate(*i));
175 
176  hit.rotate(R);
177 
178  if (match(hit)) {
179  data.push_back(hit);
180  }
181  }
182 
183  // select first hit in PMT
184 
185  sort(data.begin(), data.end(), compare);
186 
187  JDataW0_t::iterator __end = unique(data.begin(), data.end(), equal_to<JDAQPMTIdentifier>());
188 
189 
190  double E_GeV = parameters.E_GeV;
191  /*
192  if (track->getE() > 0.1) {
193  E_GeV = track->getE();
194  }
195  */
196 
197  DEBUG("line1z: "
198  << FIXED(12,3) << tz.getX() << ' '
199  << FIXED(12,3) << tz.getY() << ' '
200  << FIXED(12,3) << tz.getZ() << ' '
201  << FIXED(12,3) << tz.getT() << ' '
202  << FIXED(12,1) << E_GeV << ' '
203  << FIXED( 8,3) << track->getQ() << endl);
204 
205  // move track to bary-center of hits
206  /*
207  int N = 0;
208  double z = 0.0;
209 
210  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
211  N += 1;
212  z += hit->getZ();
213  }
214 
215  if (N != 0) {
216  tz.setZ(z/N, getSpeedOfLight());
217  }
218  */
219 
220  double Q = 0.0;
221 
222  for (JDataW0_t::const_iterator hit = data.begin(); hit != __end; ++hit) {
223 
224  const double x = hit->getX() - tz.getX();
225  const double y = hit->getY() - tz.getY();
226  const double z = hit->getZ() - tz.getZ();
227  const double R = sqrt(x*x + y*y);
228 
229  const double t1 = tz.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
230 
231  JDirection3D u(hit->getDX(), hit->getDY(), hit->getDZ()); // PMT orientation
232 
233  u.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
234 
235  const double theta = u.getTheta();
236  const double phi = fabs(u.getPhi()); // rotational symmetry of Cherenkov cone
237 
238  //const double E = gWater.getE(E_GeV, z); // correct for energy loss
239  const double E = E_GeV;
240  const double dt = T_ns.constrain(hit->getT() - t1);
241 
242  JMuonPDF_t::result_type H1 = pdf.calculate(E, R, theta, phi, dt);
243  JMuonPDF_t::result_type H0(hit->getR() * 1e-9, 0.0, T_ns);
244 
245  H1 += H0; // signal + background
246 
247  const double chi2 = H1.getChi2() - H0.getChi2(); // likelihood ratio
248 
249  DEBUG("hit: "
250  << setw(10) << hit->getModuleID() << ':' << setw( 2) << setfill('0') << hit->getPMTAddress() << setfill(' ') << ' '
251  << FIXED(12,1) << E << ' '
252  << FIXED( 9,1) << R << ' '
253  << FIXED( 6,4) << theta << ' '
254  << FIXED( 6,4) << phi << ' '
255  << FIXED( 8,3) << dt << ' '
256  << FIXED(12,3) << chi2 << endl);
257 
258  Q += getQuality(chi2);
259  }
260 
261  DEBUG("quality: " << FIXED(8,3) << Q << ' ' << distance(data.begin(), __end) << endl);
262 
263 
264  double Y = 0.0;
265 
266  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
267 
268  JPosition3D pos(i->getPosition());
269 
270  pos.transform(R, tz.getPosition());
271 
272  if (pos.getX() <= parameters.roadWidth_m) {
273 
274  JModule module = *i;
275 
276  module.transform(R, tz.getPosition());
277 
278  for (JModule::iterator pmt = module.begin(); pmt != module.end(); ++pmt) {
279 
280  const double R = sqrt(pmt->getX()*pmt->getX() + pmt->getY()*pmt->getY());
281  const double theta = pmt->getTheta();
282  const double phi = fabs(pmt->getPhi());
283  const double y = npe.calculate(1.0, R, theta, phi); // MIP
284 
285  DEBUG("PMT: "
286  << setw(10) << module.getID() << ':' << setw( 2) << setfill('0') << distance(module.begin(),pmt) << setfill(' ') << ' '
287  << FIXED(9,1) << R << ' '
288  << FIXED(6,4) << theta << ' '
289  << FIXED(6,4) << phi << ' '
290  << SCIENTIFIC(7,2) << y << endl);
291 
292  Y += y;
293  }
294  }
295  }
296 
297  DEBUG("npe: " << FIXED(7,1) << Y << endl);
298  }
299  }
300  STATUS(endl);
301 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
Q(UTCMax_s-UTCMin_s)-livetime_s
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:33
Data structure for a composite optical module.
Definition: JModule.hh:68
then usage E
Definition: JMuonPostfit.sh:35
double getQuality(const double chi2, const int NDF)
Get quality of fit.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:102
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
Router for direct addressing of module data in detector data structure.
Rotation matrix.
Definition: JRotation3D.hh:111
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
General purpose class for parallel reading of objects from a single file or multiple files...
Template specialisation of class JModel to match hit with muon trajectory along z-axis.
Definition: JFit/JModel.hh:34
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
void transform(const JRotation3D &R, const JVector3D &pos)
Transform position.
Definition: JPosition3D.hh:331
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
Rotation around Z-axis.
Definition: JRotation3D.hh:85
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
JAxis3D & rotate(const JRotation3D &R)
Rotate axis.
Definition: JAxis3D.hh:225
Detector file.
Definition: JHead.hh:224
JDirection3D getDirection(const Vec &dir)
Get direction.
JFunction1D_t::result_type result_type
Definition: JPDF_t.hh:145
Acoustic event fit.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
int getID() const
Get identifier.
Definition: JObjectID.hh:50
JPosition3D getPosition(const Vec &pos)
Get position.
File router for fast addressing of summary data.
int debug
debug level
Definition: JSirene.cc:66
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:43
Auxiliary data structure for muon PDF.
Definition: JPDF_t.hh:135
#define FATAL(A)
Definition: JMessage.hh:67
Auxiliary class for a hit with background rate value.
Definition: JHitW0.hh:21
static const int JSTART_LENGTH_METRES
distance between first and last hits in metres from JStart.cc
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
void transform(const JRotation3D &R, const JVector3D &pos)
Transformation of geometry (see method JGEOMETRY3D::JPosition3D::transform(const JRotation3D&amp;, const JVector3D&amp;)).
Definition: JModule.hh:346
const double getInverseSpeedOfLight()
Get inverse speed of light.
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:27
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do set_variable DETECTOR_TXT $WORKDIR detector
double u[N+1]
Definition: JPolint.hh:755
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
JPosition3D & rotate(const JRotation3D &R)
Rotate.
Definition: JPosition3D.hh:186
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:484
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62