#include #include #include #include #include #include #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/JProperties.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" namespace { using JGEOMETRY2D::JPosition2D; /** * Marker with text. */ struct JPoint_t : public JPosition2D { /** * Constructor. * * \param title title * \param pos position * \param marker marker attributes * \param text text attributes * \param offset text ofset */ JPoint_t(const std::string& title, const JPosition2D& pos, const TAttMarker& marker, const TAttText& text, const double offset) : JPosition2D(pos), offset(offset * cos(text.GetTextAngle()), offset * sin(text.GetTextAngle())) { static_cast(this->marker) = marker; static_cast (this->text) = text; this->text.SetTitle(title.c_str()); this->text.SetTextAngle(0.0); // overwrite text attributes this->text.SetTextColor(marker.GetMarkerColor()); // } /** * Add offset. * * \param x x-position * \param y y-position */ void add(const double x, const double y) { static_cast(*this).add(JPosition2D(x,y)); configure(); } /** * Subtract offset. * * \param x x-position * \param y y-position */ void sub(const double x, const double y) { static_cast(*this).sub(JPosition2D(x,y)); configure(); } /** * Draw. */ void Draw() { marker.Draw(); text .Draw(); } /** * Configure. */ void configure() { this->marker.SetX(this->getX()); this->marker.SetY(this->getY()); this->text .SetX(this->getX() + offset.getX()); this->text .SetY(this->getY() + offset.getY()); } TMarker marker; TText text; private: JPosition2D offset; }; /** * Container of markers. */ struct JGraph_t : public std::vector { /** * Add offset. * * \param x x-position * \param y y-position */ void add(const double x, const double y) { for (iterator i = begin(); i != end(); ++i) { i->add(x, y); } } /** * Subtract offset. * * \param x x-position * \param y y-position */ void sub(const double x, const double y) { for (iterator i = begin(); i != end(); ++i) { i->sub(x, y); } } /** * Draw. */ void Draw() { for (iterator i = begin(); i != end(); ++i) { i->Draw(); } } }; } /** * \file * Auxiliary program to draw the footprint of detector(s). * \author mdejong */ int main(int argc, char**argv) { using namespace std; using namespace JPP; typedef JContainer< vector > tripods_container; typedef JContainer< vector > transmitters_container; vector detectorFile; vector tripodsFile; vector transmittersFile; string outputFile; JCanvas canvas; JCircle2D focus; string legend; double markerSize; double textSize; bool drawCircle; bool all; bool batch; int debug; try { JProperties properties; properties["focus"] = focus; JParser<> zap("Auxiliary program to draw the footprint of detector(s)."); zap['w'] = make_field(canvas, "size of canvas x [pixels]") = JCanvas(500, 500); zap['@'] = make_field(properties, "properties") = JPARSER::initialised(); zap['a'] = make_field(detectorFile, "detector file") = JPARSER::initialised(); zap['T'] = make_field(tripodsFile, "tripod file") = JPARSER::initialised(); zap['Y'] = make_field(transmittersFile, "transmitter file") = JPARSER::initialised(); zap['o'] = make_field(outputFile, "graphics output") = ""; zap['L'] = make_field(legend, "position legend e.g. TR") = "", "TL", "TR", "BR", "BL"; zap['S'] = make_field(markerSize, "marker size") = 1.0; zap['s'] = make_field(textSize, "text size") = 0.02; zap['C'] = make_field(drawCircle, "draw smallest enclosing cicrle"); zap['A'] = make_field(all, "draw all modules"); zap['B'] = make_field(batch, "batch processing"); zap['d'] = make_field(debug) = 1; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } if (detectorFile.empty() && tripodsFile.empty()) { FATAL("No detector elements." << endl); } gROOT->SetBatch(batch); gErrorIgnoreLevel = kWarning; TApplication* tp = new TApplication("user", NULL, NULL); TCanvas* cv = new TCanvas("detector", "", canvas.x, canvas.y); JSinglePointer gStyle(new JStyle("gplot", cv->GetWw(), cv->GetWh())); gROOT->SetStyle("gplot"); gROOT->ForceStyle(); cv->SetFillStyle(4000); cv->SetFillColor(kWhite); cv->Divide(1, 1); cv->cd(1); JMarkerAttributes::getInstance().setMarkerSize(markerSize); vector text_attributes = { TAttText(kHAlignLeft + kVAlignBottom, 0.25*PI, kBlack, 62, textSize), TAttText(kHAlignRight + kVAlignBottom, 0.75*PI, kBlack, 62, textSize), TAttText(kHAlignRight + kVAlignTop, 1.25*PI, kBlack, 62, textSize), TAttText(kHAlignLeft + kVAlignTop, 1.75*PI, kBlack, 62, textSize) }; map data; // graphics data JUTMPosition position; // UTM position double offset = 0.0; // text offset for (size_t i = 0; i != detectorFile.size(); ++i) { JDetector detector; try { load(detectorFile[i], detector); } catch(const JException& error) { FATAL(error); } position = detector.getUTMPosition(); const TAttMarker& marker = *JMarkerAttributes::getInstance().next(); const TAttText& text = text_attributes[(0+i)%text_attributes.size()]; vector& buffer = data[getFilename(detectorFile[i])]; if (offset == 0.0) { offset = 3.0 * textSize * JCircle2D(detector.begin(), detector.end()).getRadius(); } set counter; for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { if (module->getFloor() == 0) { buffer.push_back(JPoint_t(MAKE_STRING(module->getString()), JPosition2D(module->getX(), module->getY()), marker, text, offset)); counter.insert(module->getString()); } } for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { const bool status = (counter.count(module->getString()) == 0); if (status || all) { buffer.push_back(JPoint_t((status ? MAKE_STRING(module->getString()) : ""), JPosition2D(module->getX(), module->getY()), marker, text, offset)); counter.insert(module->getString()); } } } for (size_t i = 0; i != tripodsFile.size(); ++i) { tripods_container tripods; tripods.load(tripodsFile[i].c_str()); const TAttMarker& marker = *JMarkerAttributes::getInstance().next(); const TAttText& text = text_attributes[(0+i)%text_attributes.size()]; vector& buffer = data[getFilename(tripodsFile[i])]; if (offset == 0.0) { offset = 3.0 * textSize * JCircle2D(tripods.begin(), tripods.end()).getRadius(); } for (tripods_container::iterator i = tripods.begin(); i != tripods.end(); ++i) { buffer.push_back(JPoint_t(MAKE_STRING(i->getID()), JPosition2D(i->getUTMEast() - position.getUTMEast(), i->getUTMNorth() - position.getUTMNorth()), marker, text, offset)); } } if (!transmittersFile.empty()) { if (transmittersFile.size() == detectorFile.size()) { for (size_t i = 0; i != transmittersFile.size(); ++i) { transmitters_container transmitters; transmitters.load(transmittersFile[i].c_str()); JDetector detector; load(detectorFile[i], detector); const TAttMarker& marker = *JMarkerAttributes::getInstance().next(); const TAttText& text = text_attributes[(2+i)%text_attributes.size()]; vector& buffer = data[getFilename(transmittersFile[i])]; if (offset == 0.0) { offset = 3.0 * textSize * JCircle2D(transmitters.begin(), transmitters.end()).getRadius(); } for (transmitters_container::iterator i = transmitters.begin(); i != transmitters.end(); ++i) { const JPosition3D pos = detector.getModule(i->getLocation()).getPosition(); buffer.push_back(JPoint_t(MAKE_STRING(i->getID()), JPosition2D(pos.getX() + i->getX(), pos.getY() + i->getY()), marker, text, offset)); } } } else { FATAL("Tranmitter files and detector files should match one-to-one." << endl); } } vector buffer; for (map::iterator i = data.begin(); i != data.end(); ++i) { copy(i->second.begin(), i->second.end(), back_inserter(buffer)); } JCircle2D circle(buffer.begin(), buffer.end()); // enclosing circle if (focus.getRadius() > 0.0) { circle = focus; } NOTICE("focus = " << FIXED(12,3) << circle.getX() << ' ' << FIXED(12,3) << circle.getY() << ' ' << FIXED(9,3) << circle.getRadius() << endl); // center for (map::iterator i = data.begin(); i != data.end(); ++i) { i->second.sub(circle.getX(), circle.getY()); } circle.sub(circle.getPosition()); // draw cv->cd(1); const Double_t xmin = circle.getX() - 1.15 * circle.getRadius(); const Double_t xmax = circle.getX() + 1.15 * circle.getRadius(); const Double_t ymin = circle.getY() - 1.15 * circle.getRadius(); const Double_t ymax = circle.getY() + 1.15 * circle.getRadius(); TH2D h2("h2", "", 200, xmin, xmax, 200, ymin, ymax); h2.GetXaxis()->SetTitle("x [m]"); h2.GetYaxis()->SetTitle("y [m]"); h2.GetXaxis()->CenterTitle(true); h2.GetYaxis()->CenterTitle(true); h2.SetStats(kFALSE); h2.Draw("AXIS"); TEllipse ellipse(circle.getX(), circle.getY(), circle.getRadius()); if (drawCircle) { ellipse.Draw(); } for (map::iterator i = data.begin(); i != data.end(); ++i) { i->second.Draw(); } if (legend != "") { Ssiz_t height = data.size(); Ssiz_t width = 1; for (map::const_iterator i = data.begin(); i != data.end(); ++i) { width = max(width, (Ssiz_t) i->first.size()); } TLegend* lg = getLegend(width, height, legend); lg->SetTextSize(textSize); for (map::const_iterator i = data.begin(); i != data.end(); ++i) { if (!i->second.empty()) { lg->AddEntry(&i->second[0].marker, i->first.c_str(), "P"); } } lg->Draw(); } cv->Update(); if (outputFile != "") { cv->SaveAs(outputFile.c_str()); } if (!batch) { tp->Run(); } }