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

Program to perform JFIT::JRegressor<JEnergy, JGandalf> fit using I/O of JFIT::JEvt data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <algorithm>
#include <limits>
#include <map>
#include "TMath.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/definitions/fitparameters.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorSubset.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JFrame.hh"
#include "JTrigger/JTimeslice.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JHitR1.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JTriggerParameters.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JParallelFileScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JFit/JLine1Z.hh"
#include "JFit/JEnergy.hh"
#include "JFit/JEnergyRegressor.hh"
#include "JFit/JFitToolkit.hh"
#include "JFit/JNPEHit.hh"
#include "JFit/JMEstimator.hh"
#include "JFit/JTimeRange.hh"
#include "JFit/JModel.hh"
#include "JReconstruction/JEvt.hh"
#include "JReconstruction/JEvtToolkit.hh"
#include "JReconstruction/JMuonEnergyParameters_t.hh"
#include "JReconstruction/JEnergyCorrection.hh"
#include "JPhysics/JConstants.hh"
#include "JTools/JRange.hh"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JFunctionalMap_t.hh"
#include "JPhysics/JPDFTable.hh"
#include "JPhysics/JPDFToolkit.hh"
#include "JPhysics/JNPETable.hh"
#include "JLang/JComparator.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 perform JFIT::JRegressor<JEnergy, JGandalf> fit using I/O of JFIT::JEvt data.

Author
mdejong

Definition in file JMuonEnergy.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 120 of file JMuonEnergy.cc.

121 {
122  using namespace std;
123  using namespace JPP;
124  using namespace KM3NETDAQ;
125 
126  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
127  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
129 
130  JParallelFileScanner_t inputFile;
132  JLimit_t& numberOfEvents = inputFile.getLimit();
133  string detectorFile;
134  string pdfFile;
135  JEnergyCorrection correct;
137  int debug;
138 
139  try {
140 
141  JParser<> zap("Program to perform fit of muon energy to data.");
142 
143  zap['f'] = make_field(inputFile);
144  zap['o'] = make_field(outputFile) = "energy.root";
145  zap['a'] = make_field(detectorFile);
146  zap['n'] = make_field(numberOfEvents) = JLimit::max();
147  zap['P'] = make_field(pdfFile);
148  zap['E'] = make_field(correct) = JEnergyCorrection(); // off
150  zap['d'] = make_field(debug) = 1;
151 
152  zap(argc, argv);
153  }
154  catch(const exception& error) {
155  FATAL(error.what() << endl);
156  }
157 
158 
159  cout.tie(&cerr);
160 
161 
163 
164  try {
165  load(detectorFile, detector);
166  }
167  catch(const JException& error) {
168  FATAL(error);
169  }
170 
171  const JModuleRouter router(detector);
172 
173  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
174 
175  typedef vector<JHitL0> JDataL0_t;
176  const JBuildL0<JHitL0> buildL0;
177 
178 
179  typedef JRegressor<JEnergy> JRegressor_t;
180 
182  JRegressor_t::T_ns.setRange(parameters.TMin_ns, parameters.TMax_ns);
183 
184  JRegressor_t fit(pdfFile);
185 
186  fit.estimator.reset(getMEstimator(parameters.mestimator));
187 
188  if (!fit.estimator.is_valid()) {
189  FATAL("Invalid M-Estimator." << endl);
190  }
191 
192  if (fit.getRmax() < parameters.roadWidth_m) {
193  parameters.roadWidth_m = fit.getRmax();
194  }
195 
196 
197  outputFile.open();
198 
199  outputFile.put(JMeta(argc, argv));
200 
201  while (inputFile.hasNext()) {
202 
203  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
204 
205  multi_pointer_type ps = inputFile.next();
206 
207  const JDAQEvent* tev = ps;
208  const JEvt* in = ps;
209 
210  summary.update(*tev);
211 
212  // select start values
213 
214  JEvt cp = *in;
215  JEvt out;
216 
217  if (parameters.reprocess) {
219  }
220 
221  cp.select(parameters.numberOfPrefits, qualitySorter);
222 
223  if (!cp.empty()) {
224  cp.select(JHistory::is_event(cp.begin()->getHistory()));
225  }
226 
227 
228  JDataL0_t dataL0;
229 
230  buildL0(*tev, router, true, back_inserter(dataL0));
231 
232 
233  for (JEvt::const_iterator track = cp.begin(); track != cp.end(); ++track) {
234 
235  const JRotation3D R (getDirection(*track));
236  const JLine1Z tz(getPosition (*track).rotate(R), track->getT());
237  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, JRegressor_t::T_ns);
238 
239 
241 
242  for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
243 
244  JHitR1 hit(*i);
245 
246  hit.rotate(R);
247 
248  if (match(hit)) {
249  top.insert(i->getPMTIdentifier());
250  }
251  }
252 
253 
254  double EMin_log = parameters.EMin_log;
255  double EMax_log = parameters.EMax_log;
256  JRange<double> Z_m;
257 
258  if (track->hasW(JSTART_LENGTH_METRES) && track->getW(JSTART_LENGTH_METRES) > 0.0) {
259 
260  Z_m.setLowerLimit(parameters.ZMin_m);
261 
262  //EMin_log = log10(track->getW(JSTART_LENGTH_METRES) * gWater.getA());
263  }
264 
265 
266  vector<JNPEHit> data;
267 
268  const JDetectorSubset_t subdetector(detector, getAxis(*track), parameters.roadWidth_m, Z_m);
269 
270  for (JDetector::const_iterator module = subdetector.begin(); module != subdetector.end(); ++module) {
271 
272  if (summary.hasSummaryFrame(module->getID())) {
273 
274  const JDAQSummaryFrame& frame = summary.getSummaryFrame(module->getID());
275 
276  for (size_t i = 0; i != module->size(); ++i) {
277 
278  if (getDAQStatus(frame, *module, i) &&
279  getPMTStatus(frame, *module, i) &&
280  !module->getPMT(i).has(PMT_DISABLE)) {
281 
282  const double rate_Hz = frame.getRate(i);
283  const size_t count = top.count(JDAQPMTIdentifier(module->getID(), i));
284 
285  data.push_back(JNPEHit(fit.getNPE(module->getPMT(i), rate_Hz), count));
286  }
287  }
288  }
289  }
290 
291  const int NDF = distance(data.begin(), data.end()) - 1;
292 
293  if (NDF >= 0) {
294 
295 
296  // 5-point search between given limits
297 
298  const int N = 5;
299 
300  JResult result[N];
301 
302  for (int i = 0; i != N; ++i) {
303  result[i].x = EMin_log + i * (EMax_log - EMin_log) / (N-1);
304  }
305 
306  map<double, double> buffer;
307 
308  do {
309 
310  int j = 0;
311 
312  for (int i = 0; i != N; ++i) {
313 
314  if (!result[i]) {
315 
316  const JEnergy x = result[i].x;
317  const double chi2 = fit(x, data.begin(), data.end());
318 
319  result[i].chi2 = chi2;
320  buffer[chi2] = x.getE();
321  }
322 
323  if (result[i].chi2 < result[j].chi2) {
324  j = i;
325  }
326  }
327 
328  for (int i = 0; i != N; ++i) {
329  DEBUG(' ' << FIXED(5,2) << result[i].x << ' ' << FIXED(9,3) << result[i].chi2);
330  }
331  DEBUG(endl);
332 
333 
334  // squeeze range
335 
336  switch (j) {
337 
338  case 0:
339  result[0] = JResult(2*result[0].x - result[2].x);
340  result[4] = result[2];
341  result[2] = result[0];
342  break;
343 
344  case 1:
345  result[4] = result[2];
346  result[2] = result[1];
347  break;
348 
349  case 2:
350  result[0] = result[1];
351  result[4] = result[3];
352  break;
353 
354  case 3:
355  result[0] = result[2];
356  result[2] = result[3];
357  break;
358 
359  case 4:
360  result[0] = result[2];
361  result[2] = result[4];
362  result[4] = JResult(2*result[4].x - result[2].x);
363  break;
364  }
365 
366  result[1] = JResult(0.5 * (result[0].x + result[2].x));
367  result[3] = JResult(0.5 * (result[2].x + result[4].x));
368 
369  } while (result[4].x - result[0].x > parameters.resolution);
370 
371 
372  if (result[1].chi2 != result[3].chi2) {
373 
374  result[2].x += 0.25 * (result[3].x - result[1].x) * (result[1].chi2 - result[3].chi2) / (result[1].chi2 + result[3].chi2 - 2*result[2].chi2);
375  result[2].chi2 = fit(result[2].x, data.begin(), data.end());
376  }
377 
378  const double chi2 = result[2].chi2;
379  const double E = result[2].x.getE();
380 
381  // calculate additional variables
382 
383  double Emin = numeric_limits<double>::max();
384  double Emax = numeric_limits<double>::lowest();
385 
386  for (map<double, double>::const_iterator i = buffer.begin(); i != buffer.end() && i->first <= chi2 + 1.0; ++i) {
387  if (i->second < Emin) { Emin = i->second; }
388  if (i->second > Emax) { Emax = i->second; }
389  }
390 
391  const double mu_range = gWater(E); // range of a muon with the reconstructed energy
392 
393  double noise_likelihood = 0.0; // log-likelihood of every hit being from K40
394  int number_of_hits = 0; // number of hits selected for JEnergy
395 
396  for (vector<JNPEHit>::const_iterator i = data.begin(); i != data.end(); ++i) {
397  noise_likelihood += log10(i->getP()); // probability of getting the observed multiplicity with K40
398  number_of_hits += i->getN(); // hit multiplicity
399  }
400 
401 
402  out.push_back(JFit(*track).add(JMUONENERGY));
403 
404  // set corrected energy
405 
406  out.rbegin()->setE(correct(E));
407 
408  // set additional values
409 
410  out.rbegin()->setW(track->getW());
411  out.rbegin()->setW(JENERGY_ENERGY, E);
412  out.rbegin()->setW(JENERGY_CHI2, chi2);
413  out.rbegin()->setW(JENERGY_MUON_RANGE_METRES, mu_range);
414  out.rbegin()->setW(JENERGY_NOISE_LIKELIHOOD, noise_likelihood);
415  out.rbegin()->setW(JENERGY_NDF, NDF);
416  out.rbegin()->setW(JENERGY_NUMBER_OF_HITS, number_of_hits);
417  out.rbegin()->setW(JENERGY_MINIMAL_ENERGY, Emin);
418  out.rbegin()->setW(JENERGY_MAXIMAL_ENERGY, Emax);
419  }
420  }
421 
422  // apply default sorter
423 
424  sort(out.begin(), out.end(), qualitySorter);
425 
426  copy(in->begin(), in->end(), back_inserter(out));
427 
428  outputFile.put(out);
429  }
430  STATUS(endl);
431 
433 
434  io >> outputFile;
435 
436  outputFile.close();
437 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
double getE() const
Get energy.
Definition: JEnergy.hh:170
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
double getRate(const int tdc, const double factor=1.0) const
Get count rate.
static const int JENERGY_ENERGY
uncorrected energy [GeV] from JEnergy.cc
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Regressor function object for JEnergy fit.
Template specialisation of L0 builder for JHitL0 data type.
Definition: JBuildL0.hh:102
static const int JENERGY_MAXIMAL_ENERGY
maximal energy [GeV] from JEnergy.cc
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
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
Auxiliary class to test history.
Definition: JHistory.hh:152
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:323
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
JFit & add(const int type)
Add event to history.
string outputFile
static const int JENERGY_CHI2
chi2 from JEnergy.cc
static const int JENERGY_NOISE_LIKELIHOOD
log likelihood of every hit being K40 from JEnergy.cc
static const int JENERGY_NDF
number of degrees of freedom from JEnergy.cc
static const int JENERGY_MINIMAL_ENERGY
minimal energy [GeV] from JEnergy.cc
Type list.
Definition: JTypeList.hh:22
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Detector file.
Definition: JHead.hh:196
JDirection3D getDirection(const Vec &dir)
Get direction.
Auxiliary class for correction of energy determined by JEnergy.cc.
Acoustic event fit.
static const int JENERGY_NUMBER_OF_HITS
number of hits from JEnergy.cc
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Data storage class for rate measurements of all PMTs in one module.
Auxiliary class to test history.
Definition: JHistory.hh:104
JAxis3D getAxis(const Trk &track)
Get axis.
return result
Definition: JPolint.hh:727
static const int PMT_DISABLE
KM3NeT Data Definitions v2.0.0-14-gbeccebb https://git.km3net.de/common/km3net-dataformat.
Definition: pmt_status.hh:12
JPosition3D getPosition(const Vec &pos)
Get position.
File router for fast addressing of summary data.
Auxiliary class for simultaneously handling light yields and response of PMT.
Definition: JNPEHit.hh:19
int debug
debug level
Definition: JSirene.cc:63
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:40
Detector subset without binary search functionality.
Data structure for fit parameters.
#define FATAL(A)
Definition: JMessage.hh:67
bool getPMTStatus(const JStatus &status)
Test status of PMT.
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.
std::vector< int > count
Definition: JAlgorithm.hh:184
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:27
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
then cp
Reduced data structure for L1 hit.
Definition: JHitR1.hh:31
Data structure for fit of energy.
Definition: JEnergy.hh:28
Object reading from a list of files.
int j
Definition: JPolint.hh:666
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
JMEstimator * getMEstimator(const int type)
Get M-Estimator.
Definition: JMEstimator.hh:166
do set_variable DETECTOR_TXT $WORKDIR detector
bool getDAQStatus(const JDAQFrameStatus &frame, const JStatus &status)
Test status of DAQ.
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 source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:38
JPosition3D & rotate(const JRotation3D &R)
Rotate.
Definition: JPosition3D.hh:186
then usage $script[input file[working directory[option]]] nWhere option can be N
Definition: JMuonPostfit.sh:37
static const int JENERGY_MUON_RANGE_METRES
range of a muon with the reconstructed energy [m] from JEnergy.cc
void setLowerLimit(argument_type x)
Set lower limit.
Definition: JRange.hh:224
bool qualitySorter(const JRECONSTRUCTION::JFit &first, const JRECONSTRUCTION::JFit &second)
Comparison of fit results.
static const int JMUONENERGY
then usage $script[input file[working directory[option]]] nWhere option can be E
Definition: JMuonPostfit.sh:37
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62