Jpp - 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 "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 "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.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 39 of file JDrawDetector2D.cc.

40 {
41  using namespace std;
42  using namespace JPP;
43 
44  typedef JContainer< vector<JTripod> > tripods_container;
45 
46  vector<string> detectorFile;
47  tripods_container tripods;
48  string outputFile;
49  JCanvas canvas;
50  bool legend;
51  double markerSize;
52  double textSize;
53  bool batch;
54  int debug;
55 
56  try {
57 
58  JParser<> zap("Auxiliary program to draw the footprint of detector(s).");
59 
60  zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(500, 500);
61  zap['a'] = make_field(detectorFile, "detector file") = JPARSER::initialised();
62  zap['T'] = make_field(tripods, "tripod data") = JPARSER::initialised();
63  zap['o'] = make_field(outputFile, "graphics output") = "";
64  zap['L'] = make_field(legend, "optional legend");
65  zap['S'] = make_field(markerSize, "marker size") = 1.0;
66  zap['s'] = make_field(textSize, "text size") = 0.02;
67  zap['B'] = make_field(batch, "batch processing");
68  zap['d'] = make_field(debug) = 1;
69 
70  zap(argc, argv);
71  }
72  catch(const exception &error) {
73  FATAL(error.what() << endl);
74  }
75 
76 
77  if (detectorFile.empty() && tripods.empty()) {
78  FATAL("No detector elements." << endl);
79  }
80 
81 
82  gROOT->SetBatch(batch);
83 
84  gErrorIgnoreLevel = kWarning;
85 
86  TApplication* tp = new TApplication("user", NULL, NULL);
87  TCanvas* cv = new TCanvas("detector", "", canvas.x, canvas.y);
88 
89  JSinglePointer<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
90 
91  gROOT->SetStyle("gplot");
92  gROOT->ForceStyle();
93 
94  cv->SetFillStyle(4000);
95  cv->SetFillColor(kWhite);
96  cv->Divide(1, 1);
97  cv->cd(1);
98 
99  JMarkerAttributes::getInstance().setMarkerSize(markerSize);
100 
101 
102  struct alignment {
103  double angle;
104  int text;
105  };
106 
107  vector<alignment> align = {
108  { 0.25*PI, kHAlignLeft + kVAlignBottom },
109  { 0.75*PI, kHAlignRight + kVAlignBottom },
110  { 1.25*PI, kHAlignRight + kVAlignTop },
111  { 1.75*PI, kHAlignLeft + kVAlignTop }
112  };
113 
114 
115  int number_of_positions = 0;
116 
117  vector<TGraph*> data;
119 
120  JCircle2D circle;
121  JUTMPosition position;
122 
123  for (vector<string>::const_iterator file_name = detectorFile.begin(); file_name != detectorFile.end(); ++file_name) {
124 
126 
127  try {
128  load(*file_name, detector);
129  }
130  catch(const JException& error) {
131  FATAL(error);
132  }
133 
134  if (detector.empty()) {
135  ERROR("Empty detector." << endl);
136  continue;
137  }
138 
139  const JCircle2D c1(detector.begin(), detector.end());
140 
141  if (c1.getRadius() > circle.getRadius()) {
142  circle = c1;
143  }
144 
145  position = detector.getUTMPosition();
146 
147  const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
148 
149  TGraph* graph = new TGraph(detector.size());
150 
151  graph->SetTitle(getFilename(*file_name).c_str());
152 
153  dynamic_cast<TAttMarker&>(*graph) = marker;
154 
155  int n = 0;
156 
157  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module, ++n) {
158 
159  number_of_positions += 1;
160 
161  graph->GetX()[n] = module->getX();
162  graph->GetY()[n] = module->getY();
163 
164  if (module->getFloor() == 1) {
165 
166  const alignment& a1 = align[distance(detectorFile.cbegin(), file_name) % align.size()];
167 
168  text.push_back(new TText(module->getX() + 3.0 * textSize * c1.getRadius()*cos(a1.angle),
169  module->getY() + 3.0 * textSize * c1.getRadius()*sin(a1.angle),
170  MAKE_CSTRING(module->getString())));
171 
172  (*text.rbegin())->SetTextSize (textSize);
173  (*text.rbegin())->SetTextAlign(a1.text);
174  (*text.rbegin())->SetTextColor(marker.GetMarkerColor());
175  }
176  }
177 
178  data.push_back(graph);
179  }
180 
181  if (!tripods.empty()) {
182 
183  if (number_of_positions == 0) {
184 
185  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
186  position += *i;
187  }
188 
189  position /= tripods.size();
190  }
191 
192  vector<JPosition3D> zbuf;
193 
194  for (tripods_container::iterator i = tripods.begin(); i != tripods.end(); ++i) {
195  zbuf.push_back(*i - position);
196  }
197 
198  circle = JCircle2D(zbuf.begin(), zbuf.end());
199 
200  TGraph* graph = new TGraph(zbuf.size());
201 
202  graph->SetTitle("tripod");
203 
204  for (size_t i = 0; i != zbuf.size(); ++i) {
205  graph->GetX()[i] = zbuf[i].getPosition().getX();
206  graph->GetY()[i] = zbuf[i].getPosition().getY();
207  }
208 
209  dynamic_cast<TAttMarker&>(*graph) = *JMarkerAttributes::getInstance().next();
210 
211  data.push_back(graph);
212  }
213 
214  if (number_of_positions == 0 && tripods.empty()) {
215  FATAL("No detector elements." << endl);
216  }
217 
218 
219  DEBUG("Detector (x,y,R): " << FIXED(12,3) << circle.getX() << ' ' << FIXED(12,3) << circle.getY() << ' ' << FIXED(9,3) << circle.getRadius() << endl);
220 
221  for (vector<TGraph*>::iterator graph = data.begin(); graph != data.end(); ++graph) {
222  for (int i = 0; i != (*graph)->GetN(); ++i) {
223  (*graph)->GetX()[i] -= circle.getX();
224  (*graph)->GetY()[i] -= circle.getY();
225  }
226  }
227 
228  for (vector<TText*>::iterator i = text.begin(); i != text.end(); ++i) {
229  (*i)->SetX((*i)->GetX() - circle.getX());
230  (*i)->SetY((*i)->GetY() - circle.getY());
231  }
232 
233  circle.sub(circle.getPosition());
234 
235  Double_t R = circle.getRadius();
236 
237  if (R <= 1.0) {
238  R = 1.0;
239  }
240 
241  Double_t xmin = circle.getX() - R;
242  Double_t xmax = circle.getX() + R;
243  Double_t ymin = circle.getY() - R;
244  Double_t ymax = circle.getY() + R;
245 
246 
247  const Double_t dx = (xmax - xmin);
248  const Double_t dy = (ymax - ymin);
249 
250  xmin -= 0.1 * dx;
251  xmax += 0.1 * dx;
252  ymin -= 0.1 * dy;
253  ymax += 0.1 * dy;
254 
255 
256  cv->cd(1);
257 
258  TH2D h2("h2", "", 1, xmin, xmax, 1, ymin, ymax);
259 
260  h2.GetXaxis()->SetTitle("x [m]");
261  h2.GetYaxis()->SetTitle("y [m]");
262 
263  h2.GetXaxis()->CenterTitle(true);
264  h2.GetYaxis()->CenterTitle(true);
265 
266  h2.SetStats(kFALSE);
267  h2.Draw();
268 
269 
270  TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius());
271 
272  ellipse.Draw();
273 
274  for (vector<TGraph*>::iterator graph = data.begin(); graph != data.end(); ++graph) {
275  (*graph)->Draw("P");
276  }
277 
278  for (vector<TText*>::iterator i = text.begin(); i != text.end(); ++i) {
279  (*i)->Draw();
280  }
281 
282 
283  if (legend) {
284 
285  const Double_t x2 = gPad->GetX2() - gStyle->GetPadRightMargin();
286  const Double_t y2 = gPad->GetY2() - gStyle->GetPadTopMargin();
287 
288  size_t length = 10;
289  Double_t font_size = 0.03; // gStyle->GetStatFontSize();
290 
291  for (vector<TGraph*>::iterator graph = data.begin(); graph != data.end(); ++graph) {
292  length = max(length, strlen((*graph)->GetTitle()));
293  }
294 
295  TLegend* legend = new TLegend(x2 - length * font_size * 0.65 - 0.02,
296  y2 - data.size() * font_size,
297  x2 - 0.01,
298  y2 - 0.01);
299 
300  //legend->SetFillStyle(4000);
301  legend->SetFillColor(0);
302  legend->SetBorderSize(0);
303  legend->SetTextSize(font_size);
304 
305  for (vector<TGraph*>::iterator graph = data.begin(); graph != data.end(); ++graph) {
306  legend->AddEntry(*graph, (*graph)->GetTitle(), "P");
307  }
308 
309  legend->Draw();
310  }
311 
312  cv->Update();
313 
314  if (outputFile != "") {
315  cv->SaveAs(outputFile.c_str());
316  }
317 
318  if (!batch) {
319  tp->Run();
320  }
321 }
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
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:30
Detector data structure.
Definition: JDetector.hh:80
#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 UTM position.
Definition: JUTMPosition.hh:36
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
#define ERROR(A)
Definition: JMessage.hh:66
static const double PI
Mathematical constants.
int debug
debug level
Definition: JSirene.cc:63
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:40
#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.
alias put_queue eval echo n
Definition: qlib.csh:19
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
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