Jpp  15.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JDrawDetector2D.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <limits>
4 #include <vector>
5 #include <map>
6 #include <set>
7 
8 #include "TROOT.h"
9 #include "TH2D.h"
10 #include "TApplication.h"
11 #include "TGraph.h"
12 #include "TEllipse.h"
13 #include "TCanvas.h"
14 #include "TMarker.h"
15 #include "TText.h"
16 #include "TStyle.h"
17 #include "TLegend.h"
18 #include "TError.h"
19 
20 #include "JLang/JSinglePointer.hh"
21 #include "JROOT/JCanvas.hh"
22 #include "JROOT/JStyle.hh"
24 #include "JROOT/JLegend.hh"
25 #include "JDetector/JDetector.hh"
27 #include "JDetector/JTripod.hh"
29 #include "JGeometry2D/JCircle2D.hh"
30 
31 #include "Jeep/JeepToolkit.hh"
32 #include "Jeep/JContainer.hh"
33 #include "Jeep/JParser.hh"
34 #include "Jeep/JMessage.hh"
35 
36 namespace {
37 
38  /**
39  * Marker with text.
40  */
41  struct JPoint_t {
42  /**
43  * Constructor.
44  *
45  * \param id identifier
46  * \param x x-position
47  * \param y y-position
48  * \param marker marker attributes
49  * \param text text attributes
50  * \param offset text ofset
51  */
52  JPoint_t(const int id,
53  const double x,
54  const double y,
55  const TAttMarker& marker,
56  const TAttText& text,
57  const double offset)
58  {
59  static_cast<TAttMarker&>(this->marker) = marker;
60  static_cast<TAttText&> (this->text) = text;
61 
62  this->marker.SetX(x);
63  this->marker.SetY(y);
64 
65  this->text.SetText(x + offset * cos(text.GetTextAngle()),
66  y + offset * sin(text.GetTextAngle()),
67  MAKE_CSTRING(id));
68 
69  // overwrite text attributes
70 
71  this->text.SetTextAngle(0.0);
72  this->text.SetTextColor(marker.GetMarkerColor());
73  }
74 
75  /**
76  * Add offset.
77  *
78  * \param x x-position
79  * \param y y-position
80  */
81  void add(const double x, const double y)
82  {
83  this->marker.SetX(this->marker.GetX() + x);
84  this->marker.SetY(this->marker.GetY() + y);
85  this->text .SetX(this->text .GetX() + x);
86  this->text .SetY(this->text .GetY() + y);
87  }
88 
89  /**
90  * Subtract offset.
91  *
92  * \param x x-position
93  * \param y y-position
94  */
95  void sub(const double x, const double y)
96  {
97  this->marker.SetX(this->marker.GetX() - x);
98  this->marker.SetY(this->marker.GetY() - y);
99  this->text .SetX(this->text .GetX() - x);
100  this->text .SetY(this->text .GetY() - y);
101  }
102 
103  /**
104  * Draw.
105  */
106  void Draw()
107  {
108  marker.Draw();
109  text .Draw();
110  }
111 
112  TMarker marker;
113  TText text;
114  };
115 
116 
117  /**
118  * Container of markers.
119  */
120  struct JGraph_t :
121  public std::vector<JPoint_t>
122  {
123  /**
124  * Add offset.
125  *
126  * \param x x-position
127  * \param y y-position
128  */
129  void add(const double x, const double y)
130  {
131  for (iterator i = begin(); i != end(); ++i) {
132  i->add(x, y);
133  }
134  }
135 
136  /**
137  * Subtract offset.
138  *
139  * \param x x-position
140  * \param y y-position
141  */
142  void sub(const double x, const double y)
143  {
144  for (iterator i = begin(); i != end(); ++i) {
145  i->sub(x, y);
146  }
147  }
148 
149  /**
150  * Draw.
151  */
152  void Draw()
153  {
154  for (iterator i = begin(); i != end(); ++i) {
155  i->Draw();
156  }
157  }
158  };
159 }
160 
161 
162 /**
163  * \file
164  * Auxiliary program to draw the footprint of detector(s).
165  * \author mdejong
166  */
167 int main(int argc, char**argv)
168 {
169  using namespace std;
170  using namespace JPP;
171 
172  typedef JContainer< vector<JTripod> > tripods_container;
173 
174  vector<string> detectorFile;
175  vector<string> tripodsFile;
176  string outputFile;
177  JCanvas canvas;
178  string legend;
179  double markerSize;
180  double textSize;
181  bool drawCircle;
182  bool batch;
183  int debug;
184 
185  try {
186 
187  JParser<> zap("Auxiliary program to draw the footprint of detector(s).");
188 
189  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(500, 500);
190  zap['a'] = make_field(detectorFile, "detector file") = JPARSER::initialised();
191  zap['T'] = make_field(tripodsFile, "tripod file") = JPARSER::initialised();
192  zap['o'] = make_field(outputFile, "graphics output") = "";
193  zap['L'] = make_field(legend, "position legend e.g. TR") = "", "TL", "TR", "BR", "BL";
194  zap['S'] = make_field(markerSize, "marker size") = 1.0;
195  zap['s'] = make_field(textSize, "text size") = 0.02;
196  zap['C'] = make_field(drawCircle, "draw smallest enclosing cicrle");
197  zap['B'] = make_field(batch, "batch processing");
198  zap['d'] = make_field(debug) = 1;
199 
200  zap(argc, argv);
201  }
202  catch(const exception &error) {
203  FATAL(error.what() << endl);
204  }
205 
206 
207  if (detectorFile.empty() && tripodsFile.empty()) {
208  FATAL("No detector elements." << endl);
209  }
210 
211 
212  gROOT->SetBatch(batch);
213 
214  gErrorIgnoreLevel = kWarning;
215 
216  TApplication* tp = new TApplication("user", NULL, NULL);
217  TCanvas* cv = new TCanvas("detector", "", canvas.x, canvas.y);
218 
219  JSinglePointer<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
220 
221  gROOT->SetStyle("gplot");
222  gROOT->ForceStyle();
223 
224  cv->SetFillStyle(4000);
225  cv->SetFillColor(kWhite);
226  cv->Divide(1, 1);
227  cv->cd(1);
228 
229  JMarkerAttributes::getInstance().setMarkerSize(markerSize);
230 
231  vector<TAttText> text_attributes = {
232  TAttText(kHAlignLeft + kVAlignBottom, 0.25*PI, kBlack, 62, textSize),
233  TAttText(kHAlignRight + kVAlignBottom, 0.75*PI, kBlack, 62, textSize),
234  TAttText(kHAlignRight + kVAlignTop, 1.25*PI, kBlack, 62, textSize),
235  TAttText(kHAlignLeft + kVAlignTop, 1.75*PI, kBlack, 62, textSize)
236  };
237 
238  map<string, JGraph_t> data; // graphics data
239  JCircle2D circle; // enclosing circle
240  JUTMPosition position; // UTM position
241 
242  for (size_t i = 0; i != detectorFile.size(); ++i) {
243 
245 
246  try {
247  load(detectorFile[i], detector);
248  }
249  catch(const JException& error) {
250  FATAL(error);
251  }
252 
253  JCircle2D c1(detector.begin(), detector.end());
254 
255  if (c1.getRadius() > circle.getRadius()) {
256  circle = c1;
257  }
258 
259  position = detector.getUTMPosition();
260 
261  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
262  const TAttText& text = text_attributes[i%text_attributes.size()];
263 
264  vector<JPoint_t>& buffer = data[getFilename(detectorFile[i])];
265 
266  set<int> counter;
267 
268  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
269 
270  if (counter.count(module->getString()) == 0) {
271 
272  buffer.push_back(JPoint_t(module->getString(),
273  module->getX(),
274  module->getY(),
275  marker,
276  text,
277  3.0 * textSize * c1.getRadius()));
278 
279  counter.insert(module->getString());
280  }
281  }
282  }
283 
284  for (size_t i = 0; i != tripodsFile.size(); ++i) {
285 
286  tripods_container tripods;
287 
288  ifstream in(tripodsFile[i].c_str());
289 
290  in >> tripods;
291 
292  in.close();
293 
294  if (detectorFile.empty()) {
295  position = getAverage(make_array(tripods.begin(), tripods.end(), &JTripod::getUTMPosition));
296  }
297 
298  JCircle2D c1(make_array(tripods.begin(), tripods.end(), &JTripod::getPosition));
299 
300  c1.sub(position.getPosition());
301 
302  if (c1.getRadius() > circle.getRadius()) {
303  circle = c1;
304  }
305 
306  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
307  const TAttText& text = text_attributes[i%text_attributes.size()];
308 
309  vector<JPoint_t>& buffer = data[getFilename(tripodsFile[i])];
310 
311  for (tripods_container::iterator i = tripods.begin(); i != tripods.end(); ++i) {
312 
313  buffer.push_back(JPoint_t(i->getID(),
314  i->getUTMEast() - position.getUTMEast(),
315  i->getUTMNorth() - position.getUTMNorth(),
316  marker,
317  text,
318  3.0 * textSize * c1.getRadius()));
319  }
320  }
321 
322 
323  DEBUG("Detector (x,y,R): " << FIXED(12,3) << circle.getX() << ' ' << FIXED(12,3) << circle.getY() << ' ' << FIXED(9,3) << circle.getRadius() << endl);
324 
325  // center
326 
327  for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
328  i->second.sub(circle.getX(), circle.getY());
329  }
330 
331  circle.sub(circle.getPosition());
332 
333 
334  // draw
335 
336  cv->cd(1);
337 
338  const Double_t xmin = circle.getX() - 1.15 * circle.getRadius();
339  const Double_t xmax = circle.getX() + 1.15 * circle.getRadius();
340  const Double_t ymin = circle.getY() - 1.15 * circle.getRadius();
341  const Double_t ymax = circle.getY() + 1.15 * circle.getRadius();
342 
343  TH2D h2("h2", "", 1, xmin, xmax, 1, ymin, ymax);
344 
345  h2.GetXaxis()->SetTitle("x [m]");
346  h2.GetYaxis()->SetTitle("y [m]");
347 
348  h2.GetXaxis()->CenterTitle(true);
349  h2.GetYaxis()->CenterTitle(true);
350 
351  h2.SetStats(kFALSE);
352  h2.Draw("AXIS");
353 
354 
355  TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius());
356 
357  if (drawCircle) {
358  ellipse.Draw();
359  }
360 
361  for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
362  i->second.Draw();
363  }
364 
365  if (legend != "") {
366 
367  Ssiz_t height = data.size();
368  Ssiz_t width = 1;
369 
370  for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
371  width = max(width, (Ssiz_t) i->first.size());
372  }
373 
374  TLegend* lg = getLegend(width, height, legend);
375 
376  lg->SetTextSize(textSize);
377 
378  for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
379  if (!i->second.empty()) {
380  lg->AddEntry(&i->second[0].marker, i->first.c_str(), "P");
381  }
382  }
383 
384  lg->Draw();
385  }
386 
387  cv->Update();
388 
389  if (outputFile != "") {
390  cv->SaveAs(outputFile.c_str());
391  }
392 
393  if (!batch) {
394  tp->Run();
395  }
396 }
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
double getRadius() const
Get radius.
Definition: JCircle2D.hh:144
int main(int argc, char *argv[])
Definition: Main.cc:15
char text[TEXT_SIZE]
Definition: elog.cc:72
const JUTMPosition & getUTMPosition() const
Get UTM position.
Definition: JUTMPosition.hh:84
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:33
Data structure for graph data.
Definition: JGraph.hh:21
Detector data structure.
Definition: JDetector.hh:89
std::iterator_traits< T >::value_type getAverage(T __begin, T __end)
Get average.
Definition: JMath.hh:497
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
Data structure for detector geometry and calibration.
Data structure for UTM position.
Definition: JUTMPosition.hh:36
JVector2D & sub(const JVector2D &vector)
Subtract vector.
Definition: JVector2D.hh:115
The template JSinglePointer class can be used to hold a pointer to an object.
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
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
const array_type< JValue_t > & make_array(const JValue_t(&array)[N])
Method to create array of values.
Definition: JVectorize.hh:54
Auxiliary methods for handling file names, type names and environment.
JPosition3D getPosition(const Vec &pos)
Get position.
static const double PI
Mathematical constants.
int debug
debug level
Definition: JSirene.cc:63
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
TCanvas * c1
Global variables to handle mouse events.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Utility class to parse command line options.
TLegend * getLegend(const Int_t width, const Int_t height, const std::string option="TR")
Get legend.
Definition: JLegend.hh:28
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
std::string getFilename(const std::string &file_name)
Get file name part, i.e. part after last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:88
do set_variable DETECTOR_TXT $WORKDIR detector
Data structure for tripod.
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:41
Wrapper class around ROOT TStyle.
Definition: JStyle.hh:20
Container I/O.
Data structure for size of TCanvas.
Definition: JCanvas.hh:26