Jpp - the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JMuonStart.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <vector>
6 #include <set>
7 #include <algorithm>
8 #include <memory>
9 
13 #include "JDAQ/JDAQEventIO.hh"
14 #include "JDAQ/JDAQTimesliceIO.hh"
16 
17 #include "JDetector/JDetector.hh"
21 
22 #include "JPhysics/JK40Rates.hh"
23 
24 #include "JTrigger/JHit.hh"
25 #include "JTrigger/JHitR1.hh"
26 #include "JTrigger/JBuildL0.hh"
28 
32 #include "JSupport/JSupport.hh"
33 #include "JSupport/JMeta.hh"
34 
35 #include "JFit/JLine1Z.hh"
36 #include "JFit/JEnergyRegressor.hh"
37 #include "JFit/JFitToolkit.hh"
38 #include "JFit/JNPEHit.hh"
39 #include "JFit/JMEstimator.hh"
40 #include "JFit/JModel.hh"
41 
42 #include "JReconstruction/JEvt.hh"
46 
47 #include "JPhysics/JConstants.hh"
48 #include "JTools/JFunction1D_t.hh"
50 
51 #include "JPhysics/JPDFTable.hh"
52 #include "JPhysics/JPDFToolkit.hh"
53 #include "JPhysics/JNPETable.hh"
54 
55 #include "JLang/JComparator.hh"
56 
57 #include "Jeep/JPrint.hh"
58 #include "Jeep/JParser.hh"
59 #include "Jeep/JMessage.hh"
60 
61 
62 /**
63  * \file
64  *
65  * Program to determine start position of muon trajectory.
66  * \author mdejong
67  */
68 int main(int argc, char **argv)
69 {
70  using namespace std;
71  using namespace JPP;
72  using namespace KM3NETDAQ;
73 
74  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
75  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
77 
78  JParallelFileScanner_t inputFile;
80  JLimit_t& numberOfEvents = inputFile.getLimit();
81  string detectorFile;
82  string pdfFile;
84  int debug;
85 
86  try {
87 
88  JParser<> zap("Program to perform PDF fit of muon trajectory to data.");
89 
90  zap['f'] = make_field(inputFile);
91  zap['o'] = make_field(outputFile) = "start.root";
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 
105  cout.tie(&cerr);
106 
107 
109 
110  try {
111  load(detectorFile, detector);
112  }
113  catch(const JException& error) {
114  FATAL(error);
115  }
116 
117  const JModuleRouter router(detector);
118 
119  typedef JHitR1 hit_type;
120 
121  typedef vector<hit_type> buffer_type;
122  const JBuildL0<hit_type> buildL0;
123 
124 
125  typedef JRegressor<JEnergy> JRegressor_t;
126 
128  JRegressor_t::T_ns.setRange(parameters.TMin_ns, parameters.TMax_ns);
129 
130  const JStart start(parameters.Pmin1, parameters.Pmin2);
131 
132  JRegressor_t fit(pdfFile);
133 
134  if (fit.getRmax() < parameters.roadWidth_m) {
135  parameters.roadWidth_m = fit.getRmax();
136  }
137 
138 
139  outputFile.open();
140 
141  outputFile.put(JMeta(argc, argv));
142 
143  while (inputFile.hasNext()) {
144 
145  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
146 
147  multi_pointer_type ps = inputFile.next();
148 
149  const JDAQEvent* tev = ps;
150  const JEvt* in = ps;
151 
152  // select start values
153 
154  JEvt cp = *in;
155  JEvt out;
156 
157  if (parameters.reprocess) {
159  }
160 
161  cp.select(parameters.numberOfPrefits, qualitySorter);
162 
163  if (!cp.empty()) {
164  cp.select(JHistory::is_event(cp.begin()->getHistory()));
165  }
166 
167 
168  buffer_type data;
169 
170  buildL0(*tev, router, true, back_inserter(data));
171 
172 
173  for (JEvt::iterator track = cp.begin(); track != cp.end(); ++track) {
174 
175  const JRotation3D R (getDirection(*track));
176  const JLine1Z tz(getPosition (*track).rotate(R), track->getT());
177  const JModel<JLine1Z> match(tz, parameters.roadWidth_m, JRegressor_t::T_ns);
178 
179 
181 
182  for (buffer_type::const_iterator i = data.begin(); i != data.end(); ++i) {
183 
184  hit_type hit(*i);
185 
186  hit.rotate(R);
187 
188  if (match(hit)) {
189  top.insert(hit.getModuleIdentifier());
190  }
191  }
192 
193  vector<JNPEHit> data;
194 
195  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
196 
197  const size_t count = top.count(i->getID());
198 
199  if (count != 0) {
200 
201  JPosition3D pos(i->getPosition());
202 
203  pos.transform(R, tz.getPosition());
204 
205  if (pos.getX() <= parameters.roadWidth_m) {
206 
207  JModule module = *i;
208 
209  module.transform(R, tz.getPosition());
210 
211  data.push_back(JNPEHit(fit.getNPE(module.begin(), module.end(), parameters.R_Hz), count));
212  }
213  }
214  }
215 
216 
217  double Zmin = 0.0; // relative start position [m]
218  double Zmax = 0.0; // relative end position [m]
219  double npe = 0.0; // npe due to mip up to barycenter of hits
220  double npe_total = 0.0; // npe due to mip up to end of track
221 
222 
223  if (!data.empty()) {
224 
225  sort(data.begin(), data.end(), make_comparator(&JNPEHit::getZ));
226 
227  vector<JNPEHit>::const_iterator track_start = start.find(data. begin(), data. end());
228  vector<JNPEHit>::const_reverse_iterator track_end = start.find(data.rbegin(), data.rend());
229 
230  if (track_start == data. end()) { track_start = data. begin(); }
231  if (track_end == data.rend()) { track_end = data.rbegin(); }
232 
233  Zmin = track_start->getZ();
234  Zmax = track_end ->getZ();
235 
236  for (vector<JNPEHit>::const_iterator i = track_start; i != track_end.base(); ++i) {
237 
238  if (i->getZ() <= 0.0) {
239  npe += i->getYA();
240  }
241 
242  npe_total += i->getYA();
243  }
244 
245  DEBUG("vertex " << FIXED(6,1) << track->getW(JCOPY_Z_M, 0.0) << endl);
246 
247  for (vector<JNPEHit>::const_iterator i = data.begin(); i != data.end(); ++i) {
248  DEBUG(FIXED(6,1) << i->getZ() << ' ' <<
249  FIXED(6,4) << i->getY0() << ' ' <<
250  FIXED(6,4) << i->getYA() << ' ' <<
251  setw(1) << i->getN() << ' ' <<
252  FIXED(6,4) << i->getP() << endl);
253  }
254 
255  DEBUG(" --> " << FIXED(6,1) << Zmin << ' ' << FIXED(6,1) << Zmax << endl);
256  }
257 
258  // move track
259 
260  track->move(Zmin, getSpeedOfLight());
261 
262  out.push_back(JFit(*track).add(JMUONSTART));
263 
264  out.rbegin()->setW(track->getW());
265  out.rbegin()->setW(JSTART_NPE_MIP, npe);
266  out.rbegin()->setW(JSTART_NPE_MIP_TOTAL, npe_total);
267  out.rbegin()->setW(JSTART_LENGTH_METRES, Zmax - Zmin);
268  }
269 
270  // apply default sorter
271 
272  sort(out.begin(), out.end(), qualitySorter);
273 
274  copy(in->begin(), in->end(), back_inserter(out));
275 
276  outputFile.put(out);
277  }
278  STATUS(endl);
279 
281 
282  io >> outputFile;
283 
284  outputFile.close();
285 }
Auxiliary methods to evaluate Poisson probabilities and chi2.
static const int JMUONSTART
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
static const int JSTART_NPE_MIP_TOTAL
number of photo-electrons along the whole track from JStart.cc
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
Auxiliary methods for PDF calculations.
JComparator< JResult_t T::*, JComparison::lt > make_comparator(JResult_t T::*member)
Helper method to create comparator between values of data member.
Definition: JComparator.hh:185
Data structure for a composite optical module.
Definition: JModule.hh:57
Auxiliary method to locate start and end point of muon trajectory.
Regressor function object for JEnergy fit.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
Recording of objects on file according a format that follows from the file name extension.
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
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
T find(T __begin, T __end) const
Get start or end point of muon trajectory.
Definition: JStart.hh:108
Basic data structure for time and time over threshold information of hit.
JFit & add(const int type)
Add event to history.
string outputFile
Data structure for detector geometry and calibration.
void transform(const JRotation3D &R, const JVector3D &pos)
Transform position.
Definition: JPosition3D.hh:331
Various implementations of functional maps.
static const int JSTART_NPE_MIP
number of photo-electrons up to the barycentre from JStart.cc
Type list.
Definition: JTypeList.hh:22
Scanning of objects from a single file according a format that follows from the extension of each fil...
Auxiliary class to extract a subset of optical modules from a detector.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
static const int JCOPY_Z_M
true vertex position along track [m] from JCopy.cc
Data structure for fit parameters.
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:196
JDirection3D getDirection(const Vec &dir)
Get direction.
Acoustic event fit.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Auxiliary class for start or end point evaluation.
Definition: JStart.hh:21
Auxiliary class to test history.
Definition: JHistory.hh:104
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
ROOT I/O of application specific meta data.
JPosition3D getPosition(const Vec &pos)
Get position.
Physics constants.
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
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
Direct access to module in detector data structure.
Reduced data structure for L1 hit.
const double getSpeedOfLight()
Get speed of light.
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
Utility class to parse command line options.
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:382
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
Object reading from a list of files.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Data regression method for JFIT::JEnergy.
Template L0 hit builder.
Definition: JBuildL0.hh:35
do set_variable DETECTOR_TXT $WORKDIR detector
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
Maximum likelihood estimator (M-estimators).
bool qualitySorter(const JRECONSTRUCTION::JFit &first, const JRECONSTRUCTION::JFit &second)
Comparison of fit results.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62