#include <DTVDriftSegmentCalibration.h>
Classes | |
class | DTVDriftSegmentCalibration |
Public Member Functions | |
void | analyze (const edm::Event &event, const edm::EventSetup &eventSetup) |
void | beginJob () |
void | beginRun (const edm::Run &run, const edm::EventSetup &setup) |
DTVDriftSegmentCalibration (const edm::ParameterSet &pset) | |
void | endJob () |
virtual | ~DTVDriftSegmentCalibration () |
Private Types | |
typedef std::map< DTChamberId, std::vector< TH1F * > > | ChamberHistosMapTH1F |
typedef std::map< DTChamberId, std::vector< TH2F * > > | ChamberHistosMapTH2F |
Private Member Functions | |
void | bookHistos (DTChamberId) |
Private Attributes | |
TFile * | rootFile_ |
DTSegmentSelector | select_ |
std::string | theCalibChamber_ |
edm::InputTag | theRecHits4DLabel_ |
ChamberHistosMapTH1F | theVDriftHistoMapTH1F_ |
ChamberHistosMapTH2F | theVDriftHistoMapTH2F_ |
Produces histograms from v-drift computation in segment fit to be used for v-drift calibration
Definition at line 24 of file DTVDriftSegmentCalibration.h.
typedef std::map<DTChamberId, std::vector<TH1F*> > DTVDriftSegmentCalibration::ChamberHistosMapTH1F [private] |
Definition at line 37 of file DTVDriftSegmentCalibration.h.
typedef std::map<DTChamberId, std::vector<TH2F*> > DTVDriftSegmentCalibration::ChamberHistosMapTH2F [private] |
Definition at line 38 of file DTVDriftSegmentCalibration.h.
Definition at line 36 of file DTVDriftSegmentCalibration.cc.
References edm::ParameterSet::getUntrackedParameter(), rootFile_, and dtT0WireCalibration_cfg::rootFileName.
: select_(pset), theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")), //writeVDriftDB_(pset.getUntrackedParameter<bool>("writeVDriftDB", false)), theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")) { LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Constructor called!"; // the root file which will contain the histos string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTVDriftHistos.root"); rootFile_ = new TFile(rootFileName.c_str(), "RECREATE"); rootFile_->cd(); }
DTVDriftSegmentCalibration::~DTVDriftSegmentCalibration | ( | ) | [virtual] |
Definition at line 56 of file DTVDriftSegmentCalibration.cc.
References rootFile_.
{ rootFile_->Close(); LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Destructor called!"; }
void DTVDriftSegmentCalibration::analyze | ( | const edm::Event & | event, |
const edm::EventSetup & | eventSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 61 of file DTVDriftSegmentCalibration.cc.
References bookHistos(), DTChamberId, HcalObjRepresent::Fill(), edm::EventSetup::get(), LogTrace, rootFile_, select_, theCalibChamber_, theRecHits4DLabel_, theVDriftHistoMapTH1F_, theVDriftHistoMapTH2F_, and GeomDet::toGlobal().
{ rootFile_->cd(); // Get the DT Geometry ESHandle<DTGeometry> dtGeom; eventSetup.get<MuonGeometryRecord>().get(dtGeom); // Get the rechit collection from the event Handle<DTRecSegment4DCollection> all4DSegments; event.getByLabel(theRecHits4DLabel_, all4DSegments); DTChamberId chosenChamberId; if(theCalibChamber_ != "All") { stringstream linestr; int selWheel, selStation, selSector; linestr << theCalibChamber_; linestr >> selWheel >> selStation >> selSector; chosenChamberId = DTChamberId(selWheel, selStation, selSector); LogVerbatim("Calibration") << " Chosen chamber: " << chosenChamberId << endl; } // Loop over segments by chamber DTRecSegment4DCollection::id_iterator chamberIdIt; for(chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt){ // Calibrate just the chosen chamber/s if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue; // Book histos if(theVDriftHistoMapTH1F_.find(*chamberIdIt) == theVDriftHistoMapTH1F_.end()){ LogTrace("Calibration") << " Booking histos for Chamber: " << *chamberIdIt; bookHistos(*chamberIdIt); } // Get the chamber from the setup const DTChamber* chamber = dtGeom->chamber(*chamberIdIt); // Get the range for the corresponding ChamberId DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt)); // Loop over the rechits of this DetUnit for(DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment){ LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition() << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition()); if( !select_(*segment, event, eventSetup) ) continue; // Fill v-drift values if( (*segment).hasPhi() ) { //if( segment->phiSegment()->ist0Valid() ){ double segmentVDrift = segment->phiSegment()->vDrift(); if( segmentVDrift != 0.00 ){ (theVDriftHistoMapTH1F_[*chamberIdIt])[0]->Fill(segmentVDrift); (theVDriftHistoMapTH2F_[*chamberIdIt])[0]->Fill(segment->localPosition().x(),segmentVDrift); (theVDriftHistoMapTH2F_[*chamberIdIt])[1]->Fill(segment->localPosition().y(),segmentVDrift); } } // Probably not meaningful if( (*segment).hasZed() ){ //if( segment->zSegment()->ist0Valid() ){ double segmentVDrift = segment->zSegment()->vDrift(); if( segmentVDrift != 0.00 ){ (theVDriftHistoMapTH1F_[*chamberIdIt])[1]->Fill(segmentVDrift); } } } // DTRecSegment4DCollection::const_iterator segment } // DTRecSegment4DCollection::id_iterator chamberIdIt } // DTVDriftSegmentCalibration::analyze
void DTVDriftSegmentCalibration::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 50 of file DTVDriftSegmentCalibration.cc.
{
TH1::SetDefaultSumw2(true);
}
void DTVDriftSegmentCalibration::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | setup | ||
) | [virtual] |
void DTVDriftSegmentCalibration::bookHistos | ( | DTChamberId | chId | ) | [private] |
Definition at line 151 of file DTVDriftSegmentCalibration.cc.
References DTChamberId::sector(), relativeConstraints::station, DTChamberId::station(), theVDriftHistoMapTH1F_, theVDriftHistoMapTH2F_, and DTChamberId::wheel().
Referenced by analyze().
{ // Compose the chamber name stringstream wheel; wheel << chId.wheel(); stringstream station; station << chId.station(); stringstream sector; sector << chId.sector(); string chHistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str(); vector<TH1F*> histosTH1F; histosTH1F.push_back(new TH1F(("hRPhiVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Phi segments", 200, -0.4, 0.4)); if(chId.station() != 4) histosTH1F.push_back(new TH1F(("hRZVDriftCorr" + chHistoName).c_str(), "v-drift corr. from Z segments", 200, -0.4, 0.4)); vector<TH2F*> histosTH2F; histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosX" + chHistoName).c_str(), "v-drift corr. vs. segment x position", 250, -125., 125., 200, -0.4, 0.4)); histosTH2F.push_back(new TH2F(("hRPhiVDriftCorrVsSegmPosY" + chHistoName).c_str(), "v-drift corr. vs. segment y position", 250, -125., 125., 200, -0.4, 0.4)); //if(chId.station() != 4) ... theVDriftHistoMapTH1F_[chId] = histosTH1F; theVDriftHistoMapTH2F_[chId] = histosTH2F; }
void DTVDriftSegmentCalibration::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 130 of file DTVDriftSegmentCalibration.cc.
References rootFile_, theVDriftHistoMapTH1F_, and theVDriftHistoMapTH2F_.
{ rootFile_->cd(); LogVerbatim("Calibration") << "[DTVDriftSegmentCalibration] Writing histos to file!" << endl; for(ChamberHistosMapTH1F::const_iterator itChHistos = theVDriftHistoMapTH1F_.begin(); itChHistos != theVDriftHistoMapTH1F_.end(); ++itChHistos){ vector<TH1F*>::const_iterator itHistTH1F = (*itChHistos).second.begin(); vector<TH1F*>::const_iterator itHistTH1F_end = (*itChHistos).second.end(); for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write(); vector<TH2F*>::const_iterator itHistTH2F = theVDriftHistoMapTH2F_[(*itChHistos).first].begin(); vector<TH2F*>::const_iterator itHistTH2F_end = theVDriftHistoMapTH2F_[(*itChHistos).first].end(); for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write(); } /*if(writeVDriftDB_){ // ... }*/ }
TFile* DTVDriftSegmentCalibration::rootFile_ [private] |
Definition at line 47 of file DTVDriftSegmentCalibration.h.
Referenced by analyze(), DTVDriftSegmentCalibration(), endJob(), and ~DTVDriftSegmentCalibration().
Definition at line 41 of file DTVDriftSegmentCalibration.h.
Referenced by analyze().
std::string DTVDriftSegmentCalibration::theCalibChamber_ [private] |
Definition at line 45 of file DTVDriftSegmentCalibration.h.
Referenced by analyze().
Definition at line 43 of file DTVDriftSegmentCalibration.h.
Referenced by analyze().
Definition at line 48 of file DTVDriftSegmentCalibration.h.
Referenced by analyze(), bookHistos(), and endJob().
Definition at line 49 of file DTVDriftSegmentCalibration.h.
Referenced by analyze(), bookHistos(), and endJob().