#include <RPCDigiValid.h>
Public Member Functions | |
RPCDigiValid (const edm::ParameterSet &ps) | |
~RPCDigiValid () | |
Protected Member Functions | |
void | analyze (const edm::Event &e, const edm::EventSetup &c) |
void | beginJob () |
void | endJob (void) |
Private Attributes | |
MonitorElement * | BxDist |
DQMStore * | dbe_ |
std::string | digiLabel |
std::string | outputFile_ |
MonitorElement * | Res |
MonitorElement * | ResWmin1 |
MonitorElement * | ResWmin2 |
MonitorElement * | ResWplu1 |
MonitorElement * | ResWplu2 |
MonitorElement * | ResWzer0 |
MonitorElement * | rzview |
MonitorElement * | StripProf |
MonitorElement * | xyview |
Definition at line 24 of file RPCDigiValid.h.
RPCDigiValid::RPCDigiValid | ( | const edm::ParameterSet & | ps | ) |
Definition at line 19 of file RPCDigiValid.cc.
References BxDist, dbe_, digiLabel, edm::ParameterSet::getUntrackedParameter(), cmsCodeRules::cppFunctionSkipper::operator, outputFile_, Res, ResWmin1, ResWmin2, ResWplu1, ResWplu2, ResWzer0, rzview, StripProf, and xyview.
:dbe_(0){ digiLabel = ps.getUntrackedParameter<std::string>("digiLabel"); outputFile_ = ps.getUntrackedParameter<string>("outputFile", "rpcDigiValidPlots.root"); dbe_ = Service<DQMStore>().operator->(); if ( dbe_ ) { dbe_->setCurrentFolder("RPCDigisV/RPCDigis"); xyview = dbe_->book2D("X Vs Y View","X Vs Y View",1000, -700., 700., 1000, -700., 700.); rzview = dbe_->book2D("R Vs Z View","X Vs Y View",1000, -1100., 1100.,1000,0., 700.); Res = dbe_->book1D("Digi SimHit difference", "Digi SimHit Difference", 300, -8, 8); ResWmin2 = dbe_->book1D("W-2 Residuals", "Residuals for Wheel -2", 300, -8, 8); ResWmin1 = dbe_->book1D("W-1 Residuals", "Residuals for Wheel -1", 300, -8, 8); ResWzer0 = dbe_->book1D("W 0 Residuals", "Residuals for Wheel 0", 300, -8, 8); ResWplu1 = dbe_->book1D("W+1 Residuals", "Residuals for Wheel +1", 300, -8, 8); ResWplu2 = dbe_->book1D("W+2 Residuals", "Residuals for Wheel +2", 300, -8, 8); BxDist = dbe_->book1D("Bunch Crossing", "Bunch Crossing", 20, -9.5, 9.5); StripProf = dbe_->book1D("Strip Profile", "Strip Profile", 100, 0, 100); } }
RPCDigiValid::~RPCDigiValid | ( | ) |
Definition at line 42 of file RPCDigiValid.cc.
{}
void RPCDigiValid::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | c | ||
) | [protected, virtual] |
Implements edm::EDAnalyzer.
Definition at line 50 of file RPCDigiValid.cc.
References BxDist, RPCRoll::centreOfStrip(), digiLabel, MonitorElement::Fill(), edm::EventSetup::get(), L1TEmulatorMonitor_cff::p, PV3DBase< T, PVType, FrameType >::perp(), RPCDetId::region(), Res, ResWmin1, ResWmin2, ResWplu1, ResWplu2, ResWzer0, RPCDetId::ring(), RPCDetId, rzview, StripProf, GeomDet::toGlobal(), ExpressReco_HICollisions_FallBack::x, PV3DBase< T, PVType, FrameType >::x(), xyview, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
{ // cout << endl <<"--- [RPCDigiQuality] Analysing Event: #Run: " << event.id().run() // << " #Event: " << event.id().event() << endl; // Get the RPC Geometry edm::ESHandle<RPCGeometry> rpcGeom; eventSetup.get<MuonGeometryRecord>().get(rpcGeom); edm::Handle<PSimHitContainer> simHit; event.getByLabel("g4SimHits", "MuonRPCHits", simHit); edm::Handle<RPCDigiCollection> rpcDigis; event.getByLabel(digiLabel, rpcDigis); // Loop on simhits PSimHitContainer::const_iterator simIt; //loop over Simhit std::map<RPCDetId, std::vector<double> > allsims; for (simIt = simHit->begin(); simIt != simHit->end(); simIt++) { RPCDetId Rsid = (RPCDetId)(*simIt).detUnitId(); const RPCRoll* soll = dynamic_cast<const RPCRoll* >( rpcGeom->roll(Rsid)); int ptype = simIt->particleType(); // std::cout <<"This is a Simhit with Parent "<<ptype<<std::endl; if (ptype == 13 || ptype == -13) { std::vector<double> buff; if (allsims.find(Rsid) != allsims.end() ){ buff= allsims[Rsid]; } buff.push_back(simIt->localPosition().x()); allsims[Rsid]=buff; } GlobalPoint p=soll->toGlobal(simIt->localPosition()); /* std::cout <<"Muon Position phi="<<p.phi() <<" R="<<p.perp() <<" z="<<p.z()<<std::endl; */ xyview->Fill(p.x(),p.y()); rzview->Fill(p.z(),p.perp()); } //loop over Digis RPCDigiCollection::DigiRangeIterator detUnitIt; for (detUnitIt=rpcDigis->begin(); detUnitIt!=rpcDigis->end();++detUnitIt){ const RPCDetId Rsid = (*detUnitIt).first; const RPCRoll* roll = dynamic_cast<const RPCRoll* >( rpcGeom->roll(Rsid)); const RPCDigiCollection::Range& range = (*detUnitIt).second; std::vector<double> sims; if (allsims.find(Rsid) != allsims.end() ){ sims = allsims[Rsid]; } int ndigi=0; for (RPCDigiCollection::const_iterator digiIt = range.first; digiIt != range.second; ++digiIt){ StripProf->Fill(digiIt->strip()); BxDist->Fill(digiIt->bx()); ndigi++; } // std::cout<<" Number of Digi "<<ndigi<<" for "<<Rsid<<std::endl; if (sims.size() == 1 && ndigi == 1){ double dis = roll->centreOfStrip(range.first->strip()).x()-sims[0]; Res->Fill(dis); if (Rsid.region() == 0 ){ if (Rsid.ring() == -2) ResWmin2->Fill(dis); else if (Rsid.ring() == -1) ResWmin1->Fill(dis); else if (Rsid.ring() == 0) ResWzer0->Fill(dis); else if (Rsid.ring() == 1) ResWplu1->Fill(dis); else if (Rsid.ring() == 2) ResWplu2->Fill(dis); } } } }
void RPCDigiValid::beginJob | ( | void | ) | [protected, virtual] |
void RPCDigiValid::endJob | ( | void | ) | [protected, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 46 of file RPCDigiValid.cc.
References dbe_, outputFile_, and DQMStore::save().
{ if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_); }
MonitorElement* RPCDigiValid::BxDist [private] |
Definition at line 46 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
DQMStore* RPCDigiValid::dbe_ [private] |
Definition at line 49 of file RPCDigiValid.h.
Referenced by endJob(), and RPCDigiValid().
std::string RPCDigiValid::digiLabel [private] |
Definition at line 51 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
std::string RPCDigiValid::outputFile_ [private] |
Definition at line 50 of file RPCDigiValid.h.
Referenced by endJob(), and RPCDigiValid().
MonitorElement* RPCDigiValid::Res [private] |
Definition at line 40 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::ResWmin1 [private] |
Definition at line 42 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::ResWmin2 [private] |
Definition at line 41 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::ResWplu1 [private] |
Definition at line 44 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::ResWplu2 [private] |
Definition at line 45 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::ResWzer0 [private] |
Definition at line 43 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::rzview [private] |
Definition at line 39 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::StripProf [private] |
Definition at line 47 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().
MonitorElement* RPCDigiValid::xyview [private] |
Definition at line 38 of file RPCDigiValid.h.
Referenced by analyze(), and RPCDigiValid().