Jpp  18.2.0
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 <vector>
#include <set>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TF1.h"
#include "TFitResult.h"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JDetector/JHydrophone.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JModuleRouter.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 "JAcoustics/JTransmission_t.hh"
#include "Jeep/JContainer.hh"
#include "Jeep/JProperties.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 51 of file JCanberra.cc.

52 {
53  using namespace std;
54  using namespace JPP;
55 
59  typedef JContainer< set<JTransmission_t> > disable_container;
60 
62  JLimit_t& numberOfEvents = inputFile.getLimit();
63  string detectorFile;
64  string outputFile;
65  JSoundVelocity V = getSoundVelocity; // default sound velocity
66  tripods_container tripods; // tripods
67  transmitters_container transmitters; // transmitters
68  hydrophones_container hydrophones; // hydrophones
69  int id; // emitter identifier
70  disable_container disable; // disable tansmissions
71  string option;
72  int numberOfEntries = 21;
73  int debug;
74 
75  try {
76 
77  JParser<> zap("Example program to plot acoustic fit results.");
78 
79  JProperties properties;
80 
81  properties.insert(gmake_property(numberOfEntries));
82 
83  zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh])");
84  zap['n'] = make_field(numberOfEvents) = JLimit::max();
85  zap['a'] = make_field(detectorFile);
86  zap['o'] = make_field(outputFile) = "canberra.root";
87  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
88  zap['T'] = make_field(tripods, "tripod data");
89  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
90  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
91  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
92  zap['E'] = make_field(id, "emitter identifier (-1 = all)") = -1;
93  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
94  zap['O'] = make_field(option, "ROOT fit option, see TH1::Fit.") = "LS";
95  zap['@'] = make_field(properties) = JPARSER::initialised();
96  zap['d'] = make_field(debug) = 2;
97 
98  zap(argc, argv);
99  }
100  catch(const exception &error) {
101  FATAL(error.what() << endl);
102  }
103 
104 
105  if (option.find('Q') == string::npos && debug < JEEP::debug_t) { option += 'Q'; }
106 
108 
109  try {
110  load(detectorFile, detector);
111  }
112  catch(const JException& error) {
113  FATAL(error);
114  }
115 
116  JHashMap<int, JLocation> receivers;
117  JHashMap<int, JEmitter> emitters;
118 
119  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
120  receivers[i->getID()] = i->getLocation();
121  }
122 
123  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
124  emitters[i->getID()] = JEmitter(i->getID(),
125  i->getUTMPosition() - detector.getUTMPosition());
126  }
127 
128  const JGeometry geometry(detector, hydrophones);
129 
130  V.set(detector.getUTMZ());
131 
132 
133  JManager<JTransmission_t, TH1D> H2(new TH1D("[%].toa", NULL, 1000, -1.0e-2, +1.0e-2));
134 
135  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
136  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
137  if (id == i->getID() || id == -1) {
138  H2[JTransmission_t(i->getID(), module->getID())];
139  }
140  }
141  }
142 
144 
145  JTreeScanner_t in(inputFile);
146 
147  JTreeScanner_t::iterator p = in.begin();
148 
149  while (inputFile.hasNext()) {
150 
151  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
152 
153  const JEvt* evt = inputFile.next();
154  const JModel model = getModel(*evt);
155 
156  DEBUG("model" << endl << model << endl);
157 
158  for ( ; p != in.end() && p-> begin()->getToA() < evt->UNIXTimeStart - 0.5; ++p) {}
159 
160  JTreeScanner_t::iterator q = p;
161 
162  for ( ; q != in.end() && q->rbegin()->getToA() <= evt->UNIXTimeStop + 0.5; ++q) {}
163 
164  DEBUG("events " << distance(p, q) << endl);
165 
166  for (JTreeScanner_t::iterator i = p; i != q; ++i) {
167 
168  if (id == i->getID() || id == -1) {
169 
170  if (emitters.has(i->getID())) {
171 
172  const JEmitter& emitter = emitters[i->getID()];
173 
174  for (JEvent::const_iterator hit = i->begin(); hit != i->end(); ++hit) {
175 
176  if (receivers.has(hit->getID())) {
177 
178  if (disable.count(JTransmission_t(i->getID(), hit->getID())) == 0 &&
179  disable.count(JTransmission_t(-1, hit->getID())) == 0) {
180 
181  const JLocation& location = receivers[hit->getID()];
182 
183  if (geometry .has(location.getString()) &&
184  model.string.has(location.getString())) {
185 
186  const JGEOMETRY::JString& string = geometry [location.getString()];
187  const JMODEL ::JString& parameters = model.string[location.getString()];
188  const JPosition3D position = string.getPosition(parameters, location.getFloor());
189 
190  const double D = emitter.getDistance(position);
191  const double Vi = V.getInverseVelocity(D, emitter.getZ(), position.getZ());
192  const double toa = hit->getToE() + D * Vi;
193 
194  H2[JTransmission_t(i->getID(), hit->getID())]->Fill(hit->getToA() - toa);
195  }
196  }
197  }
198  }
199  }
200  }
201  }
202 
203  p = q;
204  }
205  STATUS(endl);
206 
207 
209 
210  JManager<int, TH2D> HA(new TH2D("[%].mean", NULL,
212  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
213 
214  JManager<int, TH2D> HB(new TH2D("[%].sigma", NULL,
216  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
217 
218  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
219  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
220  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
221  }
222 
223  const JModuleRouter router(detector);
224 
225  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))");
226 
227  const vector<Double_t> Q = { 0.33, 0.55, 0.66 };
228 
229  for (JManager<JTransmission_t, TH1D>::iterator i = H2.begin(); i != H2.end(); ++i) {
230 
231  TH1* h1 = i->second;
232 
233  if (h1->GetEntries() >= numberOfEntries) {
234 
235  vector<Double_t> X(Q.size(), 0.0);
236 
237  h1->GetQuantiles(Q.size(), X.data(), Q.data());
238 
239  f1.SetParameter(0, h1->GetMaximum());
240  f1.SetParameter(1, X[1]);
241  f1.SetParameter(2, 0.5*(X[2] - X[0]));
242 
243  TFitResultPtr result = h1->Fit(&f1, option.c_str(), "same", X[1] - 5.0 * (X[2] - X[0]), X[1] + 5.0 * (X[2] - X[0]));
244 
245  if (result.Get() == NULL) {
246 
247  ERROR("Invalid TFitResultPtr " << h1->GetName() << endl);
248 
249  } else {
250 
251  const JLocation location = router.getModule(i->first.rx).getLocation();
252 
253  HA[i->first.tx]->Fill(geometry.getIndex(location.getString()), location.getFloor(), f1.GetParameter(1) * 1.0e3);
254  HB[i->first.tx]->Fill(geometry.getIndex(location.getString()), location.getFloor(), f1.GetParameter(2) * 1.0e3);
255  }
256  }
257  }
258 
259 
260  TFile out(outputFile.c_str(), "recreate");
261 
262  out << H2 << HA << HB;
263 
264  out.Write();
265  out.Close();
266 }
JModel getModel(const JEvt &evt)
Get model.
Utility class to parse command line options.
Definition: JParser.hh:1514
General exception.
Definition: JException.hh:24
Q(UTCMax_s-UTCMin_s)-livetime_s
debug
Definition: JMessage.hh:29
int getFloor() const
Get floor number.
Definition: JLocation.hh:145
JContainer< std::vector< JTransmitter > > transmitters_container
Definition: JSydney.cc:79
static JDetectorMechanics getMechanics
Function object to get string mechanics.
Definition: JMechanics.hh:242
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:72
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
#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.
Utility class to parse parameter values.
Definition: JProperties.hh:497
*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
#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
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
string outputFile
double getDistance(const JVector3D &pos) const
Get distance to point.
Definition: JVector3D.hh:270
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
Template definition for direct access of elements in ROOT TChain.
const JPolynome f1(1.0, 2.0, 3.0)
Function.
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:226
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition: JSydney.cc:78
Acoustic event fit.
Acoustic emitter.
Definition: JEmitter.hh:27
Logical location of module.
Definition: JLocation.hh:37
const JLocation & getLocation() const
Get location.
Definition: JLocation.hh:69
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:1989
#define ERROR(A)
Definition: JMessage.hh:66
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
JContainer< std::vector< JTripod > > tripods_container
Definition: JSydney.cc:77
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
Base class for JTreeScanner.
Definition: JTreeScanner.hh:54
JACOUSTICS::JModel::string_type string
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
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.
no fit printf nominal n $STRING awk v X
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
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 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
Acoustic transmission identifier.
bool has(const T &value) const
Test whether given value is present.
do echo Generating $dir eval D
Definition: JDrawLED.sh:53
double getZ() const
Get z position.
Definition: JVector3D.hh:115
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62