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