Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-12-24 02:19:00

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 #include "FWCore/Framework/interface/ESTransientHandle.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "DetectorDescription/Core/interface/DDFilteredView.h"
0007 #include "DetectorDescription/Core/interface/DDCompactView.h"
0008 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0009 
0010 class DDFilteredViewAnalyzer : public edm::one::EDAnalyzer<> {
0011 public:
0012   explicit DDFilteredViewAnalyzer(const edm::ParameterSet&);
0013   ~DDFilteredViewAnalyzer(void) override {}
0014 
0015   void beginJob() override {}
0016   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0017   void endJob() override {}
0018 
0019 private:
0020   std::string m_attribute;
0021   std::string m_value;
0022   bool m_shouldPrint;
0023   DDCompOp m_comp;
0024   edm::ESGetToken<DDCompactView, IdealGeometryRecord> ddToken_;
0025 };
0026 
0027 DDFilteredViewAnalyzer::DDFilteredViewAnalyzer(const edm::ParameterSet& pset) {
0028   m_attribute = pset.getParameter<std::string>("attribute");
0029   m_value = pset.getParameter<std::string>("value");
0030 
0031   m_shouldPrint = pset.getUntrackedParameter<bool>("shouldPrint", true);
0032   if (pset.getUntrackedParameter<bool>("compareNotEquals", true)) {
0033     m_comp = DDCompOp::not_equals;
0034   } else {
0035     m_comp = DDCompOp::equals;
0036   }
0037   ddToken_ = esConsumes<DDCompactView, IdealGeometryRecord>();
0038 }
0039 
0040 void DDFilteredViewAnalyzer::analyze(const edm::Event&, const edm::EventSetup& iSetup) {
0041   edm::ESTransientHandle<DDCompactView> cpv = iSetup.getTransientHandle(ddToken_);
0042 
0043   DDValue val(m_attribute, m_value, 0.0);
0044   DDSpecificsFilter filter;
0045   filter.setCriteria(val,  // name & value of a variable
0046                      m_comp);
0047   DDFilteredView fv(*cpv, filter);
0048   if (fv.firstChild()) {
0049     std::cout << "Found attribute " << m_attribute.c_str() << " with value " << m_value.c_str() << std::endl;
0050     bool dodet = true;
0051     int i = 0;
0052     while (dodet) {
0053       dodet = fv.next();
0054       if (m_shouldPrint) {
0055         std::cout << i++ << ": " << fv.logicalPart().name() << std::endl;
0056       }
0057     }
0058   } else
0059     std::cout << "No luck..." << std::endl;
0060 }
0061 
0062 DEFINE_FWK_MODULE(DDFilteredViewAnalyzer);