File indexing completed on 2024-04-06 12:15:17
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <iostream>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/EventSetup.h"
0029 #include "FWCore/Framework/interface/ESHandle.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0034 #include "Geometry/MTDNumberingBuilder/interface/GeometricTimingDet.h"
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036
0037 #include "DataFormats/ForwardDetId/interface/MTDDetId.h"
0038
0039 #include "DataFormats/Math/interface/Rounding.h"
0040
0041
0042
0043 #include "DetectorDescription/DDCMS/interface/DDTranslation.h"
0044 #include "DetectorDescription/DDCMS/interface/DDRotationMatrix.h"
0045
0046
0047
0048
0049
0050
0051 class GeometricTimingDetAnalyzer : public edm::one::EDAnalyzer<> {
0052 public:
0053 explicit GeometricTimingDetAnalyzer(const edm::ParameterSet&);
0054 ~GeometricTimingDetAnalyzer() override;
0055
0056 void beginJob() override {}
0057 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0058 void endJob() override {}
0059 void dumpGeometricTimingDet(const GeometricTimingDet* det);
0060
0061 private:
0062 edm::ESGetToken<GeometricTimingDet, IdealGeometryRecord> gtdToken_;
0063 };
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076 GeometricTimingDetAnalyzer::GeometricTimingDetAnalyzer(const edm::ParameterSet& iConfig) {
0077 gtdToken_ = esConsumes<GeometricTimingDet, IdealGeometryRecord>();
0078 }
0079
0080 GeometricTimingDetAnalyzer::~GeometricTimingDetAnalyzer() {
0081
0082
0083 }
0084
0085
0086
0087
0088
0089 using cms_rounding::roundIfNear0;
0090
0091
0092 void GeometricTimingDetAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0093 using namespace edm;
0094
0095 edm::LogInfo("GeometricTimingDetAnalyzer") << "Beginning MTD GeometricTimingDet container dump ";
0096
0097
0098
0099 auto rDD = iSetup.getTransientHandle(gtdToken_);
0100
0101 if (!rDD.isValid()) {
0102 edm::LogError("GeometricTimingDetAnalyzer") << "ESTransientHandle<GeometricTimingDet> rDD is not valid!";
0103 return;
0104 }
0105
0106 edm::LogVerbatim("GeometricTimingDetAnalyzer")
0107 << "\n Top node is " << rDD.product() << " " << rDD.product()->name() << std::endl;
0108
0109 const auto& top = rDD.product();
0110 dumpGeometricTimingDet(top);
0111
0112 edm::LogVerbatim("GeometricTimingDetAnalyzer") << " SubDetectors and layers:";
0113
0114 std::vector<const GeometricTimingDet*> det;
0115
0116 det = rDD->components();
0117 for (const auto& it : det) {
0118 dumpGeometricTimingDet(it);
0119
0120 std::vector<const GeometricTimingDet*> layer;
0121 layer = it->components();
0122 for (const auto& lay : layer) {
0123 dumpGeometricTimingDet(lay);
0124 }
0125 }
0126 det.clear();
0127
0128 edm::LogVerbatim("GeometricTimingDetAnalyzer") << " Modules:";
0129
0130 det = rDD->deepComponents();
0131 for (const auto& it : det) {
0132 dumpGeometricTimingDet(it);
0133 }
0134 }
0135
0136 void GeometricTimingDetAnalyzer::dumpGeometricTimingDet(const GeometricTimingDet* det) {
0137 const GeometricTimingDet::Translation& trans = det->translation();
0138
0139 const GeometricTimingDet::RotationMatrix& res = det->rotation();
0140
0141 DD3Vector x, y, z;
0142 res.GetComponents(x, y, z);
0143
0144 MTDDetId thisDet(det->geographicalID());
0145
0146 auto fround = [&](double in) {
0147 std::stringstream ss;
0148 ss << std::fixed << std::setw(14) << roundIfNear0(in);
0149 return ss.str();
0150 };
0151
0152 edm::LogVerbatim("GeometricTimingDetAnalyzer").log([&](auto& log) {
0153 log << "\n---------------------------------------------------------------------------------------\n";
0154 log << "Module = " << det->name() << " rawId = " << det->geographicalID().rawId()
0155 << " Sub/side/RR = " << thisDet.mtdSubDetector() << " " << thisDet.mtdSide() << " " << thisDet.mtdRR() << "\n\n"
0156 << " type = " << det->type() << "\n"
0157 << " shape = " << det->shape() << "\n\n"
0158 << " radLength " << det->radLength() << "\n"
0159 << " xi " << det->xi() << "\n"
0160 << " PixelROCRows " << det->pixROCRows() << "\n"
0161 << " PixROCCols " << det->pixROCCols() << "\n"
0162 << " PixelROC_X " << det->pixROCx() << "\n"
0163 << " PixelROC_Y " << det->pixROCy() << "\n"
0164 << " TrackerStereoDetectors " << (det->stereo() ? "true" : "false") << "\n"
0165 << " SiliconAPVNumber " << det->siliconAPVNum() << "\n\n"
0166 << " Siblings numbers = ";
0167 std::vector<int> nv = det->navType();
0168 for (auto sib : nv) {
0169 log << sib << ", ";
0170 }
0171 log << "\n"
0172 << " And Contains SubDetectors: " << det->components().size() << "\n"
0173 << " And Contains Daughters: " << det->deepComponents().size() << "\n"
0174 << " Is leaf: " << det->isLeaf() << "\n\n"
0175 << " Translation = " << fround(trans.X()) << " " << fround(trans.Y()) << " " << fround(trans.Z()) << "\n"
0176 << " Rotation = " << fround(x.X()) << " " << fround(x.Y()) << " " << fround(x.Z()) << " " << fround(y.X())
0177 << " " << fround(y.Y()) << " " << fround(y.Z()) << " " << fround(z.X()) << " " << fround(z.Y()) << " "
0178 << fround(z.Z()) << "\n"
0179 << " Phi = " << fround(det->phi()) << " Rho = " << fround(det->rho()) << "\n";
0180 log << "\n---------------------------------------------------------------------------------------\n";
0181 });
0182
0183 edm::LogVerbatim("MTDUnitTest") << det->geographicalID().rawId() << fround(trans.X()) << fround(trans.Y())
0184 << fround(trans.Z()) << fround(x.X()) << fround(x.Y()) << fround(x.Z())
0185 << fround(y.X()) << fround(y.Y()) << fround(y.Z()) << fround(z.X()) << fround(z.Y())
0186 << fround(z.Z());
0187
0188 DD3Vector colx(x.X(), x.Y(), x.Z());
0189 DD3Vector coly(y.X(), y.Y(), y.Z());
0190 DD3Vector colz(z.X(), z.Y(), z.Z());
0191 DDRotationMatrix result(colx, coly, colz);
0192 DD3Vector cx, cy, cz;
0193 result.GetComponents(cx, cy, cz);
0194 if (cx.Cross(cy).Dot(cz) < 0.5) {
0195 edm::LogInfo("GeometricTimingDetAnalyzer")
0196 << "Left Handed Rotation Matrix detected; making it right handed: " << det->name();
0197 }
0198 }
0199
0200
0201 DEFINE_FWK_MODULE(GeometricTimingDetAnalyzer);