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

Example program to monitor acoustic events. More...

#include <iostream>
#include <iomanip>
#include <limits>
#include <map>
#include <algorithm>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JTools/JHashCollection.hh"
#include "JTools/JHashMap.hh"
#include "JTools/JRange.hh"
#include "JTools/JQuantile.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JManager.hh"
#include "JLang/JComparator.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JEvent.hh"
#include "JAcoustics/JSupport.hh"
#include "Jeep/JContainer.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 monitor acoustic events.

Author
mdejong

Definition in file JAcousticsMonitor.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 78 of file JAcousticsMonitor.cc.

79 {
80  using namespace std;
81  using namespace JPP;
82 
83  typedef JContainer< vector<JTripod> > tripods_container;
84  typedef JContainer< vector<JTransmitter> > transmitters_container;
85 
86  JMultipleFileScanner<> inputFile;
87  JLimit_t& numberOfEvents = inputFile.getLimit();
88  string outputFile;
89  string detectorFile;
90  tripods_container tripods; // tripods
91  transmitters_container transmitters; // transmitters
92  double Q;
93  int debug;
94 
95  try {
96 
97  JParser<> zap("Example program to monitor acoustic events.");
98 
99  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
100  zap['n'] = make_field(numberOfEvents) = JLimit::max();
101  zap['o'] = make_field(outputFile) = "monitor.root";
102  zap['a'] = make_field(detectorFile);
103  zap['T'] = make_field(tripods, "tripod data") = JPARSER::initialised();
104  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
105  zap['Q'] = make_field(Q) = 0.9;
106  zap['d'] = make_field(debug) = 2;
107 
108  zap(argc, argv);
109  }
110  catch(const exception &error) {
111  FATAL(error.what() << endl);
112  }
113 
115 
116  try {
117  load(detectorFile, detector);
118  }
119  catch(const JException& error) {
120  FATAL(error);
121  }
122 
123  sort(detector.begin(), detector.end(), make_comparator(&JModule::getLocation));
124 
125  JHashMap<int, JModule_t> receivers;
127 
128  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
129  receivers[i->getID()] = JModule_t(i->getLocation(), i->getPosition());
130  }
131 
132  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
133  emitters[i->getID()] = (i->getUTMPosition() - detector.getUTMPosition());
134  }
135 
136  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
137  try {
138  emitters[i->getID()] = (i->getPosition() + detector.getModule(i->getLocation()).getPosition());
139  }
140  catch(const exception&) {
141  continue; // if no module available, discard transmitter
142  }
143  }
144 
145  const JHashCollection<int> string(make_array(detector.begin(), detector.end(), &JModule::getString));
146  const JRange <int> floor (make_array(detector.begin(), detector.end(), &JModule::getFloor));
147 
148  JManager<int, TH1D> HA(new TH1D("H[%].size", NULL, detector.size() + 1, -0.5, detector.size() + 0.5));
149  JManager<int, TH1D> HB(new TH1D("H[%].overlays", NULL, 201, -0.5, 200.5));
150  JManager<int, TH1D> HC(new TH1D("H[%].time", NULL, 800, 0.0, 4.0));
151  JManager<int, TH1D> HD(new TH1D("H[%].rms", NULL, 500, 0.0, 1.0e-1));
152  JManager<int, TH1D> HE(new TH1D("H[%].quantile", NULL, 500, 0.0, 1.0e-1));
153  JManager<int, TH1D> HQ(new TH1D("H[%].quality", NULL, 200, 0.0, 8.0));
154  JManager<int, TH2D> HR(new TH2D("H[%].QR", NULL, 40, 1.5, 3.5, 40, 3.0, 7.0));
155  JManager<int, TH1D> H1(new TH1D("H[%].multiplicity", NULL, 101, -0.5, 100.5));
156  JManager<int, TH2D> H2(new TH2D("H[%].event", NULL,
157  string.size(), -0.5, string.size() - 0.5,
158  floor.getUpperLimit() + 1, - 0.5, floor.getUpperLimit() + 0.5));
159 
160  for (Int_t i = 1; i <= H2->GetXaxis()->GetNbins(); ++i) {
161  H2->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(string.at(i-1)));
162  }
163 
164  for (Int_t i = 1; i <= H2->GetYaxis()->GetNbins(); ++i) {
165  H2->GetYaxis()->SetBinLabel(i, MAKE_CSTRING(i-1));
166  }
167 
168  JManager<int, TH2D> H3((TH2D*) H2->Clone("H[%].doubles"));
169 
170 
172 
173  counter_type counter = 0;
174 
175  for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
176 
177  JTreeScanner_t in(*i);
178 
180 
181  for (JTreeScanner_t::iterator event = in.begin(); event != in.end() && counter != inputFile.getLimit(); ++event, ++counter) {
182 
183  if (counter%1000 == 0) {
184  STATUS("event " << setw(8) << counter << '\r'); DEBUG(endl);
185  }
186 
187  const JEvent evt = *event;
188 
189  JQuantile Q1("", true);
190 
191  for (JEvent::const_iterator i = evt.begin(); i != evt.end(); ++i) {
192  Q1.put(i->getToE());
193  }
194 
195  HA[evt.getID()]->Fill((double) evt.size());
196  HB[evt.getID()]->Fill((double) evt.getOverlays());
197 
198  if (toe.has(evt.getID())) {
199  HC[evt.getID()]->Fill(log10(evt.begin()->getToE() - toe[evt.getID()]));
200  }
201 
202  HD[evt.getID()]->Fill(Q1.getSTDev());
203  HE[evt.getID()]->Fill(Q1.getQuantile(Q, JQuantile::symmetric_t));
204 
205  toe[evt.getID()] = evt.begin()->getToE();
206 
207  TH1D* hq = HQ[evt.getID()];
208  TH2D* hr = HR[evt.getID()];
209  TH1D* h1 = H1[evt.getID()];
210  TH2D* h2 = H2[evt.getID()];
211  TH2D* h3 = H3[evt.getID()];
212 
213  for (JEvent::const_iterator i = evt.begin(); i != evt.end(); ++i) {
214  hq->Fill(log10(i->getQ()));
215  }
216 
217  if (emitters.has(evt.getID())) {
218 
219  const JPosition3D& p0 = emitters[evt.getID()];
220 
221  for (JEvent::const_iterator i = evt.begin(); i != evt.end(); ++i) {
222  if (receivers.has(i->getID())) {
223  HR->Fill(log10(p0.getDistance(receivers[i->getID()])), log10(i->getQ()));
224  hr->Fill(log10(p0.getDistance(receivers[i->getID()])), log10(i->getQ()));
225  }
226  }
227  }
228 
229  map<int, set<double> > buffer;
230 
231  for (JEvent::const_iterator i = evt.begin(); i != evt.end(); ++i) {
232  buffer[i->getID()].insert(i->getQ());
233  }
234 
235  for (map<int, set<double> >::const_iterator i = buffer.begin(); i != buffer.end(); ++i) {
236 
237  h1->Fill((double) i->second.size());
238 
239  if (receivers.has(i->first)) {
240 
241  const JLocation& location = receivers[i->first];
242 
243  h2->Fill((double) string.getIndex(location.getString()), (double) location.getFloor(), 1.0);
244 
245  if (i->second.size() >= 2u) {
246  h3->Fill((double) string.getIndex(location.getString()), (double) location.getFloor(), 1.0);
247  }
248  }
249  }
250  }
251  }
252  STATUS(endl);
253 
254 
255  TFile out(outputFile.c_str(), "recreate");
256 
257  out << HA << HB << HC << HD << HE << HQ << HR << *HR << H1 << H2 << H3;
258 
259  out.Write();
260  out.Close();
261 }
Utility class to parse command line options.
Definition: JParser.hh:1517
General exception.
Definition: JException.hh:23
int getOverlays() const
Get overlays.
Q(UTCMax_s-UTCMin_s)-livetime_s
JComparator< JResult_t T::*, JComparison::lt > make_comparator(JResult_t T::*member)
Helper method to create comparator between values of data member.
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
General purpose class for hash map of unique keys.
Definition: JHashMap.hh:72
Auxiliary data structure for running average, standard deviation and quantiles.
Definition: JQuantile.hh:43
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
then fatal Number of tripods
Definition: JFootprint.sh:45
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:83
Long64_t counter_type
Type definition for counter.
string outputFile
double getDistance(const JVector3D &pos) const
Get distance to point.
Definition: JVector3D.hh:270
Template definition for direct access of elements in ROOT TChain.
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:226
Logical location of module.
Definition: JLocation.hh:37
Auxiliary wrapper for I/O of container with optional comment (see JComment).
Definition: JContainer.hh:39
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
set_variable E_E log10(E_{fit}/E_{#mu})"
const array_type< JValue_t > & make_array(const JValue_t(&array)[N])
Method to create array of values.
Definition: JVectorize.hh:54
JPosition3D getPosition(const Vec &pos)
Get position.
then awk string
int getIndex()
Get index for user I/O manipulation.
Definition: JManip.hh:26
#define FATAL(A)
Definition: JMessage.hh:67
Base class for JTreeScanner.
Definition: JTreeScanner.hh:54
int getString() const
Get string number.
Definition: JLocation.hh:134
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Acoustic event.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do set_variable DETECTOR_TXT $WORKDIR detector
int getID() const
Get identifier.
double u[N+1]
Definition: JPolint.hh:776
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 JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
bool has(const T &value) const
Test whether given value is present.
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62