Jpp  15.0.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
examples/JDynamics/JBallarat.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <limits>
5 #include <map>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TGraph.h"
11 
12 #include "JROOT/JGraph.hh"
13 #include "JROOT/JRootToolkit.hh"
14 #include "JROOT/JManager.hh"
15 
16 #include "JDetector/JDetector.hh"
20 
22 #include "JCompass/JEvt.hh"
23 #include "JCompass/JSupport.hh"
24 
25 #include "JDynamics/JDynamics.hh"
26 
27 #include "Jeep/JPrint.hh"
28 #include "Jeep/JParser.hh"
29 #include "Jeep/JMessage.hh"
30 
31 
32 namespace {
33 
34  using JUTC::JUTCTimeRange;
35 
36  /**
37  * Auxiliary method to get UTC time range of calibration data.
38  *
39  * \param __begin begin of calibration data
40  * \param __end end of calibration data
41  * \return UTC time range
42  */
43  template<class T>
44  static inline JUTCTimeRange getUTCTimeRange(T __begin, T __end)
45  {
46  JUTCTimeRange range = JUTCTimeRange::DEFAULT_RANGE;
47 
48  for (T i = __begin; i != __end; ++i) {
49  if (!i->second.empty()) {
50  range.combine(JUTCTimeRange(i->second.getXmin(), i->second.getXmax()));
51  }
52  }
53 
54  return range;
55  }
56 }
57 
58 
59 /**
60  * \file
61  *
62  * Example program to plot orientation data.
63  * \author mdejong
64  */
65 int main(int argc, char **argv)
66 {
67  using namespace std;
68  using namespace JPP;
69 
71  JLimit_t& numberOfEvents = inputFile.getLimit();
72  string detectorFile;
73  string outputFile;
74  double Tmax_s;
75  int debug;
76 
77  try {
78 
79  JParser<> zap("Example program to plot orientation data.");
80 
81  zap['f'] = make_field(inputFile, "output of JBallarat[.sh]");
82  zap['n'] = make_field(numberOfEvents) = JLimit::max();
83  zap['a'] = make_field(detectorFile);
84  zap['o'] = make_field(outputFile) = "gold.root";
85  zap['T'] = make_field(Tmax_s);
86  zap['d'] = make_field(debug) = 2;
87 
88  zap(argc, argv);
89  }
90  catch(const exception &error) {
91  FATAL(error.what() << endl);
92  }
93 
94 
96 
97  try {
98  load(detectorFile, detector);
99  }
100  catch(const JException& error) {
101  FATAL(error);
102  }
103 
104  const JModuleRouter router(detector);
105 
106  if (!hasDetectorAddressMap(detector.getID())) {
107  FATAL("No detector address map for detector identier " << detector.getID() << endl);
108  }
109 
110  const JDetectorAddressMap& demo = getDetectorAddressMap(detector.getID());
111 
112 
113  JDynamics dynamics(detector, Tmax_s);
114 
115  STATUS("loading input from file(s) " << inputFile << "... " << flush);
116 
117  dynamics.load(inputFile);
118 
119  STATUS("OK" << endl);
120 
121 
123 
124  for (JDynamics::JOrientation::const_iterator module = dynamics.orientation.begin(); module != dynamics.orientation.end(); ++module) {
125 
126  if (router.getModule(module->first).getFloor() != 0) {
127 
128  STATUS("Creating graphics for input data " << setw(10) << module->first << "... " << flush);
129 
130  JQuaternion3D Q0 = getRotation(getModule(demo.get(module->first)), router.getModule(module->first));
131 
132  Q0.conjugate();
133 
134  JGraph_t& z0 = Z0[module->first];
135 
136  for (JDynamics::JOrientation::function_type::const_iterator i = module->second.begin(); i != module->second.end(); ++i) {
137 
138  const JQuaternion3D Q1 = router.getModule(module->first).getQuaternion() * i->getY();
139  const JQuaternion3D Q = Q1 * Q0;
140 
142 
143  z0.put(i->getX(), (q.twist.getD() > 0.0 ? +2.0 : -2.0) * acos(q.twist.getA()));
144  }
145 
146  STATUS("OK" << endl);
147  }
148  }
149 
150  map<int, TH1D*> H0;
151 
152  for (JDynamics::JOrientation::const_iterator module = dynamics.orientation.begin(); module != dynamics.orientation.end(); ++module) {
153 
154  STATUS("Creating graphics for compiled data " << setw(10) << module->first << "... " << flush);
155 
156  if (router.getModule(module->first).getFloor() != 0 && !module->second.empty()) {
157 
158  JQuaternion3D Q0 = getRotation(getModule(demo.get(module->first)), router.getModule(module->first));
159 
160  Q0.conjugate();
161 
162  const Double_t xmin = module->second.getXmin();
163  const Double_t xmax = module->second.getXmax();
164  const Int_t nx = (Int_t) ((xmax - xmin) / Tmax_s);
165 
166  TH1D* h0 = H0[module->first] = new TH1D(MAKE_CSTRING("H[" << module->first << "].twist"), NULL, nx, xmin, xmax);
167 
168  for (Int_t i = 1; i <= h0->GetXaxis()->GetNbins(); ++i) {
169 
170  const Double_t x = h0->GetXaxis()->GetBinCenter(i);
171  const JQuaternion3D Q = module->second(x) * Q0;
172 
174 
175  h0->SetBinContent(i, (q.twist.getD() > 0.0 ? +2.0 : -2.0) * acos(q.twist.getA()));
176  }
177  }
178 
179  STATUS("OK" << endl);
180  }
181 
182 
183  const Double_t xmin = getUTCTimeRange(dynamics.orientation.begin(), dynamics.orientation.end()).getLowerLimit();
184  const Double_t xmax = getUTCTimeRange(dynamics.orientation.begin(), dynamics.orientation.end()).getUpperLimit();
185  const Int_t nx = (Int_t) ((xmax - xmin) / Tmax_s);
186 
187  JManager<int, TH1D> X0(new TH1D("X[%].twist", NULL, nx, xmin, xmax));
188 
189  for (Int_t i = 1; i <= X0->GetXaxis()->GetNbins(); ++i) {
190 
191  STATUS("Creating graphics for detector data " << FIXED(5,1) << (double) (i * 100) / (double) X0->GetXaxis()->GetNbins() << "%... " << flush << '\r');
192 
193  const Double_t x = X0->GetXaxis()->GetBinCenter(i);
194  const JDetector& buffer = dynamics.update(x);
195 
196  for (JDetector::const_iterator module = buffer.begin(); module != buffer.end(); ++module) {
197 
198  if (module->getFloor() != 0) {
199 
200  TH1D* x0 = X0[module->getID()];
201 
202  const JQuaternion3D Q = getRotation(router.getModule(module->getID()), *module);
203 
205 
206  x0->SetBinContent(i, (q.twist.getD() > 0.0 ? +2.0 : -2.0) * acos(q.twist.getA()));
207  }
208  }
209  }
210  STATUS(endl);
211 
212 
213  TFile out(outputFile.c_str(), "recreate");
214 
215  for (map<int, TH1D*>::const_iterator i = H0.begin(); i != H0.end(); ++i) {
216  out << *i->second;
217  }
218 
219  out << X0;
220 
221  for (map<int, JGraph_t>::const_iterator i = Z0.begin(); i != Z0.end(); ++i) {
222  out << JGraph(i->second, MAKE_CSTRING("G[" << i->first << "].twist"));
223  }
224 
225  out.Write();
226  out.Close();
227 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
JDAQUTCTimeRange getUTCTimeRange()
Get UTC time range.
Q(UTCMax_s-UTCMin_s)-livetime_s
int main(int argc, char *argv[])
Definition: Main.cc:15
void put(const Double_t x, const Double_t y)
Put data.
Definition: JGraph.hh:28
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for graph data.
Definition: JGraph.hh:21
collection_type::const_iterator const_iterator
Definition: JPolfit.hh:157
#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.
bool hasDetectorAddressMap(const int id)
Check if detector address map is available.
static JRotation getRotation
Function object to get rotation matrix to go from first to second module.
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
Lookup table for PMT addresses in detector.
range_type & combine(const range_type &range)
Combine ranges.
Definition: JRange.hh:433
Dynamic ROOT object management.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
Auxiliary data structure to build TGraph.
Definition: JGraph.hh:42
string outputFile
Data structure for detector geometry and calibration.
ROOT TTree parameter settings.
const JQuaternion3D & getQuaternion() const
Get quaternion.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Detector specific mapping between logical positions and readout channels of PMTs in optical modules...
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:196
JQuaternion3D twist
rotation around parallel axis
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
do set_variable OUTPUT_DIRECTORY $WORKDIR T
Dynamic detector calibration.
int debug
debug level
Definition: JSirene.cc:63
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
General purpose messaging.
data_type::const_iterator const_iterator
Definition: JDynamics.hh:123
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
double getD() const
Get d value.
Data structure for unit quaternion in three dimensions.
Direct access to module in detector data structure.
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
UTC time range [s].
Dynamic detector calibration.
Definition: JDynamics.hh:61
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
static const JVector3D JVector3Z_t(0, 0, 1)
unit z-vector
double getA() const
Get a value.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Compass event data types.
do set_variable DETECTOR_TXT $WORKDIR detector
Auxiliary data structure for decomposition of quaternion in twist and swing quaternions.
const JModule & getModule(const JDetector &detector, const JModuleLocation &location)
find module with a given string and floor number
then usage $script< detector file >< tripodfile >< stage > input file nInput files correspond to the output of JAcousticsEventBuilder[.sh] nFirst stage eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY eval JPrintDetector a $DETECTOR O CAN source JAcoustics sh $DETECTOR_ID typeset A CONFIGURATION for key in Tmax_s
JQuaternion3D & conjugate()
Conjugate quaternion.