Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JCompareDetector.cc File Reference

Auxiliary program to find differences between two detector files. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JGeometry3D/JQuaternion3D.hh"
#include "JTools/JRange.hh"
#include "JTools/JConstants.hh"
#include "JTools/JQuantile.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.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

Auxiliary program to find differences between two detector files.


A ROOT outout file with histograms is optionally produced.

Author
mjongen

Definition in file JCompareDetector.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 56 of file JCompareDetector.cc.

57 {
58  using namespace std;
59  using namespace JPP;
60 
61  string detectorFile_a;
62  string detectorFile_b;
63  string outputFile;
64  double precision;
65  int debug;
66 
67  try {
68 
69  JParser<> zap("Auxiliary program to find differences between two detector files.");
70 
71  zap['a'] = make_field(detectorFile_a);
72  zap['b'] = make_field(detectorFile_b);
73  zap['o'] = make_field(outputFile) = "";
74  zap['p'] = make_field(precision) = 0.001; // [ns,m]
75  zap['d'] = make_field(debug) = 3;
76 
77  zap(argc, argv);
78  }
79  catch(const exception &error) {
80  FATAL(error.what() << endl);
81  }
82 
83 
84  JDetector detector_a;
85  JDetector detector_b;
86 
87  try {
88  load(detectorFile_a, detector_a);
89  }
90  catch(const JException& error) {
91  FATAL(error);
92  }
93 
94  try {
95  load(detectorFile_b, detector_b);
96  }
97  catch(const JException& error) {
98  FATAL(error);
99  }
100 
101 
102  bool is_equal = true;
103 
104  const JModuleRouter module_router_a(detector_a);
105  const JModuleRouter module_router_b(detector_b);
106 
107 
108  // compare detector IDs
109 
110  if (detector_a.getID() != detector_b.getID()) {
111 
112  DEBUG("* Different detector identifiers "
113  << setw(5) << detector_a.getID() << " (A) and " << endl
114  << setw(5) << detector_b.getID() << " (B)." << endl
115  << endl);
116 
117  is_equal = false;
118  }
119 
120 
121  // check whether the same modules are present in the file
122 
123  for (JDetector::iterator module = detector_a.begin(); module != detector_a.end(); ++module) {
124 
125  if (!module_router_b.hasModule(module->getID())) {
126 
127  DEBUG("* Module " << setw(10) << module->getID() << " is in A " << getLabel(*module) << " but not in B" << endl);
128 
129  is_equal = false;
130  }
131  }
132 
133  for (JDetector::iterator module = detector_b.begin(); module != detector_b.end(); ++module) {
134 
135  if (!module_router_a.hasModule(module->getID())) {
136 
137  DEBUG("* Module " << setw(10) << module->getID() << " is in B " << getLabel(*module) << " but not in A" << endl);
138 
139  is_equal = false;
140  }
141  }
142  DEBUG(endl);
143 
144 
145  // compare the modules that the files have in common
146 
147  DEBUG("Comparing module by module." << endl);
148 
149  for (JDetector::iterator module_a = detector_a.begin(); module_a != detector_a.end(); ++module_a) {
150 
151  if (!module_router_b.hasModule(module_a->getID())) {
152 
153  continue;
154 
155  is_equal = false;
156  }
157 
158  const JModule* module_b = &module_router_b.getModule(module_a->getID());
159 
160  DEBUG(" Module " << setw(10) << module_a->getID());
161 
162  // string and floor numbers
163 
164  if (module_a->getLocation() == module_b->getLocation()) {
165 
166  DEBUG(" " << getLabel(*module_a) << endl);
167 
168  } else {
169 
170  DEBUG(endl);
171  DEBUG(" * different location: "
172  << getLabel(*module_a) << " (A), "
173  << getLabel(*module_b) << " (B)" << endl);
174 
175  is_equal = false;
176  }
177 
178  // average DOM positions
179 
180  if (module_a->getPosition().getDistance(module_b->getPosition()) > precision) {
181 
182  DEBUG(" * different position: "
183  << module_a->getPosition() << " (A), "
184  << module_b->getPosition() << " (B)"
185  << ", B - A " << JPosition3D(module_b->getPosition()-module_a->getPosition()) << endl);
186 
187  is_equal = false;
188  }
189 
190  // number of PMTs
191 
192  if (module_a->size() != module_b->size()) {
193 
194  DEBUG(" * different number of PMTs: "
195  << module_a->size() << " (A), "
196  << module_b->size() << " (B)" << endl);
197 
198  is_equal = false;
199  }
200 
201  // average t0
202 
203  const JQuantile q = compareT0(*module_a, *module_b);
204 
205  if (fabs(q.getMean()) > precision) {
206 
207  DEBUG(" * different average t0: "
208  << ", B - A " << q.getMean()
209  << endl);
210 
211  is_equal = false;
212  }
213 
214  // comparing by PMT
215 
216  // t0
217 
218  for (unsigned int pmt = 0; pmt != module_a->size() && pmt != module_b->size(); ++pmt) {
219 
220  const JPMT& pmt_a = module_a->getPMT(pmt);
221  const JPMT& pmt_b = module_b->getPMT(pmt);
222 
223  if (fabs(pmt_a.getT0() - pmt_b.getT0()) > precision) {
224 
225  DEBUG(" * different T0 PMT " << pmt << ": "
226  << pmt_a.getT0() << " (A" << FILL(2,'0') << pmt << "), "
227  << pmt_b.getT0() << " (B" << FILL(2,'0') << pmt << ")"
228  << ", B - A " << pmt_b.getT0() - pmt_a.getT0()
229  << endl);
230 
231  is_equal = false;
232  }
233  }
234 
235  // positions
236 
237  for (unsigned int pmt = 0; pmt != module_a->size() && pmt != module_b->size(); ++pmt) {
238 
239  const JPMT& pmt_a = module_a->getPMT(pmt);
240  const JPMT& pmt_b = module_b->getPMT(pmt);
241 
242  // PMT positions
243 
244  if (pmt_a.getPosition().getDistance(pmt_b.getPosition()) > precision) {
245 
246  DEBUG(" * different PMT position: "
247  << pmt_a.getPosition() << " (A" << FILL(2,'0') << pmt << "), "
248  << pmt_b.getPosition() << " (B" << FILL(2,'0') << pmt << ")"
249  << ", B - A " << JPosition3D(pmt_b.getPosition() - pmt_a.getPosition()) << endl);
250 
251  is_equal = false;
252  }
253  }
254 
255  // status
256 
257  for (unsigned int pmt = 0; pmt != module_a->size() && pmt != module_b->size(); ++pmt) {
258 
259  const JPMT& pmt_a = module_a->getPMT(pmt);
260  const JPMT& pmt_b = module_b->getPMT(pmt);
261 
262  if (pmt_a.getStatus() != pmt_b.getStatus()) {
263 
264  DEBUG(" * different status PMT " << pmt << ": "
265  << pmt_a.getStatus() << " (A" << FILL(2,'0') << pmt << "), "
266  << pmt_b.getStatus() << " (B" << FILL(2,'0') << pmt << ")"
267  << endl);
268 
269  is_equal = false;
270  }
271  }
272  }
273  DEBUG(endl);
274 
275 
276  if (outputFile != "") {
277 
278  typedef JRange<int> JRange_t;
279 
280  const JRange_t string = combine(JRange_t(detector_a.begin(), detector_a.end(), &JModule::getString),
281  JRange_t(detector_b.begin(), detector_b.end(), &JModule::getString));
282  const JRange_t floor = combine(JRange_t(detector_a.begin(), detector_a.end(), &JModule::getFloor),
283  JRange_t(detector_b.begin(), detector_b.end(), &JModule::getFloor));
284 
285  TFile out(outputFile.c_str(), "recreate");
286 
287  TH2D M2("M2", NULL,
288  string.getLength() + 1,
289  string.getLowerLimit() - 0.5,
290  string.getUpperLimit() + 0.5,
291  floor.getLength() + 1,
292  floor.getLowerLimit() - 0.5,
293  floor.getUpperLimit() + 0.5);
294 
295  TH2D* X2 = (TH2D*) M2.Clone("X2" );
296  TH2D* Y2 = (TH2D*) M2.Clone("Y2" );
297  TH2D* Z2 = (TH2D*) M2.Clone("Z2" );
298  TH2D* T2 = (TH2D*) M2.Clone("T2" );
299  TH2D* RM2 = (TH2D*) M2.Clone("RM2");
300  TH2D* R2 = (TH2D*) M2.Clone("R2" );
301  TH2D* QA = (TH2D*) M2.Clone("QA" );
302  TH2D* QB = (TH2D*) M2.Clone("QB" );
303  TH2D* QC = (TH2D*) M2.Clone("QC" );
304  TH2D* QD = (TH2D*) M2.Clone("QD" );
305 
306  for( JDetector::iterator module = detector_a.begin(); module != detector_a.end(); ++module) {
307  if( !module_router_b.hasModule(module->getID()) ) {
308  M2.Fill(module->getString(), module->getFloor(), -1.0);
309  }
310  }
311 
312  for( JDetector::iterator module = detector_b.begin(); module != detector_b.end(); ++module) {
313  if( !module_router_a.hasModule(module->getID()) ) {
314  M2.Fill(module->getString(), module->getFloor(), +1.0);
315  }
316  }
317 
318  for( JDetector::iterator module_a = detector_a.begin(); module_a != detector_a.end(); ++module_a) {
319 
320  if (!module_router_b.hasModule(module_a->getID())) {
321  continue;
322  }
323 
324  const JModule* module_b = &module_router_b.getModule(module_a->getID());
325 
326  double phi = 0.0;
327 
328  for (unsigned int pmt = 0; pmt != module_a->size() && pmt != module_b->size(); ++pmt) {
329 
330  double x = module_a->getPMT(pmt).getPhi() - module_b->getPMT(pmt).getPhi();
331 
332  while (x > +PI) { x -= PI; }
333  while (x < -PI) { x += PI; }
334 
335  phi += x;
336  }
337 
338  phi /= min(module_a->size(),
339  module_b->size());
340 
341  JQuaternion3D Q = getRotation(*module_a, *module_b);
342 
343  JQuantile q = compareT0(*module_a, *module_b);
344 
345  X2 ->Fill(module_a->getString(), module_a->getFloor(), module_a->getX() - module_b->getX());
346  Y2 ->Fill(module_a->getString(), module_a->getFloor(), module_a->getY() - module_b->getY());
347  Z2 ->Fill(module_a->getString(), module_a->getFloor(), module_a->getZ() - module_b->getZ());
348  T2 ->Fill(module_a->getString(), module_a->getFloor(), q.getMean());
349  RM2->Fill(module_a->getString(), module_a->getFloor(), q.getRMS());
350  R2 ->Fill(module_a->getString(), module_a->getFloor(), phi * 180.0/PI);
351  QA ->Fill(module_a->getString(), module_a->getFloor(), Q.getA());
352  QB ->Fill(module_a->getString(), module_a->getFloor(), Q.getB());
353  QC ->Fill(module_a->getString(), module_a->getFloor(), Q.getC());
354  QD ->Fill(module_a->getString(), module_a->getFloor(), Q.getD());
355  }
356 
357  out.Write();
358  out.Close();
359  }
360 
361  ASSERT(is_equal);
362 
363  return 0;
364 }
Utility class to parse command line options.
Definition: JParser.hh:1493
General exception.
Definition: JException.hh:23
const JStatus & getStatus() const
Get status.
Definition: JStatus.hh:80
Data structure for a composite optical module.
Definition: JModule.hh:50
Quantile calculator.
Definition: JQuantile.hh:83
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition: JLocation.hh:246
Detector data structure.
Definition: JDetector.hh:80
Router for direct addressing of module data in detector data structure.
static const double PI
Constants.
Definition: JConstants.hh:20
static JRotation getRotation
Function object to get rotation matrix to go from first to second module.
double getRMS() const
Get RMS.
Definition: JQuantile.hh:345
double getPhi() const
Get phi angle.
Definition: JVersor3D.hh:141
string outputFile
double getDistance(const JVector3D &pos) const
Get distance to point.
Definition: JVector3D.hh:269
double getMean() const
Get mean value.
Definition: JQuantile.hh:331
esac $JPP_DIR examples JDetector JTransitTime o $OUTPUT_FILE n N $NPE T $TTS_NS d $DEBUG for HISTOGRAM in tts tt2 pmt
Definition: JTransitTime.sh:36
#define ASSERT(A,...)
Assert macro.
Definition: JMessage.hh:90
const JLocation & getLocation() const
Get location.
Definition: JLocation.hh:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
Data structure for PMT geometry and calibration.
Definition: JPMT.hh:47
double getY() const
Get y position.
Definition: JVector3D.hh:103
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:129
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:174
Auxiliary data structure for sequence of same character.
Definition: JPrint.hh:361
JRange< T, JComparator_t > combine(const JRange< T, JComparator_t > &first, const JRange< T, JComparator_t > &second)
Combine ranges.
Definition: JRange.hh:699
#define FATAL(A)
Definition: JMessage.hh:67
Data structure for quaternion in three dimensions.
JRange< Double_t > JRange_t
Definition: JFitToT.hh:34
double getX() const
Get x position.
Definition: JVector3D.hh:93
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
double getZ() const
Get z position.
Definition: JVector3D.hh:114
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
double getT0() const
Get time offset.