00001 #include "Validation/MuonRPCDigis/interface/RPCDigiValid.h"
00002 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00003 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
00004 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00005 #include "DataFormats/RPCDigi/interface/RPCDigiCollection.h"
00006 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00007 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
00008 #include "DataFormats/RPCDigi/interface/RPCDigiCollection.h"
00009
00010 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
00011 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
00012 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00013 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00014 #include "DQMServices/Core/interface/DQMStore.h"
00015
00016 using namespace std;
00017 using namespace edm;
00018
00019 RPCDigiValid::RPCDigiValid(const ParameterSet& ps):dbe_(0){
00020
00021 digiLabel = ps.getUntrackedParameter<std::string>("digiLabel");
00022 outputFile_ = ps.getUntrackedParameter<string>("outputFile", "rpcDigiValidPlots.root");
00023 dbe_ = Service<DQMStore>().operator->();
00024
00025 if ( dbe_ ) {
00026 dbe_->setCurrentFolder("RPCDigisV/RPCDigis");
00027
00028 xyview = dbe_->book2D("X_Vs_Y_View","X_Vs_Y_View",1000, -760., 760., 1000, -760., 760.);
00029
00030 xyvDplu1 = dbe_->book2D("Dplu1_XvsY","Dplu1_XvsY",1000, -760., 760., 1000, -760., 760.);
00031 xyvDplu2 = dbe_->book2D("Dplu2_XvsY","Dplu2_XvsY",1000, -760., 760., 1000, -760., 760.);
00032 xyvDplu3 = dbe_->book2D("Dplu3_XvsY","Dplu3_XvsY",1000, -760., 760., 1000, -760., 760.);
00033 xyvDmin1 = dbe_->book2D("Dmin1_XvsY","Dmin1_XvsY",1000, -760., 760., 1000, -760., 760.);
00034 xyvDmin2 = dbe_->book2D("Dmin2_XvsY","Dmin2_XvsY",1000, -760., 760., 1000, -760., 760.);
00035 xyvDmin3 = dbe_->book2D("Dmin3_XvsY","Dmin3_XvsY",1000, -760., 760., 1000, -760., 760.);
00036
00037 rzview = dbe_->book2D("R_Vs_Z_View","R_Vs_Z_View",1000, -1100., 1100.,1000,0., 800.);
00038 Res = dbe_->book1D("Digi_SimHit_difference", "Digi_SimHit_difference", 300, -8, 8);
00039 ResWmin2 = dbe_->book1D("W_Min2_Residuals", "W_Min2_Residuals", 300, -8, 8);
00040 ResWmin1 = dbe_->book1D("W_Min1_Residuals", "W_Min1_Residuals", 300, -8, 8);
00041 ResWzer0 = dbe_->book1D("W_Zer0_Residuals", "W_Zer0_Residuals", 300, -8, 8);
00042 ResWplu1 = dbe_->book1D("W_Plu1_Residuals", "W_Plu1_Residuals", 300, -8, 8);
00043 ResWplu2 = dbe_->book1D("W_Plu2_Residuals", "W_Plu2_Residuals", 300, -8, 8);
00044
00045 BxDist = dbe_->book1D("Bunch_Crossing", "Bunch_Crossing", 20, -10., 10.);
00046 StripProf = dbe_->book1D("Strip_Profile", "Strip_Profile", 100, 0, 100);
00047
00048 ResDmin1 = dbe_->book1D("Disk_Min1_Residuals", "Disk_Min1_Residuals", 300, -8, 8);
00049 ResDmin2 = dbe_->book1D("Disk_Min2_Residuals", "Disk_Min2_Residuals", 300, -8, 8);
00050 ResDmin3 = dbe_->book1D("Disk_Min3_Residuals", "Disk_Min3_Residuals", 300, -8, 8);
00051 ResDplu1 = dbe_->book1D("Disk_Plu1_Residuals", "Disk_Plu1_Residuals", 300, -8, 8);
00052 ResDplu2 = dbe_->book1D("Disk_Plu2_Residuals", "Disk_Plu2_Residuals", 300, -8, 8);
00053 ResDplu3 = dbe_->book1D("Disk_Plu3_Residuals", "Disk_Plu3_Residuals", 300, -8, 8);
00054
00055 }
00056 }
00057
00058 RPCDigiValid::~RPCDigiValid(){}
00059
00060 void RPCDigiValid::beginJob(){}
00061
00062 void RPCDigiValid::endJob() {
00063 if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00064 }
00065
00066 void RPCDigiValid::analyze(const Event& event, const EventSetup& eventSetup){
00067
00068
00069
00070
00071
00072 edm::ESHandle<RPCGeometry> rpcGeom;
00073 eventSetup.get<MuonGeometryRecord>().get(rpcGeom);
00074
00075 edm::Handle<PSimHitContainer> simHit;
00076 event.getByLabel("g4SimHits", "MuonRPCHits", simHit);
00077
00078 edm::Handle<RPCDigiCollection> rpcDigis;
00079 event.getByLabel(digiLabel, rpcDigis);
00080
00081
00082 PSimHitContainer::const_iterator simIt;
00083
00084
00085 std::map<RPCDetId, std::vector<double> > allsims;
00086
00087 for (simIt = simHit->begin(); simIt != simHit->end(); simIt++) {
00088 RPCDetId Rsid = (RPCDetId)(*simIt).detUnitId();
00089 const RPCRoll* soll = dynamic_cast<const RPCRoll* >( rpcGeom->roll(Rsid));
00090 int ptype = simIt->particleType();
00091
00092
00093 if (ptype == 13 || ptype == -13) {
00094 std::vector<double> buff;
00095 if (allsims.find(Rsid) != allsims.end() ){
00096 buff= allsims[Rsid];
00097 }
00098 buff.push_back(simIt->localPosition().x());
00099 allsims[Rsid]=buff;
00100 }
00101 GlobalPoint p=soll->toGlobal(simIt->localPosition());
00102
00103
00104
00105
00106
00107
00108 double sim_x = p.x();
00109 double sim_y = p.y();
00110
00111 xyview->Fill(sim_x, sim_y);
00112
00113 if (Rsid.region() == 0 ){
00114 }
00115 else if (Rsid.region() == (+1)){
00116 if (Rsid.station() == 1)
00117 xyvDplu1->Fill(sim_x, sim_y);
00118 else if (Rsid.station() == 2)
00119 xyvDplu2->Fill(sim_x, sim_y);
00120 else if (Rsid.station() == 3)
00121 xyvDplu3->Fill(sim_x, sim_y);
00122 }
00123 else if (Rsid.region() == (-1)){
00124 if (Rsid.station() == 1)
00125 xyvDmin1->Fill(sim_x, sim_y);
00126 else if (Rsid.station() == 2)
00127 xyvDmin2->Fill(sim_x, sim_y);
00128 else if (Rsid.station() == 3)
00129 xyvDmin3->Fill(sim_x, sim_y);
00130 }
00131
00132
00133
00134 rzview->Fill(p.z(),p.perp());
00135
00136 }
00137
00138 RPCDigiCollection::DigiRangeIterator detUnitIt;
00139 for (detUnitIt=rpcDigis->begin(); detUnitIt!=rpcDigis->end();++detUnitIt){
00140 const RPCDetId Rsid = (*detUnitIt).first;
00141 const RPCRoll* roll = dynamic_cast<const RPCRoll* >( rpcGeom->roll(Rsid));
00142 const RPCDigiCollection::Range& range = (*detUnitIt).second;
00143 std::vector<double> sims;
00144 if (allsims.find(Rsid) != allsims.end() ){
00145 sims = allsims[Rsid];
00146 }
00147 int ndigi=0;
00148 for (RPCDigiCollection::const_iterator digiIt = range.first;
00149 digiIt != range.second; ++digiIt){
00150 StripProf->Fill(digiIt->strip());
00151 BxDist->Fill(digiIt->bx());
00152 ndigi++;
00153 }
00154
00155
00156
00157
00158 if (sims.size() == 1 && ndigi == 1){
00159 double dis = roll->centreOfStrip(range.first->strip()).x()-sims[0];
00160 Res->Fill(dis);
00161
00162 if (Rsid.region() == 0 ){
00163 if (Rsid.ring() == -2)
00164 ResWmin2->Fill(dis);
00165 else if (Rsid.ring() == -1)
00166 ResWmin1->Fill(dis);
00167 else if (Rsid.ring() == 0)
00168 ResWzer0->Fill(dis);
00169 else if (Rsid.ring() == 1)
00170 ResWplu1->Fill(dis);
00171 else if (Rsid.ring() == 2)
00172 ResWplu2->Fill(dis);
00173 }
00174 else if (Rsid.region() == (+1)){
00175 if (Rsid.station() == 1)
00176 ResDplu1->Fill(dis);
00177 else if (Rsid.station() == 2)
00178 ResDplu2->Fill(dis);
00179 else if (Rsid.station() == 3)
00180 ResDplu3->Fill(dis);
00181 }
00182 else if (Rsid.region() == (-1)){
00183 if (Rsid.station() == 1)
00184 ResDmin1->Fill(dis);
00185 else if (Rsid.station() == 2)
00186 ResDmin2->Fill(dis);
00187 else if (Rsid.station() == 3)
00188 ResDmin3->Fill(dis);
00189 }
00190 }
00191 }
00192 }