Jpp  17.3.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
examples/JReconstruction/JMuonStart.cc File Reference

Program to test JMuonStart.cc. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <map>
#include <algorithm>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/online/JDAQ.hh"
#include "km3net-dataformat/online/JDAQPMTIdentifier.hh"
#include "km3net-dataformat/online/JDAQSummaryFrame.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDetector/JPMTRouter.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JTriggerToolkit.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JTriggeredFileScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JFit/JLine1Z.hh"
#include "JFit/JNPEHit.hh"
#include "JFit/JModel.hh"
#include "JPhysics/KM3NeT.hh"
#include "JPhysics/JGeane.hh"
#include "JReconstruction/JHitW0.hh"
#include "JReconstruction/JEvt.hh"
#include "JReconstruction/JEvtToolkit.hh"
#include "JReconstruction/JMuonStartParameters_t.hh"
#include "JReconstruction/JStart.hh"
#include "JROOT/JManager.hh"
#include "JTools/JAbstractHistogram.hh"
#include "JLang/JPredicate.hh"
#include "JLang/JComparator.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 test JMuonStart.cc.

Author
mdejong

Definition in file examples/JReconstruction/JMuonStart.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 117 of file examples/JReconstruction/JMuonStart.cc.

118 {
119  using namespace std;
120  using namespace JPP;
121  using namespace KM3NETDAQ;
122 
123  typedef JTriggeredFileScanner<JEvt, JSingleFileScanner> JTriggeredFileScanner_t;
124  typedef JTriggeredFileScanner_t::multi_pointer_type multi_pointer_type;
125 
126  typedef JAbstractHistogram<double> JHistogram_t;
127 
128  JTriggeredFileScanner_t inputFile;
129  JLimit_t& numberOfEvents = inputFile.getLimit();
130  string detectorFile;
131  string outputFile;
133  JK40Rates rates_Hz;
134  JHistogram_t histogram;
135  int debug;
136 
137  try {
138 
139  parameters.numberOfPrefits = 1;
140 
141  JParser<> zap("Program to test JMuonStart.");
142 
143  zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)");
144  zap['a'] = make_field(detectorFile);
145  zap['n'] = make_field(numberOfEvents) = JLimit::max();
146  zap['o'] = make_field(outputFile) = "start.root";
148  zap['B'] = make_field(rates_Hz) = KM3NET::getK40Rates();
149  zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t(200, -250.0, +250.0);
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 
160 
161  try {
162  load(detectorFile, detector);
163  }
164  catch(const JException& error) {
165  FATAL(error);
166  }
167 
168  const JModuleRouter moduleRouter(detector);
169  const JPMTRouter pmtRouter (detector);
170 
171  JSummaryFileRouter summary(inputFile, parameters.R_Hz);
172 
173  const JTimeRange T_ns (parameters.TMin_ns, parameters.TMax_ns);
174  const JStart start(parameters.Pmin1, parameters.Pmin2, parameters.Nmax2);
175 
176  typedef JHitL0 hit_type;
177 
178  typedef vector<hit_type> buffer_type;
179  const JBuildL0<hit_type> buildL0;
180 
181 
182  JHead head;
183 
184  try {
185  head = getHeader(inputFile);
186  }
187  catch(const JException& error) {
188  FATAL(error);
189  }
190 
191  const JPosition3D center = get<JPosition3D>(head);
192 
193  NOTICE("center: " << center << endl);
194 
195  JManager<string, TH1D> H1(new TH1D("[%].1L", NULL,
196  histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit()));
197  JManager<string, TH2D> H2(new TH2D("[%].2D", NULL,
198  histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit(),
199  histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit()));
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  JDAQEvent* tev = ps;
208  JEvt* evt = ps;
209  Evt* event = ps;
210 
211  summary.update(*tev);
212 
213  vector<Trk>::iterator muon = event->mc_trks.end();
214 
215  for (vector<Trk>::iterator track = event->mc_trks.begin(); track != event->mc_trks.end(); ++track) {
216 
217  if (is_muon(*track)) {
218  if (muon == event->mc_trks.end() || track->E > muon->E) {
219  muon = track;
220  }
221  }
222  }
223 
224  if (muon == event->mc_trks.end()) {
225  continue;
226  }
227 
228  if (muon->len == 0.0) {
229  muon->len = gWater(muon->E);
230  }
231 
232 
233  typedef JRange<double> JRange_t;
234 
235  JRange_t za(JRange_t::DEFAULT_RANGE); // Monte Carlo
236  JRange_t zb(JRange_t::DEFAULT_RANGE); // original result
237  JRange_t zd(JRange_t::DEFAULT_RANGE); // new
238 
239  {
240  const JRotation3D R(getDirection(*muon));
241 
242  JTrack3D ta = getTrack(*muon);
243 
244  ta.add(center);
245  ta.rotate(R);
246 
247  za.setLowerLimit(ta.getZ());
248  za.setLength(fabs(muon->len));
249  }
250 
251 
252  buffer_type dataL0;
253 
254  buildL0(*tev, moduleRouter, true, back_inserter(dataL0));
255 
256 
257  if (!evt->empty()) {
258 
259  sort(evt->begin(), evt->end(), qualitySorter);
260 
261  JEvt::const_iterator track = evt->begin();
262 
263  const JRotation3D R (getDirection(*track));
264  const JLine1Z tz(getPosition (*track).rotate(R), track->getT());
265  const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, T_ns);
266 
267  if (track->has(JMUONSTART)) {
268  if (track->hasW(JSTART_LENGTH_METRES) && track->getW(JSTART_LENGTH_METRES) > 0.0) {
269  zb = JRange_t(tz.getZ(), tz.getZ() + track->getW(JSTART_LENGTH_METRES));
270  }
271  }
272 
273  {
275 
276  for (buffer_type::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
277 
278  hit_type hit(*i);
279 
280  hit.rotate(R);
281 
282  if (match(hit)) {
283  top[hit.getModuleID()].insert(hit.getPMTAddress());
284  }
285  }
286 
287  struct JHit_t {
288 
289  double getZ() const { return z; }
290  double getP() const { return p; }
291 
292  double z;
293  double p;
294  };
295 
296  vector<JHit_t> data;
297 
298  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
299 
300  if (summary.hasSummaryFrame(module->getID())) {
301 
302  const JDAQSummaryFrame& frame = summary.getSummaryFrame(module->getID());
303 
304  JPosition3D pos(module->getPosition());
305 
306  pos.transform(R, tz.getPosition());
307 
308  if (pos.getX() <= parameters.roadWidth_m) {
309 
310  const double z = pos.getZ() - pos.getX() / getTanThetaC();
311  const double p = getProbability(*module, frame, rates_Hz.getMultiplesRates(), T_ns.getLength(), top[module->getID()]);
312 
313  data.push_back({ z, p });
314  }
315  }
316  }
317 
318  if (!data.empty()) {
319 
320  sort(data.begin(), data.end(), make_comparator(&JHit_t::getZ));
321 
322  vector<JHit_t>::const_iterator track_start = start.find(data. begin(), data. end());
323  vector<JHit_t>::const_reverse_iterator track_end = start.find(data.rbegin(), data.rend());
324 
325  if (track_start == data. end()) { track_start = data. begin(); }
326  if (track_end == data.rend()) { track_end = data.rbegin(); }
327 
328  zd = JRange_t(track_start->getZ(), track_end->getZ());
329 
330  zd.add(tz.getZ());
331  }
332  }
333  }
334 
335  H1["B"]->Fill(zb.getLength() - za.getLength());
336  H1["D"]->Fill(zd.getLength() - za.getLength());
337 
338  H2["B"]->Fill(zb.getLowerLimit() - za.getLowerLimit(), zb.getUpperLimit() - za.getUpperLimit());
339  H2["D"]->Fill(zd.getLowerLimit() - za.getLowerLimit(), zd.getUpperLimit() - za.getUpperLimit());
340  }
341  STATUS(endl);
342 
343 
344  TFile out(outputFile.c_str(), "recreate");
345 
346  out << H1 << H2;
347 
348  out.Write();
349  out.Close();
350 }
static const int JMUONSTART
Router for direct addressing of PMT data in detector data structure.
Definition: JPMTRouter.hh:35
Utility class to parse command line options.
Definition: JParser.hh:1517
General exception.
Definition: JException.hh:23
JComparator< JResult_t T::*, JComparison::lt > make_comparator(JResult_t T::*member)
Helper method to create comparator between values of data member.
JTrack3E getTrack(const Trk &track)
Get track.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
Router for direct addressing of module data in detector data structure.
Rotation matrix.
Definition: JRotation3D.hh:111
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
*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
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:83
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
const JK40Rates & getK40Rates()
Get K40 rates.
Definition: KM3NeT.hh:36
Simple data structure for histogram binning.
string outputFile
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.
JTime & add(const JTime &value)
Addition operator.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Data structure for fit parameters.
Type definition of range.
Definition: JHead.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
JAxis3D & rotate(const JRotation3D &R)
Rotate axis.
Definition: JAxis3D.hh:225
Detector file.
Definition: JHead.hh:226
JDirection3D getDirection(const Vec &dir)
Get direction.
Acoustic event fit.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
Auxiliary class for start or end point evaluation.
Definition: JStart.hh:21
Data storage class for rate measurements of all PMTs in one module.
JPosition3D getPosition(const Vec &pos)
Get position.
#define NOTICE(A)
Definition: JMessage.hh:64
File router for fast addressing of summary data.
Monte Carlo run header.
Definition: JHead.hh:1167
double getP(const double expval, bool hit)
Get Poisson probability to observe a hit or not for given expectation value for the number of hits...
Definition: JFitToolkit.hh:41
#define FATAL(A)
Definition: JMessage.hh:67
then JCookie sh JDataQuality D $DETECTOR_ID R
Definition: JDataQuality.sh:41
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.
Data structure for L0 hit.
Definition: JHitL0.hh:27
Auxiliary class to set-up Hit.
Definition: JSirene.hh:57
double getProbability(const size_t N, const size_t M, const JK40Rates &R_Hz, const double T_ns)
Get probability due to random background.
Definition: JFitToolkit.hh:248
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
Template L0 hit builder.
Definition: JBuildL0.hh:35
do set_variable DETECTOR_TXT $WORKDIR detector
JPosition3D & rotate(const JRotation3D &R)
Rotate.
Definition: JPosition3D.hh:186
double getZ() const
Get z position.
Definition: JVector3D.hh:115
int debug
debug level
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62