CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/SimTracker/TrackerMaterialAnalysis/plugins/ListGroups.cc

Go to the documentation of this file.
00001 #include <string>
00002 #include <vector>
00003 #include <iostream>
00004 #include <iomanip>
00005 
00006 #include "FWCore/Framework/interface/EDAnalyzer.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "FWCore/Framework/interface/ESTransientHandle.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "FWCore/ParameterSet/interface/types.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 #include "DataFormats/DetId/interface/DetId.h"
00016 #include "DataFormats/Math/interface/Vector3D.h"
00017 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00018 #include "DetectorDescription/Core/interface/DDCompactView.h"
00019 #include "DetectorDescription/Core/interface/DDMaterial.h"
00020 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00021 
00022 static
00023 bool dddGetStringRaw(const DDFilteredView & view, const std::string & name, std::string & value) {
00024   std::vector<const DDsvalues_type *> result;
00025   view.specificsV(result);
00026   for (std::vector<const DDsvalues_type *>::iterator it = result.begin(); it != result.end(); ++it)   {
00027     DDValue parameter(name);
00028     if (DDfetch(*it, parameter)) {
00029       if (parameter.strings().size() == 1) {
00030         value = parameter.strings().front();
00031         return true;
00032       } else {
00033         throw cms::Exception("Configuration")<< " ERROR: multiple " << name << " tags encountered";
00034         return false;
00035       }
00036     }
00037   }
00038   return false;
00039 }
00040 
00041 static inline
00042 double dddGetDouble(const std::string & s, DDFilteredView const & view) {
00043   std::string value;
00044   if (dddGetStringRaw(view, s, value))
00045     return double(::atof(value.c_str()));
00046   else
00047     return NAN;
00048 }
00049 
00050 static inline
00051 std::string dddGetString(const std::string & s, DDFilteredView const & view) {
00052   std::string value;
00053   if (dddGetStringRaw(view, s, value))
00054     return value;
00055   else
00056     return std::string();
00057 }
00058 
00059 static inline 
00060 std::ostream & operator<<(std::ostream & out, const math::XYZVector & v) {
00061   out << std::fixed << std::setprecision(3);
00062   return out << "(" << v.rho() << ", " << v.z() << ", " << v.phi() << ")";
00063 }
00064 
00065 class ListGroups : public edm::EDAnalyzer
00066 {
00067 public:
00068   ListGroups(const edm::ParameterSet &);
00069   virtual ~ListGroups();
00070 
00071 private:
00072   void analyze(const edm::Event &, const edm::EventSetup &);
00073   void beginJob() {}
00074   void endJob();
00075 };
00076 
00077 ListGroups::ListGroups(const edm::ParameterSet &) {
00078 }
00079 
00080 ListGroups::~ListGroups() {
00081 }
00082 
00083 void
00084 ListGroups::analyze(const edm::Event& evt, const edm::EventSetup& setup) {
00085   edm::ESTransientHandle<DDCompactView> hDdd;
00086   setup.get<IdealGeometryRecord>().get( hDdd );
00087   DDFilteredView fv(*hDdd);
00088 
00089   DDSpecificsFilter filter;
00090   filter.setCriteria(DDValue("TrackingMaterialGroup", ""), DDSpecificsFilter::not_equals);
00091   fv.addFilter(filter);
00092 
00093   while (fv.next()) {
00094     // print the group name and full hierarchy of all items
00095     std::cout << dddGetString("TrackingMaterialGroup", fv) << '\t';
00096 
00097     // start from 2 to skip the leading /OCMS[0]/CMSE[1] part
00098     const DDGeoHistory & history = fv.geoHistory();
00099     std::cout << '/';
00100     for (unsigned int h = 2; h < history.size(); ++h)
00101       std::cout << '/' << history[h].logicalPart().name().name() << '[' << history[h].copyno() << ']';
00102 
00103     // DD3Vector and DDTranslation are the same type as math::XYZVector
00104     math::XYZVector position = fv.translation() / 10.;  // mm -> cm
00105     std::cout << "\t" << position << std::endl;
00106   };
00107   std::cout << std::endl;
00108 }
00109 
00110 void
00111 ListGroups::endJob() {
00112 }
00113 
00114 //-------------------------------------------------------------------------
00115 // define as a plugin
00116 #include "FWCore/PluginManager/interface/ModuleDef.h"
00117 #include "FWCore/Framework/interface/MakerMacros.h"
00118 DEFINE_FWK_MODULE(ListGroups);