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

Auxiliary program to draw the footprint of detector(s). More...

#include <string>
#include <iostream>
#include <limits>
#include <vector>
#include <map>
#include <set>
#include "TROOT.h"
#include "TH2D.h"
#include "TApplication.h"
#include "TGraph.h"
#include "TEllipse.h"
#include "TCanvas.h"
#include "TMarker.h"
#include "TText.h"
#include "TStyle.h"
#include "TLegend.h"
#include "TError.h"
#include "JLang/JSinglePointer.hh"
#include "JROOT/JCanvas.hh"
#include "JROOT/JStyle.hh"
#include "JROOT/JMarkerAttributes.hh"
#include "JROOT/JLegend.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JGeometry2D/JPosition2D.hh"
#include "JGeometry2D/JCircle2D.hh"
#include "Jeep/JeepToolkit.hh"
#include "Jeep/JContainer.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

Auxiliary program to draw the footprint of detector(s).

Author
mdejong

Definition in file JDrawDetector2D.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 177 of file JDrawDetector2D.cc.

178 {
179  using namespace std;
180  using namespace JPP;
181 
182  typedef JContainer< vector<JTripod> > tripods_container;
183  typedef JContainer< vector<JTransmitter> > transmitters_container;
184 
185  vector<string> detectorFile;
186  vector<string> tripodsFile;
187  vector<string> transmittersFile;
188  string outputFile;
189  JCanvas canvas;
190  string legend;
191  double markerSize;
192  double textSize;
193  bool drawCircle;
194  bool batch;
195  int debug;
196 
197  try {
198 
199  JParser<> zap("Auxiliary program to draw the footprint of detector(s).");
200 
201  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(500, 500);
202  zap['a'] = make_field(detectorFile, "detector file") = JPARSER::initialised();
203  zap['T'] = make_field(tripodsFile, "tripod file") = JPARSER::initialised();
204  zap['Y'] = make_field(transmittersFile, "transmitter file") = JPARSER::initialised();
205  zap['o'] = make_field(outputFile, "graphics output") = "";
206  zap['L'] = make_field(legend, "position legend e.g. TR") = "", "TL", "TR", "BR", "BL";
207  zap['S'] = make_field(markerSize, "marker size") = 1.0;
208  zap['s'] = make_field(textSize, "text size") = 0.02;
209  zap['C'] = make_field(drawCircle, "draw smallest enclosing cicrle");
210  zap['B'] = make_field(batch, "batch processing");
211  zap['d'] = make_field(debug) = 1;
212 
213  zap(argc, argv);
214  }
215  catch(const exception &error) {
216  FATAL(error.what() << endl);
217  }
218 
219 
220  if (detectorFile.empty() && tripodsFile.empty()) {
221  FATAL("No detector elements." << endl);
222  }
223 
224 
225  gROOT->SetBatch(batch);
226 
227  gErrorIgnoreLevel = kWarning;
228 
229  TApplication* tp = new TApplication("user", NULL, NULL);
230  TCanvas* cv = new TCanvas("detector", "", canvas.x, canvas.y);
231 
232  JSinglePointer<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
233 
234  gROOT->SetStyle("gplot");
235  gROOT->ForceStyle();
236 
237  cv->SetFillStyle(4000);
238  cv->SetFillColor(kWhite);
239  cv->Divide(1, 1);
240  cv->cd(1);
241 
242  JMarkerAttributes::getInstance().setMarkerSize(markerSize);
243 
244  vector<TAttText> text_attributes = {
245  TAttText(kHAlignLeft + kVAlignBottom, 0.25*PI, kBlack, 62, textSize),
246  TAttText(kHAlignRight + kVAlignBottom, 0.75*PI, kBlack, 62, textSize),
247  TAttText(kHAlignRight + kVAlignTop, 1.25*PI, kBlack, 62, textSize),
248  TAttText(kHAlignLeft + kVAlignTop, 1.75*PI, kBlack, 62, textSize)
249  };
250 
251  map<string, JGraph_t> data; // graphics data
252  JUTMPosition position; // UTM position
253  double offset = 0.0; // text offset
254 
255  for (size_t i = 0; i != detectorFile.size(); ++i) {
256 
258 
259  try {
260  load(detectorFile[i], detector);
261  }
262  catch(const JException& error) {
263  FATAL(error);
264  }
265 
266  position = detector.getUTMPosition();
267 
268  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
269  const TAttText& text = text_attributes[(0+i)%text_attributes.size()];
270  vector<JPoint_t>& buffer = data[getFilename(detectorFile[i])];
271 
272  if (offset == 0.0) {
273  offset = 3.0 * textSize * JCircle2D(detector.begin(), detector.end()).getRadius();
274  }
275 
276  set<int> counter;
277 
278  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
279 
280  if (counter.count(module->getString()) == 0) {
281 
282  buffer.push_back(JPoint_t(MAKE_STRING(module->getString()),
283  JPosition2D(module->getX(), module->getY()),
284  marker,
285  text,
286  offset));
287 
288  counter.insert(module->getString());
289  }
290  }
291  }
292 
293  for (size_t i = 0; i != tripodsFile.size(); ++i) {
294 
295  tripods_container tripods;
296 
297  tripods.load(tripodsFile[i].c_str());
298 
299  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
300  const TAttText& text = text_attributes[(0+i)%text_attributes.size()];
301  vector<JPoint_t>& buffer = data[getFilename(tripodsFile[i])];
302 
303  if (offset == 0.0) {
304  offset = 3.0 * textSize * JCircle2D(tripods.begin(), tripods.end()).getRadius();
305  }
306 
307  for (tripods_container::iterator i = tripods.begin(); i != tripods.end(); ++i) {
308 
309  buffer.push_back(JPoint_t(MAKE_STRING(i->getID()),
310  JPosition2D(i->getUTMEast() - position.getUTMEast(), i->getUTMNorth() - position.getUTMNorth()),
311  marker,
312  text,
313  offset));
314  }
315  }
316 
317  if (!transmittersFile.empty()) {
318 
319  if (transmittersFile.size() == detectorFile.size()) {
320 
321  for (size_t i = 0; i != transmittersFile.size(); ++i) {
322 
323  transmitters_container transmitters;
324 
325  transmitters.load(transmittersFile[i].c_str());
326 
328 
329  load(detectorFile[i], detector);
330 
331  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
332  const TAttText& text = text_attributes[(2+i)%text_attributes.size()];
333  vector<JPoint_t>& buffer = data[getFilename(transmittersFile[i])];
334 
335  if (offset == 0.0) {
336  offset = 3.0 * textSize * JCircle2D(transmitters.begin(), transmitters.end()).getRadius();
337  }
338 
339  for (transmitters_container::iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
340 
341  const JPosition3D pos = detector.getModule(i->getLocation()).getPosition();
342 
343  buffer.push_back(JPoint_t(MAKE_STRING(i->getID()),
344  JPosition2D(pos.getX() + i->getX(), pos.getY() + i->getY()),
345  marker,
346  text,
347  offset));
348  }
349  }
350 
351  } else {
352 
353  FATAL("Tranmitter files and detector files should match one-to-one." << endl);
354  }
355  }
356 
357  vector<JPosition2D> buffer;
358 
359  for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
360  copy(i->second.begin(), i->second.end(), back_inserter(buffer));
361  }
362 
363  JCircle2D circle(buffer.begin(), buffer.end()); // enclosing circle
364 
365  DEBUG("Detector (x,y,R): " << FIXED(12,3) << circle.getX() << ' ' << FIXED(12,3) << circle.getY() << ' ' << FIXED(9,3) << circle.getRadius() << endl);
366 
367  // center
368 
369  for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
370  i->second.sub(circle.getX(), circle.getY());
371  }
372 
373  circle.sub(circle.getPosition());
374 
375 
376  // draw
377 
378  cv->cd(1);
379 
380  const Double_t xmin = circle.getX() - 1.15 * circle.getRadius();
381  const Double_t xmax = circle.getX() + 1.15 * circle.getRadius();
382  const Double_t ymin = circle.getY() - 1.15 * circle.getRadius();
383  const Double_t ymax = circle.getY() + 1.15 * circle.getRadius();
384 
385  TH2D h2("h2", "", 1, xmin, xmax, 1, ymin, ymax);
386 
387  h2.GetXaxis()->SetTitle("x [m]");
388  h2.GetYaxis()->SetTitle("y [m]");
389 
390  h2.GetXaxis()->CenterTitle(true);
391  h2.GetYaxis()->CenterTitle(true);
392 
393  h2.SetStats(kFALSE);
394  h2.Draw("AXIS");
395 
396 
397  TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius());
398 
399  if (drawCircle) {
400  ellipse.Draw();
401  }
402 
403  for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
404  i->second.Draw();
405  }
406 
407  if (legend != "") {
408 
409  Ssiz_t height = data.size();
410  Ssiz_t width = 1;
411 
412  for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
413  width = max(width, (Ssiz_t) i->first.size());
414  }
415 
416  TLegend* lg = getLegend(width, height, legend);
417 
418  lg->SetTextSize(textSize);
419 
420  for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
421  if (!i->second.empty()) {
422  lg->AddEntry(&i->second[0].marker, i->first.c_str(), "P");
423  }
424  }
425 
426  lg->Draw();
427  }
428 
429  cv->Update();
430 
431  if (outputFile != "") {
432  cv->SaveAs(outputFile.c_str());
433  }
434 
435  if (!batch) {
436  tp->Run();
437  }
438 }
const double xmax
Definition: JQuadrature.cc:24
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
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
Detector data structure.
Definition: JDetector.hh:89
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:66
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
Data structure for UTM position.
Definition: JUTMPosition.hh:36
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
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
TLegend * getLegend(const Int_t width, const Int_t height, const std::string option, const Double_t factor=1.0)
Get legend.
Definition: JLegend.hh:29
JPosition3D getPosition(const Vec &pos)
Get position.
static const double PI
Mathematical constants.
double getY() const
Get y position.
Definition: JVector3D.hh:104
int debug
debug level
Definition: JSirene.cc:66
#define FATAL(A)
Definition: JMessage.hh:67
const JModule & getModule(const JModuleAddress &address) const
Get module parameters.
Definition: JDetector.hh:270
const double xmin
Definition: JQuadrature.cc:23
Data structure for position in two dimensions.
Definition: JPosition2D.hh:31
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
double getX() const
Get x position.
Definition: JVector3D.hh:94
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
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:128
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
do set_variable DETECTOR_TXT $WORKDIR detector
Wrapper class around ROOT TStyle.
Definition: JStyle.hh:20
Data structure for size of TCanvas.
Definition: JCanvas.hh:26
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62