CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoTracker/DeDx/src/DeDxEstimatorProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    DeDxEstimatorProducer
00004 // Class:      DeDxEstimatorProducer
00005 // 
00013 //
00014 // Original Author:  andrea
00015 //         Created:  Thu May 31 14:09:02 CEST 2007
00016 //    Code Updates:  loic Quertenmont (querten)
00017 //         Created:  Thu May 10 14:09:02 CEST 2008
00018 // $Id: DeDxEstimatorProducer.cc,v 1.28 2010/05/31 13:07:25 querten Exp $
00019 //
00020 //
00021 
00022 
00023 // system include files
00024 #include <memory>
00025 #include "DataFormats/Common/interface/ValueMap.h"
00026 
00027 #include "RecoTracker/DeDx/interface/DeDxEstimatorProducer.h"
00028 #include "DataFormats/TrackReco/interface/DeDxData.h"
00029 #include "DataFormats/TrackReco/interface/TrackDeDxHits.h"
00030 #include "DataFormats/TrackReco/interface/DeDxHit.h"
00031 #include "DataFormats/TrackReco/interface/Track.h"
00032 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00033 
00034 #include "RecoTracker/DeDx/interface/GenericAverageDeDxEstimator.h"
00035 #include "RecoTracker/DeDx/interface/TruncatedAverageDeDxEstimator.h"
00036 #include "RecoTracker/DeDx/interface/MedianDeDxEstimator.h"
00037 #include "RecoTracker/DeDx/interface/UnbinnedFitDeDxEstimator.h"
00038 
00039 #include "FWCore/Framework/interface/ESHandle.h"
00040 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00041 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00042 
00043 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00044 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00045 
00046 
00047 using namespace reco;
00048 using namespace std;
00049 using namespace edm;
00050 
00051 DeDxEstimatorProducer::DeDxEstimatorProducer(const edm::ParameterSet& iConfig)
00052 {
00053 
00054    produces<ValueMap<DeDxData> >();
00055 
00056 
00057    string estimatorName = iConfig.getParameter<string>("estimator");
00058    if(estimatorName == "median")      m_estimator = new MedianDeDxEstimator(-2.);
00059    if(estimatorName == "generic")     m_estimator = new GenericAverageDeDxEstimator  (iConfig.getParameter<double>("exponent"));
00060    if(estimatorName == "truncated")   m_estimator = new TruncatedAverageDeDxEstimator(iConfig.getParameter<double>("fraction"));
00061    if(estimatorName == "unbinnedFit") m_estimator = new UnbinnedFitDeDxEstimator();
00062 
00063    MaxNrStrips         = iConfig.getUntrackedParameter<unsigned>("maxNrStrips"        ,  255);
00064    MinTrackHits        = iConfig.getUntrackedParameter<unsigned>("MinTrackHits"       ,  4);
00065 
00066    m_tracksTag = iConfig.getParameter<edm::InputTag>("tracks");
00067    m_trajTrackAssociationTag   = iConfig.getParameter<edm::InputTag>("trajectoryTrackAssociation");
00068 
00069    usePixel = iConfig.getParameter<bool>("UsePixel"); 
00070    useStrip = iConfig.getParameter<bool>("UseStrip");
00071    MeVperADCPixel = iConfig.getParameter<double>("MeVperADCPixel"); 
00072    MeVperADCStrip = iConfig.getParameter<double>("MeVperADCStrip"); 
00073 
00074    useCalibration = iConfig.getParameter<bool>("UseCalibration");
00075    m_calibrationPath = iConfig.getParameter<string>("calibrationPath");
00076 
00077    if(!usePixel && !useStrip)
00078    edm::LogWarning("DeDxHitsProducer") << "Pixel Hits AND Strip Hits will not be used to estimate dEdx --> BUG, Please Update the config file";
00079 }
00080 
00081 
00082 DeDxEstimatorProducer::~DeDxEstimatorProducer()
00083 {
00084   delete m_estimator;
00085 }
00086 
00087 
00088 
00089 // ------------ method called once each job just before starting event loop  ------------
00090 void  DeDxEstimatorProducer::beginRun(edm::Run & run, const edm::EventSetup& iSetup)
00091 {
00092    if(MODsColl.size()!=0)return;
00093 
00094 
00095    edm::ESHandle<TrackerGeometry> tkGeom;
00096    iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
00097 
00098    vector<GeomDet*> Det = tkGeom->dets();
00099    for(unsigned int i=0;i<Det.size();i++){
00100       DetId  Detid  = Det[i]->geographicalId();
00101 
00102        StripGeomDetUnit* StripDetUnit = dynamic_cast<StripGeomDetUnit*> (Det[i]);
00103        PixelGeomDetUnit* PixelDetUnit = dynamic_cast<PixelGeomDetUnit*> (Det[i]);
00104 
00105        double Thick=-1, Dist=-1, Norma=-1;
00106        if(StripDetUnit){
00107           Dist      = StripDetUnit->position().mag();
00108           Thick     = StripDetUnit->surface().bounds().thickness();
00109           Norma     = MeVperADCStrip/Thick;
00110        }else if(PixelDetUnit){
00111           Dist      = PixelDetUnit->position().mag();
00112           Thick     = PixelDetUnit->surface().bounds().thickness();
00113           Norma     = MeVperADCPixel/Thick;
00114        }
00115 
00116        stModInfo* MOD       = new stModInfo;
00117        MOD->DetId           = Detid.rawId();
00118        MOD->Thickness       = Thick;
00119        MOD->Distance        = Dist;
00120        MOD->Normalization   = Norma;
00121        MOD->Gain            = 1;
00122        MODsColl[MOD->DetId] = MOD;
00123    }
00124 
00125    MakeCalibrationMap();
00126 }
00127 
00128 // ------------ method called once each job just after ending the event loop  ------------
00129 void  DeDxEstimatorProducer::endJob(){
00130    MODsColl.clear();
00131 }
00132 
00133 
00134 
00135 
00136 void DeDxEstimatorProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00137 {
00138   auto_ptr<ValueMap<DeDxData> > trackDeDxEstimateAssociation(new ValueMap<DeDxData> );  
00139   ValueMap<DeDxData>::Filler filler(*trackDeDxEstimateAssociation);
00140 
00141   Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
00142   iEvent.getByLabel(m_trajTrackAssociationTag, trajTrackAssociationHandle);
00143   const TrajTrackAssociationCollection TrajToTrackMap = *trajTrackAssociationHandle.product();
00144 
00145   edm::Handle<reco::TrackCollection> trackCollectionHandle;
00146   iEvent.getByLabel(m_tracksTag,trackCollectionHandle);
00147 
00148   size_t n =  TrajToTrackMap.size();
00149   std::vector<DeDxData> dedxEstimate(n);
00150 
00151   //assume trajectory collection size is equal to track collection size and that order is kept
00152   int j=0;
00153   for(TrajTrackAssociationCollection::const_iterator cit=TrajToTrackMap.begin(); cit!=TrajToTrackMap.end(); cit++,j++){
00154      
00155      const edm::Ref<std::vector<Trajectory> > traj = cit->key;
00156      const reco::TrackRef track = cit->val;
00157 
00158      DeDxHitCollection dedxHits;
00159      vector<DeDxTools::RawHits> hits; 
00160 //     DeDxTools::trajectoryRawHits(traj, hits, usePixel, useStrip);
00161 
00162     const vector<TrajectoryMeasurement> & measurements = traj->measurements();
00163     for(vector<TrajectoryMeasurement>::const_iterator it = measurements.begin(); it!=measurements.end(); it++){
00164       TrajectoryStateOnSurface trajState=it->updatedState();
00165       if( !trajState.isValid()) continue;
00166      
00167       const TrackingRecHit * recHit=(*it->recHit()).hit();
00168       LocalVector trackDirection = trajState.localDirection();
00169       double cosine = trackDirection.z()/trackDirection.mag();
00170               
00171        if(const SiStripMatchedRecHit2D* matchedHit=dynamic_cast<const SiStripMatchedRecHit2D*>(recHit)){
00172            if(!useStrip) continue;
00173            DeDxTools::RawHits mono,stereo; 
00174            mono.trajectoryMeasurement = &(*it);
00175            stereo.trajectoryMeasurement = &(*it);
00176            mono.angleCosine = cosine; 
00177            stereo.angleCosine = cosine;
00178 
00179            mono.charge = getCharge((matchedHit->monoHit()->cluster()).get(),mono.NSaturating);
00180            stereo.charge = getCharge((matchedHit->stereoHit()->cluster()).get(),stereo.NSaturating);
00181 
00182            mono.detId= matchedHit->monoHit()->geographicalId();
00183            stereo.detId= matchedHit->stereoHit()->geographicalId();
00184 
00185            hits.push_back(mono);
00186            hits.push_back(stereo);
00187         }else if(const ProjectedSiStripRecHit2D* projectedHit=dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit)) {
00188            if(!useStrip) continue;
00189            const SiStripRecHit2D* singleHit=&(projectedHit->originalHit());
00190            DeDxTools::RawHits mono;
00191 
00192            mono.trajectoryMeasurement = &(*it);
00193            mono.angleCosine = cosine;
00194            mono.charge = getCharge((singleHit->cluster()).get(),mono.NSaturating);
00195            mono.detId= singleHit->geographicalId();
00196            hits.push_back(mono);
00197         }else if(const SiStripRecHit2D* singleHit=dynamic_cast<const SiStripRecHit2D*>(recHit)){
00198            if(!useStrip) continue;
00199            DeDxTools::RawHits mono;
00200                
00201            mono.trajectoryMeasurement = &(*it);
00202            mono.angleCosine = cosine; 
00203            mono.charge = getCharge((singleHit->cluster()).get(),mono.NSaturating);
00204            mono.detId= singleHit->geographicalId();
00205            hits.push_back(mono); 
00206         }else if(const SiStripRecHit1D* single1DHit=dynamic_cast<const SiStripRecHit1D*>(recHit)){
00207            if(!useStrip) continue;
00208            DeDxTools::RawHits mono;
00209                
00210            mono.trajectoryMeasurement = &(*it);
00211            mono.angleCosine = cosine; 
00212            mono.charge = getCharge((single1DHit->cluster()).get(),mono.NSaturating);
00213            mono.detId= single1DHit->geographicalId();
00214            hits.push_back(mono); 
00215         }else if(const SiPixelRecHit* pixelHit=dynamic_cast<const SiPixelRecHit*>(recHit)){
00216            if(!usePixel) continue;
00217 
00218            DeDxTools::RawHits pixel;
00219 
00220            pixel.trajectoryMeasurement = &(*it);
00221            pixel.angleCosine = cosine; 
00222            pixel.charge = pixelHit->cluster()->charge();
00223            pixel.NSaturating=-1;
00224            pixel.detId= pixelHit->geographicalId();
00225            hits.push_back(pixel);
00226        } 
00227     }
00229 
00230      int NClusterSaturating = 0; 
00231      for(size_t i=0; i < hits.size(); i++)
00232      {
00233          stModInfo* MOD = MODsColl[hits[i].detId];
00234          float pathLen = MOD->Thickness/std::abs(hits[i].angleCosine);
00235          float charge  = MOD->Normalization*hits[i].charge*std::abs(hits[i].angleCosine);
00236          dedxHits.push_back( DeDxHit( charge, MOD->Distance, pathLen, hits[i].detId) );
00237 
00238          if(hits[i].NSaturating>0)NClusterSaturating++;
00239      }
00240 
00241      sort(dedxHits.begin(),dedxHits.end(),less<DeDxHit>());   
00242      std::pair<float,float> val_and_error = m_estimator->dedx(dedxHits);
00243 
00244      //WARNING: Since the dEdX Error is not properly computed for the moment
00245      //It was decided to store the number of saturating cluster in that dataformat
00246      val_and_error.second = NClusterSaturating; 
00247 
00248      dedxEstimate[j] = DeDxData(val_and_error.first, val_and_error.second, dedxHits.size() );
00249   }
00250   filler.insert(trackCollectionHandle, dedxEstimate.begin(), dedxEstimate.end());
00251 
00252   // really fill the association map
00253   filler.fill();
00254    // put into the event 
00255   iEvent.put(trackDeDxEstimateAssociation);   
00256 }
00257 
00258 
00259 
00260 void DeDxEstimatorProducer::MakeCalibrationMap()
00261 {
00262    if(!useCalibration)return;
00263 
00264    TChain* t1 = new TChain("SiStripCalib/APVGain");
00265    t1->Add(m_calibrationPath.c_str());
00266 
00267    unsigned int  tree_DetId;
00268    unsigned char tree_APVId;
00269    double        tree_Gain;
00270 
00271    t1->SetBranchAddress("DetId"             ,&tree_DetId      );
00272    t1->SetBranchAddress("APVId"             ,&tree_APVId      );
00273    t1->SetBranchAddress("Gain"              ,&tree_Gain       );
00274 
00275 
00276 
00277    for (unsigned int ientry = 0; ientry < t1->GetEntries(); ientry++) {
00278        t1->GetEntry(ientry);
00279        stModInfo* MOD  = MODsColl[tree_DetId];
00280        MOD->Gain = tree_Gain;
00281    }
00282 
00283    delete t1;
00284 
00285 }
00286 
00287 
00288 int DeDxEstimatorProducer::getCharge(const SiStripCluster*   Cluster, int& Saturating_Strips)
00289 {
00290    const vector<uint8_t>&  Ampls       = Cluster->amplitudes();
00291    uint32_t                DetId       = Cluster->geographicalId();
00292 
00293 //   float G=1.0f;
00294 
00295    int toReturn = 0;
00296    Saturating_Strips = 0;
00297    for(unsigned int i=0;i<Ampls.size();i++){
00298       int CalibratedCharge = Ampls[i];
00299 
00300       if(useCalibration){
00301          stModInfo* MOD = MODsColl[DetId];
00302 //         G = MOD->Gain;
00303          CalibratedCharge = (int)(CalibratedCharge / MOD->Gain );
00304          if(CalibratedCharge>=1024){
00305             CalibratedCharge = 255;
00306          }else if(CalibratedCharge>=255){
00307             CalibratedCharge = 254;
00308          } 
00309       }
00310 
00311       toReturn+=CalibratedCharge;
00312       if(CalibratedCharge>=254)Saturating_Strips++;
00313    }
00314 
00315 //   printf("Charge = %i --> %i  (Gain=%f)\n", accumulate(Ampls.begin(), Ampls.end(), 0), toReturn, G);         
00316 
00317 
00318    return toReturn;
00319 }
00320 
00321 
00322 
00323 
00324 //define this as a plug-in
00325 DEFINE_FWK_MODULE(DeDxEstimatorProducer);