CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTracker/DeDx/plugins/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.2 2012/01/17 13:46:51 innocent Exp $
00019 //
00020 //
00021 
00022 
00023 // system include files
00024 #include <memory>
00025 #include "DataFormats/Common/interface/ValueMap.h"
00026 
00027 #include "RecoTracker/DeDx/plugins/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    shapetest = iConfig.getParameter<bool>("ShapeTest");
00075    useCalibration = iConfig.getParameter<bool>("UseCalibration");
00076    m_calibrationPath = iConfig.getParameter<string>("calibrationPath");
00077 
00078    if(!usePixel && !useStrip)
00079    edm::LogWarning("DeDxHitsProducer") << "Pixel Hits AND Strip Hits will not be used to estimate dEdx --> BUG, Please Update the config file";
00080 }
00081 
00082 
00083 DeDxEstimatorProducer::~DeDxEstimatorProducer()
00084 {
00085   delete m_estimator;
00086 }
00087 
00088 // ------------ method called once each job just before starting event loop  ------------
00089 void  DeDxEstimatorProducer::beginRun(edm::Run & run, const edm::EventSetup& iSetup)
00090 {
00091    if(MODsColl.size()!=0)return;
00092 
00093 
00094    edm::ESHandle<TrackerGeometry> tkGeom;
00095    iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
00096 
00097    vector<GeomDet*> Det = tkGeom->dets();
00098    for(unsigned int i=0;i<Det.size();i++){
00099       DetId  Detid  = Det[i]->geographicalId();
00100 
00101        StripGeomDetUnit* StripDetUnit = dynamic_cast<StripGeomDetUnit*> (Det[i]);
00102        PixelGeomDetUnit* PixelDetUnit = dynamic_cast<PixelGeomDetUnit*> (Det[i]);
00103 
00104        double Thick=-1, Dist=-1, Norma=-1;
00105        if(StripDetUnit){
00106           Dist      = StripDetUnit->position().mag();
00107           Thick     = StripDetUnit->surface().bounds().thickness();
00108           Norma     = MeVperADCStrip/Thick;
00109        }else if(PixelDetUnit){
00110           Dist      = PixelDetUnit->position().mag();
00111           Thick     = PixelDetUnit->surface().bounds().thickness();
00112           Norma     = MeVperADCPixel/Thick;
00113        }
00114 
00115        stModInfo* MOD       = new stModInfo;
00116        MOD->DetId           = Detid.rawId();
00117        MOD->Thickness       = Thick;
00118        MOD->Distance        = Dist;
00119        MOD->Normalization   = Norma;
00120        MOD->Gain            = 1;
00121        MODsColl[MOD->DetId] = MOD;
00122    }
00123 
00124    MakeCalibrationMap();
00125 }
00126 
00127 // ------------ method called once each job just after ending the event loop  ------------
00128 void  DeDxEstimatorProducer::endJob(){
00129    MODsColl.clear();
00130 }
00131 
00132 
00133 
00134 
00135 void DeDxEstimatorProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00136 {
00137   auto_ptr<ValueMap<DeDxData> > trackDeDxEstimateAssociation(new ValueMap<DeDxData> );  
00138   ValueMap<DeDxData>::Filler filler(*trackDeDxEstimateAssociation);
00139 
00140   Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
00141   iEvent.getByLabel(m_trajTrackAssociationTag, trajTrackAssociationHandle);
00142   const TrajTrackAssociationCollection & TrajToTrackMap = *trajTrackAssociationHandle.product();
00143 
00144   edm::Handle<reco::TrackCollection> trackCollectionHandle;
00145   iEvent.getByLabel(m_tracksTag,trackCollectionHandle);
00146 
00147   size_t n =  TrajToTrackMap.size();
00148   std::vector<DeDxData> dedxEstimate(n);
00149 
00150   //assume trajectory collection size is equal to track collection size and that order is kept
00151   int j=0;
00152   for(TrajTrackAssociationCollection::const_iterator cit=TrajToTrackMap.begin(); cit!=TrajToTrackMap.end(); cit++,j++){
00153      
00154      const edm::Ref<std::vector<Trajectory> > traj = cit->key;
00155      const reco::TrackRef track = cit->val;
00156 
00157      DeDxHitCollection dedxHits;
00158      vector<DeDxTools::RawHits> hits; 
00159 //     DeDxTools::trajectoryRawHits(traj, hits, usePixel, useStrip);
00160 
00161     const vector<TrajectoryMeasurement> & measurements = traj->measurements();
00162     for(vector<TrajectoryMeasurement>::const_iterator it = measurements.begin(); it!=measurements.end(); it++){
00163       TrajectoryStateOnSurface trajState=it->updatedState();
00164       if( !trajState.isValid()) continue;
00165      
00166       const TrackingRecHit * recHit=(*it->recHit()).hit();
00167       LocalVector trackDirection = trajState.localDirection();
00168       double cosine = trackDirection.z()/trackDirection.mag();
00169               
00170        if(const SiStripMatchedRecHit2D* matchedHit=dynamic_cast<const SiStripMatchedRecHit2D*>(recHit)){
00171            if(!useStrip) continue;
00172            DeDxTools::RawHits mono,stereo; 
00173            mono.trajectoryMeasurement = &(*it);
00174            stereo.trajectoryMeasurement = &(*it);
00175            mono.angleCosine = cosine; 
00176            stereo.angleCosine = cosine;
00177 
00178            mono.charge = getCharge(DeDxTools::GetCluster(matchedHit->monoHit()),mono.NSaturating,matchedHit->monoId());
00179            stereo.charge = getCharge(DeDxTools::GetCluster(matchedHit->stereoHit()),stereo.NSaturating,matchedHit->stereoId());
00180 
00181            mono.detId= matchedHit->monoId();
00182            stereo.detId= matchedHit->stereoId();
00183 
00184            if(shapetest && !(DeDxTools::shapeSelection((DeDxTools::GetCluster(matchedHit->stereoHit()))->amplitudes()))) hits.push_back(stereo);
00185            if(shapetest && !(DeDxTools::shapeSelection((DeDxTools::GetCluster(matchedHit->  monoHit()))->amplitudes()))) hits.push_back(mono);
00186         }else if(const ProjectedSiStripRecHit2D* projectedHit=dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit)) {
00187            if(!useStrip) continue;
00188            const SiStripRecHit2D* singleHit=&(projectedHit->originalHit());
00189            DeDxTools::RawHits mono;
00190 
00191            mono.trajectoryMeasurement = &(*it);
00192            mono.angleCosine = cosine;
00193            mono.charge = getCharge(DeDxTools::GetCluster(singleHit),mono.NSaturating,singleHit->geographicalId());
00194            mono.detId= singleHit->geographicalId();
00195            if(shapetest && !(DeDxTools::shapeSelection((DeDxTools::GetCluster(singleHit))->amplitudes()))) continue;
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(DeDxTools::GetCluster(singleHit),mono.NSaturating,singleHit->geographicalId());
00204            mono.detId= singleHit->geographicalId();
00205            if(shapetest && !(DeDxTools::shapeSelection((DeDxTools::GetCluster(singleHit))->amplitudes()))) continue;
00206            hits.push_back(mono); 
00207         }else if(const SiStripRecHit1D* single1DHit=dynamic_cast<const SiStripRecHit1D*>(recHit)){
00208            if(!useStrip) continue;
00209            DeDxTools::RawHits mono;
00210                
00211            mono.trajectoryMeasurement = &(*it);
00212            mono.angleCosine = cosine; 
00213            mono.charge = getCharge(DeDxTools::GetCluster(single1DHit),mono.NSaturating,single1DHit->geographicalId());
00214            mono.detId= single1DHit->geographicalId();
00215            if(shapetest && !(DeDxTools::shapeSelection((DeDxTools::GetCluster(single1DHit))->amplitudes()))) continue;
00216            hits.push_back(mono); 
00217         }else if(const SiPixelRecHit* pixelHit=dynamic_cast<const SiPixelRecHit*>(recHit)){
00218            if(!usePixel) continue;
00219 
00220            DeDxTools::RawHits pixel;
00221 
00222            pixel.trajectoryMeasurement = &(*it);
00223            pixel.angleCosine = cosine; 
00224            pixel.charge = pixelHit->cluster()->charge();
00225            pixel.NSaturating=-1;
00226            pixel.detId= pixelHit->geographicalId();
00227            hits.push_back(pixel);
00228        } 
00229     }
00231 
00232      int NClusterSaturating = 0; 
00233      for(size_t i=0; i < hits.size(); i++)
00234      {
00235          stModInfo* MOD = MODsColl[hits[i].detId];
00236          float pathLen = MOD->Thickness/std::abs(hits[i].angleCosine);
00237          float charge  = MOD->Normalization*hits[i].charge*std::abs(hits[i].angleCosine);
00238          dedxHits.push_back( DeDxHit( charge, MOD->Distance, pathLen, hits[i].detId) );
00239 
00240          if(hits[i].NSaturating>0)NClusterSaturating++;
00241      }
00242 
00243      sort(dedxHits.begin(),dedxHits.end(),less<DeDxHit>());   
00244      std::pair<float,float> val_and_error = m_estimator->dedx(dedxHits);
00245 
00246      //WARNING: Since the dEdX Error is not properly computed for the moment
00247      //It was decided to store the number of saturating cluster in that dataformat
00248      val_and_error.second = NClusterSaturating; 
00249 
00250      dedxEstimate[j] = DeDxData(val_and_error.first, val_and_error.second, dedxHits.size() );
00251   }
00252   filler.insert(trackCollectionHandle, dedxEstimate.begin(), dedxEstimate.end());
00253 
00254   // really fill the association map
00255   filler.fill();
00256    // put into the event 
00257   iEvent.put(trackDeDxEstimateAssociation);   
00258 }
00259 
00260 
00261 
00262 void DeDxEstimatorProducer::MakeCalibrationMap()
00263 {
00264    if(!useCalibration)return;
00265 
00266    TChain* t1 = new TChain("SiStripCalib/APVGain");
00267    t1->Add(m_calibrationPath.c_str());
00268 
00269    unsigned int  tree_DetId;
00270    unsigned char tree_APVId;
00271    double        tree_Gain;
00272 
00273    t1->SetBranchAddress("DetId"             ,&tree_DetId      );
00274    t1->SetBranchAddress("APVId"             ,&tree_APVId      );
00275    t1->SetBranchAddress("Gain"              ,&tree_Gain       );
00276 
00277 
00278 
00279    for (unsigned int ientry = 0; ientry < t1->GetEntries(); ientry++) {
00280        t1->GetEntry(ientry);
00281        stModInfo* MOD  = MODsColl[tree_DetId];
00282        MOD->Gain = tree_Gain;
00283    }
00284 
00285    delete t1;
00286 
00287 }
00288 
00289 
00290 int DeDxEstimatorProducer::getCharge(const SiStripCluster*   Cluster, int& Saturating_Strips,
00291                                      const uint32_t & DetId)
00292 {
00293    //const vector<uint8_t>&  Ampls       = Cluster->amplitudes();
00294    const vector<uint8_t>&  Ampls       = Cluster->amplitudes();
00295    //uint32_t                DetId       = Cluster->geographicalId();
00296 
00297 //   float G=1.0f;
00298 
00299    int toReturn = 0;
00300    Saturating_Strips = 0;
00301    for(unsigned int i=0;i<Ampls.size();i++){
00302       int CalibratedCharge = Ampls[i];
00303 
00304       if(useCalibration){
00305          stModInfo* MOD = MODsColl[DetId];
00306 //         G = MOD->Gain;
00307          CalibratedCharge = (int)(CalibratedCharge / MOD->Gain );
00308          if(CalibratedCharge>=1024){
00309             CalibratedCharge = 255;
00310          }else if(CalibratedCharge>=255){
00311             CalibratedCharge = 254;
00312          } 
00313       }
00314 
00315       toReturn+=CalibratedCharge;
00316       if(CalibratedCharge>=254)Saturating_Strips++;
00317    }
00318 
00319 //   printf("Charge = %i --> %i  (Gain=%f)\n", accumulate(Ampls.begin(), Ampls.end(), 0), toReturn, G);         
00320 
00321 
00322    return toReturn;
00323 }
00324 
00325 
00326 
00327 
00328 //define this as a plug-in
00329 DEFINE_FWK_MODULE(DeDxEstimatorProducer);