CMS 3D CMS Logo

MuonGeometrySVGTemplate.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MuonGeometrySVGTemplate
4 // Class: MuonGeometrySVGTemplate
5 //
13 //
14 // Original Author: Jim Pivarski
15 // Created: Wed Jul 14 18:31:18 CDT 2010
16 //
17 //
18 
19 // system include files
20 #include <fstream>
21 
31 
41 
42 //
43 // class decleration
44 //
45 
47 public:
48  explicit MuonGeometrySVGTemplate(const edm::ParameterSet &iConfig);
49  ~MuonGeometrySVGTemplate() override;
50 
51 private:
52  void analyze(const edm::Event &, const edm::EventSetup &iConfig) override;
53 
55  // std::string m_disk1TemplateName;
56  // std::string m_disk23TemplateName;
57  // std::string m_diskp4TemplateName;
58  // std::string m_diskm4TemplateName;
59 };
60 
61 //
62 // constructors and destructor
63 //
64 
66  : m_wheelTemplateName(iConfig.getParameter<std::string>("wheelTemplateName")) {}
67 
69 
70 // ------------ method called to for each event ------------
72  // loads ideal geometry
73  MuonAlignmentInputMethod inputMethod;
74  MuonAlignment muonAlignment(iSetup, inputMethod);
75  AlignableNavigator *alignableNavigator = muonAlignment.getAlignableNavigator();
76 
77  edm::FileInPath fip_BEGINNING("Alignment/MuonAlignment/data/wheel_template.svg_BEGINNING");
78  std::ifstream in_BEGINNING(fip_BEGINNING.fullPath().c_str());
79  edm::FileInPath fip_END("Alignment/MuonAlignment/data/wheel_template.svg_END");
80  std::ifstream in_END(fip_END.fullPath().c_str());
81 
82  const double height = 45.; // assume all chambers are 45 cm tall (local z)
83  std::ofstream out(m_wheelTemplateName.c_str());
84 
85  while (in_BEGINNING.good()) {
86  char c = (char)in_BEGINNING.get();
87  if (in_BEGINNING.good())
88  out << c;
89  }
90 
91  for (int station = 1; station <= 4; station++) {
92  int numSectors = 12;
93  if (station == 4)
94  numSectors = 14;
95  for (int sector = 1; sector <= numSectors; sector++) {
96  DTChamberId id(-2, station, sector); // wheel -2 has a +1 signConvention for x
97  Alignable *chamber = &*(alignableNavigator->alignableFromDetId(id));
98 
99  // different stations, sectors have different widths (*very* fortunate that Alignment software provides this)
100  double width = chamber->surface().width();
101 
102  // lower-left corner of chamber in the chamber's coordinates
103  double x = -width / 2.;
104  double y = -height / 2.;
105 
106  // phi position of chamber
107  align::GlobalVector direction = chamber->surface().toGlobal(LocalVector(1., 0., 0.));
108  double phi = atan2(direction.y(), direction.x());
109 
110  // we'll apply a translation to put the chamber in its place
111  double tx = chamber->surface().position().x();
112  double ty = chamber->surface().position().y();
113 
114  out << " <rect id=\"MB_" << station << "_" << sector << "\" x=\"" << x << "\" y=\"" << y << "\" width=\""
115  << width << "\" height=\"" << height << "\" transform=\"translate(" << tx << ", " << ty << ") rotate("
116  << phi * 180. / M_PI
117  << ")\" style=\"fill:#e1e1e1;fill-opacity:1;stroke:#000000;stroke-width:5.0;stroke-dasharray:1, "
118  "1;stroke-dashoffset:0\" />"
119  << std::endl;
120  }
121  }
122 
123  while (in_END.good()) {
124  char c = (char)in_END.get();
125  if (in_END.good())
126  out << c;
127  }
128 }
129 
130 //define this as a plug-in
Local3DVector LocalVector
Definition: LocalVector.h:12
MuonGeometrySVGTemplate(const edm::ParameterSet &iConfig)
std::string fullPath() const
Definition: FileInPath.cc:161
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define M_PI
void analyze(const edm::Event &, const edm::EventSetup &iConfig) override
~MuonGeometrySVGTemplate() override
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.