#include <DQM/DTMonitorModule/src/DTResolutionAnalysisTask.h>
Public Member Functions | |
void | analyze (const edm::Event &event, const edm::EventSetup &setup) |
void | beginJob (const edm::EventSetup &c) |
BeginJob. | |
void | beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) |
To reset the MEs. | |
DTResolutionAnalysisTask (const edm::ParameterSet &pset) | |
Constructor. | |
void | endJob () |
Endjob. | |
virtual | ~DTResolutionAnalysisTask () |
Destructor. | |
Private Member Functions | |
void | bookHistos (DTSuperLayerId slId) |
void | fillHistos (DTSuperLayerId slId, float distExtr, float residual) |
Private Attributes | |
edm::ESHandle< DTGeometry > | dtGeom |
std::map< DTSuperLayerId, std::vector< MonitorElement * > > | histosPerSL |
int | prescaleFactor |
int | resetCycle |
DQMStore * | theDbe |
std::string | theRecHitLabel |
std::string | theRecHits4DLabel |
Definition at line 32 of file DTResolutionAnalysisTask.h.
DTResolutionAnalysisTask::DTResolutionAnalysisTask | ( | const edm::ParameterSet & | pset | ) |
Constructor.
Definition at line 35 of file DTResolutionAnalysisTask.cc.
References lat::endl(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), prescaleFactor, resetCycle, theRecHitLabel, and theRecHits4DLabel.
00035 { 00036 00037 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] Constructor called!" << endl; 00038 00039 // the name of the 4D rec hits collection 00040 theRecHits4DLabel = pset.getParameter<string>("recHits4DLabel"); 00041 // the name of the rechits collection 00042 theRecHitLabel = pset.getParameter<string>("recHitLabel"); 00043 00044 prescaleFactor = pset.getUntrackedParameter<int>("diagnosticPrescale", 1); 00045 resetCycle = pset.getUntrackedParameter<int>("ResetCycle", -1); 00046 }
DTResolutionAnalysisTask::~DTResolutionAnalysisTask | ( | ) | [virtual] |
Destructor.
Definition at line 49 of file DTResolutionAnalysisTask.cc.
References lat::endl().
00049 { 00050 00051 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] Destructor called!" << endl; 00052 00053 }
void DTResolutionAnalysisTask::analyze | ( | const edm::Event & | event, | |
const edm::EventSetup & | setup | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 108 of file DTResolutionAnalysisTask.cc.
References edmNew::copy(), funct::cos(), dtGeom, lat::endl(), fillHistos(), DTWireId::layerId(), range, sl, DTRecSegment2D::specificRecHits(), DTLayer::specificTopology(), DTSuperLayerId::superlayer(), DTChamber::superLayer(), DTLayerId::superlayerId(), theRecHitLabel, theRecHits4DLabel, GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTTopology::wirePosition(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
00108 { 00109 00110 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] Analyze #Run: " << event.id().run() 00111 << " #Event: " << event.id().event() << endl; 00112 00113 00114 // Get the 4D segment collection from the event 00115 edm::Handle<DTRecSegment4DCollection> all4DSegments; 00116 event.getByLabel(theRecHits4DLabel, all4DSegments); 00117 00118 // Get the rechit collection from the event 00119 Handle<DTRecHitCollection> dtRecHits; 00120 event.getByLabel(theRecHitLabel, dtRecHits); 00121 00122 // Loop over all chambers containing a segment 00123 DTRecSegment4DCollection::id_iterator chamberId; 00124 for (chamberId = all4DSegments->id_begin(); 00125 chamberId != all4DSegments->id_end(); 00126 ++chamberId) { 00127 // Get the range for the corresponding ChamerId 00128 DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId); 00129 int nsegm = distance(range.first, range.second); 00130 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Chamber: " << *chamberId << " has " << nsegm 00131 //<< " 4D segments" << endl; 00132 // Get the chamber 00133 const DTChamber* chamber = dtGeom->chamber(*chamberId); 00134 00135 // Loop over the rechits of this ChamerId 00136 for (DTRecSegment4DCollection::const_iterator segment4D = range.first; 00137 segment4D!=range.second; 00138 ++segment4D) { 00139 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " == RecSegment dimension: " << (*segment4D).dimension() << endl; 00140 00141 // If Statio != 4 skip RecHits with dimension != 4 00142 // For the Station 4 consider 2D RecHits 00143 if((*chamberId).station() != 4 && (*segment4D).dimension() != 4) { 00144 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask]***Warning: RecSegment dimension is not 4 but " 00145 << (*segment4D).dimension() << "!" << endl; 00146 continue; 00147 } else if((*chamberId).station() == 4 && (*segment4D).dimension() != 2) { 00148 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask]***Warning: RecSegment dimension is not 2 but " 00149 << (*segment4D).dimension() << "!" << endl; 00150 continue; 00151 } 00152 00153 00154 // Get all 1D RecHits at step 3 within the 4D segment 00155 vector<DTRecHit1D> recHits1D_S3; 00156 00157 00158 // Get 1D RecHits at Step 3 and select only events with 00159 // 8 hits in phi and 4 hits in theta (if any) 00160 00161 if((*segment4D).hasPhi()) { // has phi component 00162 const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment(); 00163 vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits(); 00164 00165 if(phiRecHits.size() != 8) { 00166 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] Phi segments has: " << phiRecHits.size() 00167 //<< " hits" << endl; // FIXME: info output 00168 continue; 00169 } 00170 copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3)); 00171 } else { 00172 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] 4D segment has not phi component!" << endl; 00173 } 00174 00175 if((*segment4D).hasZed()) { 00176 const DTSLRecSegment2D* zSeg = (*segment4D).zSegment(); 00177 vector<DTRecHit1D> zRecHits = zSeg->specificRecHits(); 00178 if(zRecHits.size() != 4) { 00179 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] Theta segments has: " << zRecHits.size() 00180 //<< " hits, skipping" << endl; // FIXME: info output 00181 continue; 00182 } 00183 copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3)); 00184 } 00185 00186 // Loop over 1D RecHit inside 4D segment 00187 for(vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin(); 00188 recHit1D != recHits1D_S3.end(); 00189 recHit1D++) { 00190 const DTWireId wireId = (*recHit1D).wireId(); 00191 00192 // Get the layer and the wire position 00193 const DTLayer* layer = chamber->superLayer(wireId.superlayerId())->layer(wireId.layerId()); 00194 float wireX = layer->specificTopology().wirePosition(wireId.wire()); 00195 00196 // Distance of the 1D rechit from the wire 00197 float distRecHitToWire = fabs(wireX - (*recHit1D).localPosition().x()); 00198 00199 // Extrapolate the segment to the z of the wire 00200 00201 // Get wire position in chamber RF 00202 LocalPoint wirePosInLay(wireX,0,0); 00203 GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay); 00204 LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob); 00205 00206 // Segment position at Wire z in chamber local frame 00207 LocalPoint segPosAtZWire = (*segment4D).localPosition() 00208 + (*segment4D).localDirection()*wirePosInChamber.z()/cos((*segment4D).localDirection().theta()); 00209 00210 // Compute the distance of the segment from the wire 00211 int sl = wireId.superlayer(); 00212 00213 double distSegmToWire = -1; 00214 if(sl == 1 || sl == 3) { 00215 // RPhi SL 00216 distSegmToWire = fabs(wirePosInChamber.x() - segPosAtZWire.x()); 00217 } else if(sl == 2) { 00218 // RZ SL 00219 distSegmToWire = fabs(wirePosInChamber.y() - segPosAtZWire.y()); 00220 } 00221 00222 if(distSegmToWire > 2.1) 00223 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Warning: dist segment-wire: " << distSegmToWire << endl; 00224 00225 double residual = distRecHitToWire - distSegmToWire; 00226 // FIXME: Fill the histos 00227 fillHistos(wireId.superlayerId(), distSegmToWire, residual); 00228 00229 // edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Dist. segment extrapolation - wire (cm): " << distSegmToWire << endl; 00230 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Dist. RecHit - wire (cm): " << distRecHitToWire << endl; 00231 //edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Residual (cm): " << residual << endl; 00232 00233 00234 }// End of loop over 1D RecHit inside 4D segment 00235 }// End of loop over the rechits of this ChamerId 00236 } 00237 // ----------------------------------------------------------------------------- 00238 }
void DTResolutionAnalysisTask::beginJob | ( | const edm::EventSetup & | c | ) | [virtual] |
BeginJob.
Reimplemented from edm::EDAnalyzer.
Definition at line 56 of file DTResolutionAnalysisTask.cc.
References bookHistos(), muonGeometry::chambers, dtGeom, edm::EventSetup::get(), sl, DTChamberId::station(), and theDbe.
00056 { 00057 // Get the DQM needed services 00058 theDbe = edm::Service<DQMStore>().operator->(); 00059 // theDbe->setCurrentFolder("DT/02-Segments"); 00060 00061 // Get the DT Geometry 00062 setup.get<MuonGeometryRecord>().get(dtGeom); 00063 00064 00065 // Book the histograms 00066 vector<DTChamber*> chambers = dtGeom->chambers(); 00067 for(vector<DTChamber*>::const_iterator chamber = chambers.begin(); 00068 chamber != chambers.end(); ++chamber) { // Loop over all chambers 00069 DTChamberId dtChId = (*chamber)->id(); 00070 for(int sl = 1; sl <= 3; ++sl) { // Loop over SLs 00071 if(dtChId.station() == 4 && sl == 2) continue; 00072 const DTSuperLayerId dtSLId(dtChId,sl); 00073 bookHistos(dtSLId); 00074 } 00075 } 00076 00077 }
void DTResolutionAnalysisTask::beginLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, | |
edm::EventSetup const & | context | |||
) | [virtual] |
To reset the MEs.
Reimplemented from edm::EDAnalyzer.
Definition at line 82 of file DTResolutionAnalysisTask.cc.
References lat::endl(), histo, histosPerSL, i, edm::LuminosityBlock::id(), edm::LuminosityBlockID::luminosityBlock(), resetCycle, and size.
00083 { 00084 00085 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionTask]: Begin of LS transition"<<endl; 00086 00087 if(resetCycle != -1 && lumiSeg.id().luminosityBlock() % resetCycle == 0) { 00088 for(map<DTSuperLayerId, vector<MonitorElement*> > ::const_iterator histo = histosPerSL.begin(); 00089 histo != histosPerSL.end(); 00090 histo++) { 00091 int size = (*histo).second.size(); 00092 for(int i=0; i<size; i++){ 00093 (*histo).second[i]->Reset(); 00094 } 00095 } 00096 } 00097 }
void DTResolutionAnalysisTask::bookHistos | ( | DTSuperLayerId | slId | ) | [private] |
Definition at line 242 of file DTResolutionAnalysisTask.cc.
References DQMStore::book1D(), lat::endl(), combine::histos, histosPerSL, DTChamberId::sector(), DQMStore::setCurrentFolder(), DTChamberId::station(), DTSuperLayerId::superlayer(), theDbe, muonGeometry::wheel, and DTChamberId::wheel().
Referenced by beginJob().
00242 { 00243 00244 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << " Booking histos for SL: " << slId << endl; 00245 00246 // Compose the chamber name 00247 stringstream wheel; wheel << slId.wheel(); 00248 stringstream station; station << slId.station(); 00249 stringstream sector; sector << slId.sector(); 00250 stringstream superLayer; superLayer << slId.superlayer(); 00251 00252 00253 string slHistoName = 00254 "_W" + wheel.str() + 00255 "_St" + station.str() + 00256 "_Sec" + sector.str() + 00257 "_SL" + superLayer.str(); 00258 00259 theDbe->setCurrentFolder("DT/02-Segments/Wheel" + wheel.str() + 00260 "/Sector" + sector.str() + 00261 "/Station" + station.str()); 00262 // Create the monitor elements 00263 vector<MonitorElement *> histos; 00264 // Note hte order matters 00265 histos.push_back(theDbe->book1D("hResDist"+slHistoName, 00266 "Residuals on the distance from wire (rec_hit - segm_extr) (cm)", 00267 200, -0.4, 0.4)); 00268 //FIXME: 2D plot removed to reduce the # of ME 00269 // histos.push_back(theDbe->book2D("hResDistVsDist"+slHistoName, 00270 // "Residuals on the distance (cm) from wire (rec_hit - segm_extr) vs distance (cm)", 00271 // 100, 0, 2.5, 200, -0.4, 0.4)); 00272 histosPerSL[slId] = histos; 00273 }
Endjob.
Reimplemented from edm::EDAnalyzer.
Definition at line 100 of file DTResolutionAnalysisTask.cc.
References lat::endl().
00100 { 00101 00102 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTResolutionAnalysisTask") << "[DTResolutionAnalysisTask] endjob called!"<<endl; 00103 00104 }
void DTResolutionAnalysisTask::fillHistos | ( | DTSuperLayerId | slId, | |
float | distExtr, | |||
float | residual | |||
) | [private] |
Definition at line 277 of file DTResolutionAnalysisTask.cc.
References combine::histos, and histosPerSL.
Referenced by analyze().
00279 { 00280 vector<MonitorElement *> histos = histosPerSL[slId]; 00281 histos[0]->Fill(residual); 00282 //FIXME: 2D plot removed to reduce the # of ME 00283 // histos[1]->Fill(distExtr, residual); 00284 00285 }
std::map<DTSuperLayerId, std::vector<MonitorElement*> > DTResolutionAnalysisTask::histosPerSL [private] |
Definition at line 75 of file DTResolutionAnalysisTask.h.
Referenced by beginLuminosityBlock(), bookHistos(), and fillHistos().
int DTResolutionAnalysisTask::prescaleFactor [private] |
int DTResolutionAnalysisTask::resetCycle [private] |
Definition at line 61 of file DTResolutionAnalysisTask.h.
Referenced by beginLuminosityBlock(), and DTResolutionAnalysisTask().
DQMStore* DTResolutionAnalysisTask::theDbe [private] |
Definition at line 56 of file DTResolutionAnalysisTask.h.
Referenced by beginJob(), and bookHistos().
std::string DTResolutionAnalysisTask::theRecHitLabel [private] |
Definition at line 66 of file DTResolutionAnalysisTask.h.
Referenced by analyze(), and DTResolutionAnalysisTask().
std::string DTResolutionAnalysisTask::theRecHits4DLabel [private] |
Definition at line 64 of file DTResolutionAnalysisTask.h.
Referenced by analyze(), and DTResolutionAnalysisTask().