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

Example program to plot orientation data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TGraph.h"
#include "JROOT/JGraph.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JManager.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDetectorAddressMapToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JCompass/JEvt.hh"
#include "JCompass/JSupport.hh"
#include "JDynamics/JDynamics.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

Example program to plot orientation data.

Author
mdejong

Definition in file examples/JDynamics/JBallarat.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 86 of file examples/JDynamics/JBallarat.cc.

87 {
88  using namespace std;
89  using namespace JPP;
90 
92  JLimit_t& numberOfEvents = inputFile.getLimit();
93  string detectorFile;
94  string outputFile;
95  double Tmax_s;
96  int debug;
97 
98  try {
99 
100  JParser<> zap("Example program to plot orientation data.");
101 
102  zap['f'] = make_field(inputFile, "output of JBallarat[.sh]");
103  zap['n'] = make_field(numberOfEvents) = JLimit::max();
104  zap['a'] = make_field(detectorFile);
105  zap['o'] = make_field(outputFile) = "gold.root";
106  zap['T'] = make_field(Tmax_s);
107  zap['d'] = make_field(debug) = 2;
108 
109  zap(argc, argv);
110  }
111  catch(const exception &error) {
112  FATAL(error.what() << endl);
113  }
114 
115 
117 
118  try {
119  load(detectorFile, detector);
120  }
121  catch(const JException& error) {
122  FATAL(error);
123  }
124 
125  const JModuleRouter router(detector);
126 
127  if (!hasDetectorAddressMap(detector.getID())) {
128  FATAL("No detector address map for detector identier " << detector.getID() << endl);
129  }
130 
131  const JDetectorAddressMap& demo = getDetectorAddressMap(detector.getID());
132 
133 
134  JDynamics dynamics(detector, Tmax_s);
135 
136  STATUS("loading input from file(s) " << inputFile << "... " << flush);
137 
138  dynamics.load(inputFile);
139 
140  STATUS("OK" << endl);
141 
142 
144 
145  for (JDynamics::JOrientation::const_iterator module = dynamics.orientation.begin(); module != dynamics.orientation.end(); ++module) {
146 
147  if (router.getModule(module->first).getFloor() != 0) {
148 
149  STATUS("Creating graphics for input data " << setw(10) << module->first << flush << '\r');
150 
151  JQuaternion3D Q0 = getRotation(getModule(demo.get(module->first)), router.getModule(module->first));
152 
153  Q0.conjugate();
154 
155  JGraph_t& z0 = Z0[module->first];
156 
157  for (JDynamics::JOrientation::function_type::const_iterator i = module->second.begin(); i != module->second.end(); ++i) {
158 
159  const JQuaternion3D Q1 = router.getModule(module->first).getQuaternion() * i->getY();
160 
161  z0.put(i->getX(), getTwist(Q1 * Q0));
162  }
163  }
164  }
165  STATUS(endl);
166 
167  map<int, TH1D*> H1;
168 
169  for (JDynamics::JOrientation::const_iterator module = dynamics.orientation.begin(); module != dynamics.orientation.end(); ++module) {
170 
171  if (router.getModule(module->first).getFloor() != 0 && module->second.size() > 1u) {
172 
173  STATUS("Creating graphics for input data " << setw(10) << module->first << flush << '\r');
174 
175  TH1D* h1 = H1[module->first] = new TH1D(MAKE_CSTRING("U[" << module->first << "].twist"), NULL, 500, -180.0, +180.0);
176 
177  for (JDynamics::JOrientation::function_type::const_iterator i1 = module->second.begin(), i0 = i1++; i1 != module->second.end(); ++i0, ++i1) {
178 
179  JQuaternion3D Q0 = i0->getY();
180  JQuaternion3D Q1 = i1->getY();
181 
182  h1->Fill(getTwist(Q1 * Q0.conjugate()) * 180.0 / PI);
183  }
184  }
185  }
186  STATUS(endl);
187 
188  map<int, TH1D*> H0;
189 
190  for (JDynamics::JOrientation::const_iterator module = dynamics.orientation.begin(); module != dynamics.orientation.end(); ++module) {
191 
192  STATUS("Creating graphics for compiled data " << setw(10) << module->first << flush << '\r');
193 
194  if (router.getModule(module->first).getFloor() != 0 && !module->second.empty()) {
195 
196  JQuaternion3D Q0 = getRotation(getModule(demo.get(module->first)), router.getModule(module->first));
197 
198  Q0.conjugate();
199 
200  const Double_t xmin = module->second.getXmin();
201  const Double_t xmax = module->second.getXmax();
202  const Int_t nx = (Int_t) ((xmax - xmin) / Tmax_s);
203 
204  TH1D* h0 = H0[module->first] = new TH1D(MAKE_CSTRING("H[" << module->first << "].twist"), NULL, nx, xmin, xmax);
205 
206  for (Int_t i = 1; i <= h0->GetXaxis()->GetNbins(); ++i) {
207 
208  const Double_t x = h0->GetXaxis()->GetBinCenter(i);
209  const JQuaternion3D Q = module->second(x) * Q0;
210 
211  h0->SetBinContent(i, getTwist(Q));
212  }
213  }
214  }
215  STATUS(endl);
216 
217  const Double_t xmin = getUTCTimeRange(dynamics.orientation.begin(), dynamics.orientation.end()).getLowerLimit();
218  const Double_t xmax = getUTCTimeRange(dynamics.orientation.begin(), dynamics.orientation.end()).getUpperLimit();
219  const Int_t nx = (Int_t) ((xmax - xmin) / Tmax_s);
220 
221  JManager<int, TH1D> X0(new TH1D("X[%].twist", NULL, nx, xmin, xmax));
222 
223  for (Int_t i = 1; i <= X0->GetXaxis()->GetNbins(); ++i) {
224 
225  STATUS("Creating graphics for detector data " << FIXED(5,1) << (double) (i * 100) / (double) X0->GetXaxis()->GetNbins() << "%" << flush << '\r');
226 
227  const Double_t x = X0->GetXaxis()->GetBinCenter(i);
228  const JDetector& buffer = dynamics.update(x);
229 
230  for (JDetector::const_iterator module = buffer.begin(); module != buffer.end(); ++module) {
231 
232  if (module->getFloor() != 0) {
233 
234  TH1D* x0 = X0[module->getID()];
235 
236  const JQuaternion3D Q = getRotation(router.getModule(module->getID()), *module);
237 
238  x0->SetBinContent(i, getTwist(Q));
239  }
240  }
241  }
242  STATUS(endl);
243 
244 
245  TFile out(outputFile.c_str(), "recreate");
246 
247  for (map<int, TH1D*>::const_iterator i = H0.begin(); i != H0.end(); ++i) {
248  out << *i->second;
249  }
250 
251  for (map<int, TH1D*>::const_iterator i = H1.begin(); i != H1.end(); ++i) {
252  out << *i->second;
253  }
254 
255  out << X0;
256 
257  for (map<int, JGraph_t>::const_iterator i = Z0.begin(); i != Z0.end(); ++i) {
258  out << JGraph(i->second, MAKE_CSTRING("G[" << i->first << "].twist"));
259  }
260 
261  out.Write();
262  out.Close();
263 }
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
void put(const Double_t x, const Double_t y)
Put data.
Definition: JGraph.hh:28
Data structure for graph data.
Definition: JGraph.hh:21
collection_type::const_iterator const_iterator
Definition: JPolfit.hh:162
#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
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
Lookup table for PMT addresses in detector.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
Auxiliary data structure to build TGraph.
Definition: JGraph.hh:42
string outputFile
const JQuaternion3D & getQuaternion() const
Get quaternion.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:224
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
static const double PI
Mathematical constants.
int debug
debug level
Definition: JSirene.cc:63
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
data_type::const_iterator const_iterator
Definition: JDynamics.hh:124
#define FATAL(A)
Definition: JMessage.hh:67
Data structure for unit quaternion in three dimensions.
Dynamic detector calibration.
Definition: JDynamics.hh:62
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.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do set_variable DETECTOR_TXT $WORKDIR detector
double u[N+1]
Definition: JPolint.hh:755
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.