Jpp  16.0.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JNarrabri.cc File Reference

Example program to check contents of acoustic events. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <fstream>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.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 "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JMeta.hh"
#include "JSupport/JSupport.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/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 check contents of acoustic events.

Author
mdejong

Definition in file JNarrabri.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 39 of file JNarrabri.cc.

40 {
41  using namespace std;
42  using namespace JPP;
43 
45  JLimit_t& numberOfEvents = inputFile.getLimit();
46  string detectorFile;
47  string outputFile;
48  double Tmax_s;
49  int debug;
50 
51  try {
52 
53  JParser<> zap("Example program to check contents of acoustic events.");
54 
55  zap['f'] = make_field(inputFile, "output of JKatoomba[.sh]");
56  zap['n'] = make_field(numberOfEvents) = JLimit_t::max();
57  zap['a'] = make_field(detectorFile);
58  zap['o'] = make_field(outputFile) = "narrabri.root";
59  zap['T'] = make_field(Tmax_s);
60  zap['d'] = make_field(debug) = 1;
61 
62  zap(argc, argv);
63  }
64  catch(const exception &error) {
65  FATAL(error.what() << endl);
66  }
67 
68 
70 
71  try {
72  load(detectorFile, detector);
73  }
74  catch(const JException& error) {
75  FATAL(error);
76  }
77 
78 
79  JDynamics dynamics(detector, Tmax_s);
80 
81  STATUS("loading input from file(s) " << inputFile << "... " << flush);
82 
83  dynamics.load(inputFile);
84 
85  STATUS("OK" << endl);
86 
87 
90 
91  Double_t xmin = numeric_limits<Double_t>::max();
92  Double_t xmax = numeric_limits<Double_t>::lowest();
93 
94  for (JDynamics::JPosition::const_iterator string = dynamics.position.begin(); string != dynamics.position.end(); ++string) {
95 
96  if (!string->second.empty()) {
97  xmin = min(xmin, string->second.getXmin());
98  xmax = max(xmax, string->second.getXmax());
99  }
100  }
101 
102  const JFormat_t format(4, 0, std::ios_base::fmtflags(), '0');
103 
104  JManager<int, TH2D> H2(new TH2D ("[%].tilt", NULL, 300, -3.0, +3.0, 300, -3.0, +3.0), '%', format);
105  JManager<int, TH1D> HT(new TH1D ("[%].time", NULL, 200, 0.0, 4.0), '%', format);
106 
107  JManager<int, TH1D> HO(new TH1D ("H[%].orientation", NULL, 1000, xmin, xmax), '%', format);
108  JManager<int, TH1D> HA(new TH1D ("H[%].amplitude", NULL, 1000, xmin, xmax), '%', format);
109 
110  for (JDynamics::JPosition::const_iterator string = dynamics.position.begin(); string != dynamics.position.end(); ++string) {
111 
112  TH1D* ho = HO[string->first];
113  TH1D* ha = HA[string->first];
114 
115  for (Int_t i = 1; i <= HO->GetXaxis()->GetNbins(); ++i) {
116 
117  const Double_t x = HO->GetXaxis()->GetBinCenter(i);
118 
119  try {
120 
121  const JMODEL::JString tilt = string->second(x);
122 
123  ho->SetBinContent(i, tilt.getAngle());
124  ha->SetBinContent(i, tilt.getLength());
125  }
126  catch(const exception& error) {
127  ERROR("Error string " << setw(4) << string->first << " at " << FIXED(20,5) << x << ' ' << error.what() << endl);
128  }
129  }
130  }
131 
132 
133  for (JDynamics::JPosition::const_iterator string = dynamics.position.begin(); string != dynamics.position.end(); ++string) {
134 
135  if (string->second.size() > 1) {
136 
137  TH2D* h2 = H2[string->first];
138  TH1D* hc = HT[string->first];
139 
140  for (typename JDynamics::JPosition::function_type::const_iterator q = string->second.begin(), p = q++; q != string->second.end(); ++p, ++q) {
141 
142  const double t1 = q->getX() - p->getX();
143 
144  hc->Fill(log10(t1));
145 
146  if (t1 > 0 && t1 < Tmax_s) {
147  h2->Fill(600.0e3 * (q->getY().tx - p->getY().tx) / t1, // [mrad/10 min]
148  600.0e3 * (q->getY().ty - p->getY().ty) / t1); // [mrad/10 min]
149  }
150  }
151 
152  for (typename JDynamics::JPosition::function_type::const_iterator i = string->second.begin(); i != string->second.end(); ++i) {
153  ZO[string->first].put(i->getX(), i->getY().getAngle());
154  ZA[string->first].put(i->getX(), i->getY().getLength());
155  }
156  }
157  }
158 
159 
160  TFile out(outputFile.c_str(), "recreate");
161 
162  out << H2 << HT << HO << HA;
163 
164  for (map<int, JGraph_t>::const_iterator i = ZO.begin(); i != ZO.end(); ++i) {
165  out << JGraph(i->second, MAKE_CSTRING("G[" << i->first << "].orientation"));
166  }
167 
168  for (map<int, JGraph_t>::const_iterator i = ZA.begin(); i != ZA.end(); ++i) {
169  out << JGraph(i->second, MAKE_CSTRING("G[" << i->first << "].amplitude"));
170  }
171 
172  out.Write();
173  out.Close();
174 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
double getAngle() const
Get angle.
collection_type::const_iterator const_iterator
Definition: JPolfit.hh:162
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
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_type::const_iterator const_iterator
Definition: JDynamics.hh:299
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
set_variable E_E log10(E_{fit}/E_{#mu})"
#define ERROR(A)
Definition: JMessage.hh:66
int debug
debug level
Definition: JSirene.cc:63
#define FATAL(A)
Definition: JMessage.hh:67
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
Data structure for format specifications.
Definition: JManip.hh:522
double getLength() const
Get length.
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