Jpp 19.3.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
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#include <memory>
8
9#include "TROOT.h"
10#include "TH2D.h"
11#include "TApplication.h"
12#include "TGraph.h"
13#include "TEllipse.h"
14#include "TCanvas.h"
15#include "TRootCanvas.h"
16#include "TMarker.h"
17#include "TText.h"
18#include "TStyle.h"
19#include "TLegend.h"
20#include "TError.h"
21
22#include "JROOT/JCanvas.hh"
23#include "JROOT/JStyle.hh"
25#include "JROOT/JLegend.hh"
28#include "JDetector/JTripod.hh"
32
33#include "Jeep/JeepToolkit.hh"
34#include "Jeep/JContainer.hh"
35#include "Jeep/JProperties.hh"
36#include "Jeep/JParser.hh"
37#include "Jeep/JMessage.hh"
38
39namespace {
40
42
43 /**
44 * Marker with text.
45 */
46 struct JPoint_t :
47 public JPosition2D
48 {
49 /**
50 * Constructor.
51 *
52 * \param title title
53 * \param pos position
54 * \param marker marker attributes
55 * \param text text attributes
56 * \param offset text ofset
57 */
58 JPoint_t(const std::string& title,
59 const JPosition2D& pos,
60 const TAttMarker& marker,
61 const TAttText& text,
62 const double offset) :
63 JPosition2D(pos),
64 offset(offset * cos(text.GetTextAngle()),
65 offset * sin(text.GetTextAngle()))
66 {
67 static_cast<TAttMarker&>(this->marker) = marker;
68 static_cast<TAttText&> (this->text) = text;
69
70 this->text.SetTitle(title.c_str());
71 this->text.SetTextAngle(0.0); // overwrite text attributes
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 static_cast<JPosition2D&>(*this).add(JPosition2D(x,y));
84
85 configure();
86 }
87
88 /**
89 * Subtract offset.
90 *
91 * \param x x-position
92 * \param y y-position
93 */
94 void sub(const double x, const double y)
95 {
96 static_cast<JPosition2D&>(*this).sub(JPosition2D(x,y));
97
98 configure();
99 }
100
101 /**
102 * Draw.
103 */
104 void Draw()
105 {
106 marker.Draw();
107 text .Draw();
108 }
109
110 /**
111 * Configure.
112 */
113 void configure()
114 {
115 this->marker.SetX(this->getX());
116 this->marker.SetY(this->getY());
117 this->text .SetX(this->getX() + offset.getX());
118 this->text .SetY(this->getY() + offset.getY());
119 }
120
121 TMarker marker;
122 TText text;
123
124 private:
125 JPosition2D offset;
126 };
127
128
129 /**
130 * Container of markers.
131 */
132 struct JGraph_t :
133 public std::vector<JPoint_t>
134 {
135 /**
136 * Add offset.
137 *
138 * \param x x-position
139 * \param y y-position
140 */
141 void add(const double x, const double y)
142 {
143 for (iterator i = begin(); i != end(); ++i) {
144 i->add(x, y);
145 }
146 }
147
148 /**
149 * Subtract offset.
150 *
151 * \param x x-position
152 * \param y y-position
153 */
154 void sub(const double x, const double y)
155 {
156 for (iterator i = begin(); i != end(); ++i) {
157 i->sub(x, y);
158 }
159 }
160
161 /**
162 * Draw.
163 */
164 void Draw()
165 {
166 for (iterator i = begin(); i != end(); ++i) {
167 i->Draw();
168 }
169 }
170 };
171}
172
173
174/**
175 * \file
176 * Auxiliary program to draw the footprint of detector(s).
177 * \author mdejong
178 */
179int main(int argc, char**argv)
180{
181 using namespace std;
182 using namespace JPP;
183
186
187 vector<string> detectorFile;
188 vector<string> tripodsFile;
189 vector<string> transmittersFile;
190 string outputFile;
191 JCanvas canvas;
192 JCircle2D focus;
193 string legend;
194 double markerSize;
195 double textSize;
196 bool drawCircle;
197 bool all;
198 bool batch;
199 int debug;
200
201 try {
202
203 JProperties properties;
204
205 properties["focus"] = focus;
206
207 JParser<> zap("Auxiliary program to draw the footprint of detector(s).");
208
209 zap['w'] = make_field(canvas, "size of canvas <nx>x<ny> [pixels]") = JCanvas(500, 500);
210 zap['@'] = make_field(properties, "properties") = JPARSER::initialised();
211 zap['a'] = make_field(detectorFile, "detector file") = JPARSER::initialised();
212 zap['T'] = make_field(tripodsFile, "tripod file") = JPARSER::initialised();
213 zap['Y'] = make_field(transmittersFile, "transmitter file") = JPARSER::initialised();
214 zap['o'] = make_field(outputFile, "graphics output") = "";
215 zap['L'] = make_field(legend, "position legend e.g. TR") = "", "TL", "TR", "BR", "BL";
216 zap['S'] = make_field(markerSize, "marker size") = 1.0;
217 zap['s'] = make_field(textSize, "text size") = 0.02;
218 zap['C'] = make_field(drawCircle, "draw smallest enclosing cicrle");
219 zap['A'] = make_field(all, "draw all modules");
220 zap['B'] = make_field(batch, "batch processing");
221 zap['d'] = make_field(debug) = 1;
222
223 zap(argc, argv);
224 }
225 catch(const exception &error) {
226 FATAL(error.what() << endl);
227 }
228
229
230 if (detectorFile.empty() && tripodsFile.empty()) {
231 FATAL("No detector elements." << endl);
232 }
233
234
235 gROOT->SetBatch(batch);
236
237 gErrorIgnoreLevel = kWarning;
238
239 TApplication* tp = new TApplication("user", NULL, NULL);
240 TCanvas* cv = new TCanvas("detector", "", canvas.x, canvas.y);
241
242 if (!batch) {
243 ((TRootCanvas *) cv->GetCanvasImp())->Connect("CloseWindow()", "TApplication", tp, "Terminate()");
244 }
245
246 unique_ptr<TStyle> gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh()));
247
248 gROOT->SetStyle("gplot");
249 gROOT->ForceStyle();
250
251 cv->SetFillStyle(4000);
252 cv->SetFillColor(kWhite);
253 cv->Divide(1, 1);
254 cv->cd(1);
255
256 JMarkerAttributes::getInstance().setMarkerSize(markerSize);
257
258 vector<TAttText> text_attributes = {
259 TAttText(kHAlignLeft + kVAlignBottom, 0.25*PI, kBlack, 62, textSize),
260 TAttText(kHAlignRight + kVAlignBottom, 0.75*PI, kBlack, 62, textSize),
261 TAttText(kHAlignRight + kVAlignTop, 1.25*PI, kBlack, 62, textSize),
262 TAttText(kHAlignLeft + kVAlignTop, 1.75*PI, kBlack, 62, textSize)
263 };
264
265 map<string, JGraph_t> data; // graphics data
266 JUTMPosition position; // UTM position
267 double offset = 0.0; // text offset
268
269 for (size_t i = 0; i != detectorFile.size(); ++i) {
270
272
273 try {
274 load(detectorFile[i], detector);
275 }
276 catch(const JException& error) {
277 FATAL(error);
278 }
279
280 position = detector.getUTMPosition();
281
282 const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
283 const TAttText& text = text_attributes[(0+i)%text_attributes.size()];
284 vector<JPoint_t>& buffer = data[getFilename(detectorFile[i])];
285
286 if (offset == 0.0) {
287 offset = 3.0 * textSize * JCircle2D(detector.begin(), detector.end()).getRadius();
288 }
289
290 set<int> counter;
291
292 for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
293
294 if (module->getFloor() == 0) {
295
296 buffer.push_back(JPoint_t(MAKE_STRING(module->getString()),
297 JPosition2D(module->getX(), module->getY()),
298 marker,
299 text,
300 offset));
301
302 counter.insert(module->getString());
303 }
304 }
305
306 for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
307
308 const bool status = (counter.count(module->getString()) == 0);
309
310 if (status || all) {
311
312 buffer.push_back(JPoint_t((status ? MAKE_STRING(module->getString()) : ""),
313 JPosition2D(module->getX(), module->getY()),
314 marker,
315 text,
316 offset));
317
318 counter.insert(module->getString());
319 }
320 }
321 }
322
323 for (size_t i = 0; i != tripodsFile.size(); ++i) {
324
325 tripods_container tripods;
326
327 tripods.load(tripodsFile[i].c_str());
328
329 const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
330 const TAttText& text = text_attributes[(0+i)%text_attributes.size()];
331 vector<JPoint_t>& buffer = data[getFilename(tripodsFile[i])];
332
333 if (offset == 0.0) {
334 offset = 3.0 * textSize * JCircle2D(tripods.begin(), tripods.end()).getRadius();
335 }
336
337 for (tripods_container::iterator i = tripods.begin(); i != tripods.end(); ++i) {
338
339 buffer.push_back(JPoint_t(MAKE_STRING(i->getID()),
340 JPosition2D(i->getUTMEast() - position.getUTMEast(), i->getUTMNorth() - position.getUTMNorth()),
341 marker,
342 text,
343 offset));
344 }
345 }
346
347 if (!transmittersFile.empty()) {
348
349 if (transmittersFile.size() == detectorFile.size()) {
350
351 for (size_t i = 0; i != transmittersFile.size(); ++i) {
352
353 transmitters_container transmitters;
354
355 transmitters.load(transmittersFile[i].c_str());
356
358
359 load(detectorFile[i], detector);
360
361 const TAttMarker& marker = *JMarkerAttributes::getInstance().next();
362 const TAttText& text = text_attributes[(2+i)%text_attributes.size()];
363 vector<JPoint_t>& buffer = data[getFilename(transmittersFile[i])];
364
365 if (offset == 0.0) {
366 offset = 3.0 * textSize * JCircle2D(transmitters.begin(), transmitters.end()).getRadius();
367 }
368
369 for (transmitters_container::iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
370
371 const JPosition3D pos = detector.getModule(i->getLocation()).getPosition();
372
373 buffer.push_back(JPoint_t(MAKE_STRING(i->getID()),
374 JPosition2D(pos.getX() + i->getX(), pos.getY() + i->getY()),
375 marker,
376 text,
377 offset));
378 }
379 }
380
381 } else {
382
383 FATAL("Tranmitter files and detector files should match one-to-one." << endl);
384 }
385 }
386
387 vector<JPosition2D> buffer;
388
389 for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
390 copy(i->second.begin(), i->second.end(), back_inserter(buffer));
391 }
392
393 JCircle2D circle(buffer.begin(), buffer.end()); // enclosing circle
394
395 if (focus.getRadius() > 0.0) {
396 circle = focus;
397 }
398
399 NOTICE("focus = " << FIXED(12,3) << circle.getX() << ' ' << FIXED(12,3) << circle.getY() << ' ' << FIXED(9,3) << circle.getRadius() << endl);
400
401 // center
402
403 for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
404 i->second.sub(circle.getX(), circle.getY());
405 }
406
407 circle.sub(circle.getPosition());
408
409
410 // draw
411
412 cv->cd(1);
413
414 const Double_t xmin = circle.getX() - 1.15 * circle.getRadius();
415 const Double_t xmax = circle.getX() + 1.15 * circle.getRadius();
416 const Double_t ymin = circle.getY() - 1.15 * circle.getRadius();
417 const Double_t ymax = circle.getY() + 1.15 * circle.getRadius();
418
419 TH2D h2("h2", "", 200, xmin, xmax, 200, ymin, ymax);
420
421 h2.GetXaxis()->SetTitle("x [m]");
422 h2.GetYaxis()->SetTitle("y [m]");
423
424 h2.GetXaxis()->CenterTitle(true);
425 h2.GetYaxis()->CenterTitle(true);
426
427 h2.SetStats(kFALSE);
428 h2.Draw("AXIS");
429
430
431 TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius());
432
433 if (drawCircle) {
434 ellipse.Draw();
435 }
436
437 for (map<string, JGraph_t>::iterator i = data.begin(); i != data.end(); ++i) {
438 i->second.Draw();
439 }
440
441 if (legend != "") {
442
443 Ssiz_t height = data.size();
444 Ssiz_t width = 1;
445
446 for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
447 width = max(width, (Ssiz_t) i->first.size());
448 }
449
450 TLegend* lg = getLegend(width, height, legend);
451
452 lg->SetTextSize(textSize);
453
454 for (map<string, JGraph_t>::const_iterator i = data.begin(); i != data.end(); ++i) {
455 if (!i->second.empty()) {
456 lg->AddEntry(&i->second[0].marker, i->first.c_str(), "P");
457 }
458 }
459
460 lg->Draw();
461 }
462
463 cv->Update();
464
465 if (outputFile != "") {
466 cv->SaveAs(outputFile.c_str());
467 }
468
469 if (!batch) {
470 tp->Run();
471 }
472}
Container I/O.
string outputFile
Data structure for detector geometry and calibration.
int main(int argc, char **argv)
General purpose messaging.
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
Utility class to parse parameter values.
Data structure for transmitter.
Data structure for tripod.
Auxiliary methods for handling file names, type names and environment.
Detector data structure.
Definition JDetector.hh:96
Utility class to parse parameter values.
Data structure for circle in two dimensions.
Definition JCircle2D.hh:35
double getRadius() const
Get radius.
Definition JCircle2D.hh:144
Data structure for position in two dimensions.
const JPosition2D & getPosition() const
Get position.
double getY() const
Get y position.
Definition JVector2D.hh:74
double getX() const
Get x position.
Definition JVector2D.hh:63
JVector2D & sub(const JVector2D &vector)
Subtract vector.
Definition JVector2D.hh:115
Data structure for position in three dimensions.
double getY() const
Get y position.
Definition JVector3D.hh:104
double getX() const
Get x position.
Definition JVector3D.hh:94
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
Data structure for size of TCanvas.
Definition JCanvas.hh:26
int y
number of pixels in Y
Definition JCanvas.hh:99
int x
number of pixels in X
Definition JCanvas.hh:98
Wrapper class around ROOT TStyle.
Definition JStyle.hh:24
Data structure for UTM position.
double getUTMNorth() const
Get UTM north.
double getUTMEast() const
Get UTM east.
char text[TEXT_SIZE]
Definition elog.cc:72
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
void configure(const T &value, const JAbstractCollection< JAbscissa_t > &bounds, JBool< false > option)
Configuration of value.
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
Auxiliary wrapper for I/O of container with optional comment (see JComment).
Definition JContainer.hh:42
void load(const char *file_name)
Load from input file.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68