CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

DTResidualCalibration Class Reference

#include <DTResidualCalibration.h>

Inheritance diagram for DTResidualCalibration:
edm::EDAnalyzer

List of all members.

Classes

class  DTResidualCalibration

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup)
void beginJob ()
void beginRun (const edm::Run &, const edm::EventSetup &)
 DTResidualCalibration (const edm::ParameterSet &pset)
 Constructor.
void endJob ()
virtual ~DTResidualCalibration ()
 Destructor.

Private Member Functions

void bookHistos (DTSuperLayerId slId)
void fillHistos (DTSuperLayerId slId, float distance, float residualOnDistance)
float segmentToWireDistance (const DTRecHit1D &recHit1D, const DTRecSegment4D &segment)

Private Attributes

const DTGeometrydtGeom_
std::map< DTSuperLayerId,
std::vector< TH1F * > > 
histoMapTH1F_
std::map< DTSuperLayerId,
std::vector< TH2F * > > 
histoMapTH2F_
std::string rootBaseDir_
TFile * rootFile_
edm::InputTag segment4DLabel_
DTSegmentSelector select_

Detailed Description

Extracts DT segment residuals

Date:
2011/02/22 18:43:20
Revision:
1.1

Definition at line 29 of file DTResidualCalibration.h.


Constructor & Destructor Documentation

Constructor.

Definition at line 36 of file DTResidualCalibration.cc.

References edm::ParameterSet::getUntrackedParameter(), rootFile_, and dtT0WireCalibration_cfg::rootFileName.

                                                                       :
  select_(pset),
  segment4DLabel_(pset.getParameter<edm::InputTag>("segment4DLabel")),
  rootBaseDir_(pset.getUntrackedParameter<std::string>("rootBaseDir","DT/Residuals")) {

  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Constructor called.";

  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName","residuals.root");
  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
  rootFile_->cd();
}
DTResidualCalibration::~DTResidualCalibration ( ) [virtual]

Destructor.

Definition at line 48 of file DTResidualCalibration.cc.

                                              {
  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Destructor called.";
}

Member Function Documentation

void DTResidualCalibration::analyze ( const edm::Event event,
const edm::EventSetup setup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 80 of file DTResidualCalibration.cc.

References DTGeometry::chamber(), DTRecHitSegmentResidual::compute(), filterCSVwithJSON::copy, dtGeom_, fillHistos(), LogTrace, rootFile_, segment4DLabel_, segmentToWireDistance(), select_, DTRecSegment2D::specificRecHits(), DTLayerId::superlayerId(), and GeomDet::toGlobal().

                                                                                     {
  rootFile_->cd();

  // Get the 4D rechits from the event
  edm::Handle<DTRecSegment4DCollection> segment4Ds;
  event.getByLabel(segment4DLabel_, segment4Ds);
 
  // Loop over segments by chamber
  DTRecSegment4DCollection::id_iterator chamberIdIt;
  for(chamberIdIt = segment4Ds->id_begin(); chamberIdIt != segment4Ds->id_end(); ++chamberIdIt){

     const DTChamber* chamber = dtGeom_->chamber(*chamberIdIt);

     // Get the range for the corresponding ChamberId
     DTRecSegment4DCollection::range range = segment4Ds->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, setup) ) continue;

        // Get all 1D RecHits at step 3 within the 4D segment
        std::vector<DTRecHit1D> recHits1D_S3;
  
        if( (*segment).hasPhi() ){
           const DTChamberRecSegment2D* phiSeg = (*segment).phiSegment();
           const std::vector<DTRecHit1D>& phiRecHits = phiSeg->specificRecHits();
           std::copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
        }

        if( (*segment).hasZed() ){
           const DTSLRecSegment2D* zSeg = (*segment).zSegment();
           const std::vector<DTRecHit1D>& zRecHits = zSeg->specificRecHits();
           std::copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
        }

        // Loop over 1D RecHit inside 4D segment
        for(std::vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin();
                                                    recHit1D != recHits1D_S3.end(); ++recHit1D) {
           const DTWireId wireId = recHit1D->wireId();

           float segmDistance = segmentToWireDistance(*recHit1D,*segment);
           if(segmDistance > 2.1) LogTrace("Calibration") << "WARNING: segment-wire distance: " << segmDistance;
           else                   LogTrace("Calibration") << "segment-wire distance: " << segmDistance;

           float residualOnDistance = DTRecHitSegmentResidual().compute(dtGeom_,*recHit1D,*segment);
           LogTrace("Calibration") << "Wire Id " << wireId << " residual on distance: " << residualOnDistance;

           fillHistos(wireId.superlayerId(), segmDistance, residualOnDistance);
        }
     }
  }

}
void DTResidualCalibration::beginJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 52 of file DTResidualCalibration.cc.

                                     {
  TH1::SetDefaultSumw2(true);
}
void DTResidualCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 56 of file DTResidualCalibration.cc.

References bookHistos(), DTGeometry::chambers(), dtGeom_, edm::EventSetup::get(), histoMapTH1F_, and edm::ESHandle< T >::product().

                                                                                  {
  
  // get the geometry
  edm::ESHandle<DTGeometry> dtGeomH; 
  setup.get<MuonGeometryRecord>().get(dtGeomH);
  dtGeom_ = dtGeomH.product();

  // Loop over all the chambers
  if(histoMapTH1F_.size() == 0){         
     std::vector<DTChamber*>::const_iterator ch_it = dtGeom_->chambers().begin();        
     std::vector<DTChamber*>::const_iterator ch_end = dtGeom_->chambers().end();         
     for (; ch_it != ch_end; ++ch_it) {          
        DTChamberId chId = (*ch_it)->id();
        std::vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();        
        std::vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();         
        // Loop over the SLs     
        for(; sl_it != sl_end; ++sl_it) { 
           DTSuperLayerId slId = (*sl_it)->id();
           bookHistos(slId);
        }
     }
  }
}
void DTResidualCalibration::bookHistos ( DTSuperLayerId  slId) [private]

Definition at line 186 of file DTResidualCalibration.cc.

References cmsMakeMELists::baseDir, histoMapTH1F_, histoMapTH2F_, rootBaseDir_, rootFile_, DTChamberId::sector(), mergeVDriftHistosByStation::sectorStr, DTChamberId::station(), mergeVDriftHistosByStation::stationStr, launcher::step, DTSuperLayerId::superlayer(), DTChamberId::wheel(), and mergeVDriftHistosByStation::wheelStr.

Referenced by beginRun().

                                                          {
  TH1AddDirectorySentry addDir;
  rootFile_->cd();

  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Booking histos for SL: " << slId;

  // Compose the chamber name
  std::stringstream wheelStr; wheelStr << slId.wheel(); 
  std::stringstream stationStr; stationStr << slId.station();   
  std::stringstream sectorStr; sectorStr << slId.sector();      
  std::stringstream superLayerStr; superLayerStr << slId.superlayer();
  // Define the step
  int step = 3;
  std::stringstream stepStr; stepStr << step;

  std::string slHistoName =
    "_STEP" + stepStr.str() +
    "_W" + wheelStr.str() +
    "_St" + stationStr.str() +
    "_Sec" + sectorStr.str() +
    "_SL" + superLayerStr.str();
  
  edm::LogVerbatim("Calibration") << "Accessing " << rootBaseDir_;
  TDirectory* baseDir = rootFile_->GetDirectory(rootBaseDir_.c_str());
  if(!baseDir) baseDir = rootFile_->mkdir(rootBaseDir_.c_str());
  edm::LogVerbatim("Calibration") << "Accessing " << ("Wheel" + wheelStr.str());
  TDirectory* wheelDir = baseDir->GetDirectory(("Wheel" + wheelStr.str()).c_str());
  if(!wheelDir) wheelDir = baseDir->mkdir(("Wheel" + wheelStr.str()).c_str());
  edm::LogVerbatim("Calibration") << "Accessing " << ("Station" + stationStr.str());
  TDirectory* stationDir = wheelDir->GetDirectory(("Station" + stationStr.str()).c_str());
  if(!stationDir) stationDir = wheelDir->mkdir(("Station" + stationStr.str()).c_str());
  edm::LogVerbatim("Calibration") << "Accessing " << ("Sector" + sectorStr.str());
  TDirectory* sectorDir = stationDir->GetDirectory(("Sector" + sectorStr.str()).c_str());
  if(!sectorDir) sectorDir = stationDir->mkdir(("Sector" + sectorStr.str()).c_str()); 

  /*std::string dirName = rootBaseDir_ + "/Wheel" + wheelStr.str() +
                                       "/Station" + stationStr.str() +
                                       "/Sector" + sectorStr.str();

  TDirectory* dir = rootFile_->GetDirectory(dirName.c_str());
  if(!dir) dir = rootFile_->mkdir(dirName.c_str());
  dir->cd();*/
  sectorDir->cd();
  // Create the monitor elements
  std::vector<TH1F*> histosTH1F;
  histosTH1F.push_back(new TH1F(("hResDist"+slHistoName).c_str(),
                                 "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
                                 200, -0.4, 0.4));
  std::vector<TH2F*> histosTH2F;
  histosTH2F.push_back(new TH2F(("hResDistVsDist"+slHistoName).c_str(),
                                 "Residuals on the dist. (cm) from wire (rec_hit - segm_extr) vs dist. (cm)",
                                 100, 0, 2.5, 200, -0.4, 0.4));
  histoMapTH1F_[slId] = histosTH1F;
  histoMapTH2F_[slId] = histosTH2F;
}
void DTResidualCalibration::endJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 165 of file DTResidualCalibration.cc.

References rootFile_.

                                  {
  
  edm::LogVerbatim("Calibration") << "[DTResidualCalibration] Writing histos to file.";
  rootFile_->cd();
  rootFile_->Write();
  rootFile_->Close();

  /*std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos = histoMapTH1F_.begin();
  std::map<DTSuperLayerId, std::vector<TH1F*> >::const_iterator itSlHistos_end = histoMapTH1F_.end(); 
  for(; itSlHistos != itSlHistos_end; ++itSlHistos){
     std::vector<TH1F*>::const_iterator itHistTH1F = (*itSlHistos).second.begin();
     std::vector<TH1F*>::const_iterator itHistTH1F_end = (*itSlHistos).second.end();
     for(; itHistTH1F != itHistTH1F_end; ++itHistTH1F) (*itHistTH1F)->Write();

     std::vector<TH2F*>::const_iterator itHistTH2F = histoMapTH2F_[(*itSlHistos).first].begin();
     std::vector<TH2F*>::const_iterator itHistTH2F_end = histoMapTH2F_[(*itSlHistos).first].end();
     for(; itHistTH2F != itHistTH2F_end; ++itHistTH2F) (*itHistTH2F)->Write();
  }*/

}
void DTResidualCalibration::fillHistos ( DTSuperLayerId  slId,
float  distance,
float  residualOnDistance 
) [private]

Definition at line 243 of file DTResidualCalibration.cc.

References histoMapTH1F_, and histoMapTH2F_.

Referenced by analyze().

                                                                 {
  std::vector<TH1F*> const& histosTH1F = histoMapTH1F_[slId];
  std::vector<TH2F*> const& histosTH2F = histoMapTH2F_[slId];                          
  histosTH1F[0]->Fill(residualOnDistance);
  histosTH2F[0]->Fill(distance, residualOnDistance);
}
float DTResidualCalibration::segmentToWireDistance ( const DTRecHit1D recHit1D,
const DTRecSegment4D segment 
) [private]

Definition at line 138 of file DTResidualCalibration.cc.

References DTGeometry::chamber(), DTSuperLayerId::chamberId(), funct::cos(), dtGeom_, DTGeometry::layer(), DTWireId::layerId(), DTRecSegment4D::localDirection(), DTRecHit1D::localPosition(), DTRecSegment4D::localPosition(), DTLayer::specificTopology(), DTSuperLayerId::superlayer(), PV3DBase< T, PVType, FrameType >::theta(), GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTRecHit1D::wireId(), DTTopology::wirePosition(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

                                                                                                           {

  // Get the layer and the wire position
  const DTWireId wireId = recHit1D.wireId();
  const DTLayer* layer = dtGeom_->layer(wireId);
  float wireX = layer->specificTopology().wirePosition(wireId.wire());
      
  // Extrapolate the segment to the z of the wire
  // Get wire position in chamber RF
  // (y and z must be those of the hit to be coherent in the transf. of RF in case of rotations of the layer alignment)
  LocalPoint wirePosInLay(wireX,recHit1D.localPosition().y(),recHit1D.localPosition().z());
  GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
  const DTChamber* chamber = dtGeom_->chamber(wireId.layerId().chamberId());
  LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
      
  // Segment position at Wire z in chamber local frame
  LocalPoint segPosAtZWire = segment.localPosition()    + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
      
  // Compute the distance of the segment from the wire
  int sl = wireId.superlayer();
  float segmDistance = -1;
  if(sl == 1 || sl == 3) segmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
  else if(sl == 2)       segmDistance =  fabs(segPosAtZWire.y() - wirePosInChamber.y());
     
  return segmDistance;
}

Member Data Documentation

Definition at line 56 of file DTResidualCalibration.h.

Referenced by analyze(), beginRun(), and segmentToWireDistance().

std::map<DTSuperLayerId, std::vector<TH1F*> > DTResidualCalibration::histoMapTH1F_ [private]

Definition at line 58 of file DTResidualCalibration.h.

Referenced by beginRun(), bookHistos(), and fillHistos().

std::map<DTSuperLayerId, std::vector<TH2F*> > DTResidualCalibration::histoMapTH2F_ [private]

Definition at line 59 of file DTResidualCalibration.h.

Referenced by bookHistos(), and fillHistos().

std::string DTResidualCalibration::rootBaseDir_ [private]

Definition at line 52 of file DTResidualCalibration.h.

Referenced by bookHistos().

Definition at line 54 of file DTResidualCalibration.h.

Referenced by analyze(), bookHistos(), DTResidualCalibration(), and endJob().

Definition at line 51 of file DTResidualCalibration.h.

Referenced by analyze().

Definition at line 50 of file DTResidualCalibration.h.

Referenced by analyze().