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

Example program to plot acoustic fit results. More...

#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JHydrophone.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JROOT/JRootToolkit.hh"
#include "JROOT/JManager.hh"
#include "JAcoustics/JSoundVelocity.hh"
#include "JAcoustics/JGeometry.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JAcousticsToolkit.hh"
#include "JAcoustics/JHit.hh"
#include "JAcoustics/JEvent.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JEvtToolkit.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 plot acoustic fit results.

Author
mdejong

Definition in file JCanberra.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 43 of file JCanberra.cc.

44 {
45  using namespace std;
46  using namespace JPP;
47 
48  typedef JContainer< vector<JTripod> > tripods_container;
49  typedef JContainer< vector<JHydrophone> > hydrophones_container;
50 
52  JLimit_t& numberOfEvents = inputFile.getLimit();
53  string detectorFile;
54  string outputFile;
55  JSoundVelocity V = getSoundVelocity; // default sound velocity
56  tripods_container tripods; // tripods
57  hydrophones_container hydrophones; // hydrophones
58  int id; // emitter identifier
59  int debug;
60 
61  try {
62 
63  JParser<> zap("Example program to plot acoustic fit results.");
64 
65  zap['f'] = make_field(inputFile, "input file (output of JKatoomba)");
66  zap['n'] = make_field(numberOfEvents) = JLimit::max();
67  zap['a'] = make_field(detectorFile);
68  zap['o'] = make_field(outputFile) = "canberra.root";
69  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
70  zap['T'] = make_field(tripods, "tripod data");
71  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
72  zap['E'] = make_field(id, "emitter identifier (-1 = all)") = -1;
73  zap['d'] = make_field(debug) = 2;
74 
75  zap(argc, argv);
76  }
77  catch(const exception &error) {
78  FATAL(error.what() << endl);
79  }
80 
81 
83 
84  try {
85  load(detectorFile, detector);
86  }
87  catch(const JException& error) {
88  FATAL(error);
89  }
90 
91  JHashMap<int, JLocation> receivers;
92  JHashMap<int, JEmitter> emitters;
93 
94  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
95  receivers[i->getID()] = i->getLocation();
96  }
97 
98  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
99  emitters[i->getID()] = JEmitter(i->getID(),
100  i->getUTMPosition() - detector.getUTMPosition());
101  }
102 
103  const JGeometry geometry(detector, hydrophones);
104 
105  V.set(detector.getUTMZ());
106 
107 
108  JManager<int, TH1D> H2(new TH1D("[%].toa", NULL, 1000, -1.0e-2, +1.0e-2));
109 
110 
112 
113  JTreeScanner_t in(inputFile);
114 
115  JTreeScanner_t::iterator p = in.begin();
116 
117  while (inputFile.hasNext()) {
118 
119  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
120 
121  const JEvt* evt = inputFile.next();
122  const JModel model = getModel(*evt);
123 
124  DEBUG("model" << endl << model << endl);
125 
126  for ( ; p != in.end() && p-> begin()->getToA() < evt->UNIXTimeStart - 0.5; ++p) {}
127 
128  JTreeScanner_t::iterator q = p;
129 
130  for ( ; q != in.end() && q->rbegin()->getToA() <= evt->UNIXTimeStop + 0.5; ++q) {}
131 
132  DEBUG("events " << distance(p, q) << endl);
133 
134  for (JTreeScanner_t::iterator i = p; i != q; ++i) {
135 
136  if (id == i->getID() || id == -1) {
137 
138  if (emitters.has(i->getID())) {
139 
140  const JEmitter& emitter = emitters[i->getID()];
141 
142  for (JEvent::const_iterator hit = i->begin(); hit != i->end(); ++hit) {
143 
144  if (receivers.has(hit->getID())) {
145 
146  const JLocation& location = receivers[hit->getID()];
147  const JGEOMETRY::JString& string = geometry[location.getString()];
148  const JMODEL ::JString& parameters = model.string[location.getString()];
149  const JPosition3D position = string.getPosition(parameters, location.getFloor());
150 
151  const double D = emitter.getDistance(position);
152  const double Vi = V.getInverseVelocity(D, emitter.getZ(), position.getZ());
153  const double toa = hit->getToE() + D * Vi;
154 
155  H2[hit->getID()]->Fill(hit->getToA() - toa);
156  }
157  }
158  }
159  }
160  }
161 
162  p = q;
163  }
164  STATUS(endl);
165 
166  TFile out(outputFile.c_str(), "recreate");
167 
168  out << H2;
169 
170  out.Write();
171  out.Close();
172 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
do echo Generating $dir eval D
Definition: JDrawLED.sh:50
JModel getModel(const JEvt &evt)
Get model.
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
General purpose class for hash map of unique keys.
Definition: JHashMap.hh:71
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
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.
Definition: JTreeScanner.hh:91
Model for fit to acoustics data.
double UNIXTimeStop
stop time
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
static const JSoundVelocity getSoundVelocity(1541.0,-17.0e-3,-2000.0)
Function object for velocity of sound.
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:196
Acoustic event fit.
Acoustic emitter.
Definition: JEmitter.hh:27
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:1961
int debug
debug level
Definition: JSirene.cc:63
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
Implementation for velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
Base class for JTreeScanner.
Definition: JTreeScanner.hh:52
JACOUSTICS::JModel::string_type string
int getString() const
Get string number.
Definition: JLocation.hh:134
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
double UNIXTimeStart
start time
General purpose class for object reading from a list of file names.
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
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 source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:38
double getZ() const
Get z position.
Definition: JVector3D.hh:115
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62