CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTracker/FinalTrackSelectors/src/MultiTrackSelector.cc

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/src/MultiTrackSelector.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/Common/interface/ValueMap.h"
00005 #include "CondFormats/DataRecord/interface/GBRWrapperRcd.h"
00006 #include "FWCore/Framework/interface/EventSetup.h"
00007 #include "FWCore/Framework/interface/ESHandle.h" 
00008 
00009 //#include "RecoTracker/DebugTools/interface/FixTrackHitPattern.h"
00010 #include <Math/DistFunc.h>
00011 #include "TMath.h"
00012 #include "TFile.h"
00013 
00014 using reco::modules::MultiTrackSelector;
00015 
00016 MultiTrackSelector::MultiTrackSelector()
00017 {
00018   useForestFromDB_ = true;
00019   forest_ = 0;
00020   gbrVals_ = 0;
00021 }
00022 
00023 MultiTrackSelector::MultiTrackSelector( const edm::ParameterSet & cfg ) :
00024   src_( cfg.getParameter<edm::InputTag>( "src" ) ),
00025   beamspot_( cfg.getParameter<edm::InputTag>( "beamspot" ) ),
00026   useVertices_( cfg.getParameter<bool>( "useVertices" ) ),
00027   useVtxError_( cfg.getParameter<bool>( "useVtxError" ) ),
00028   vertices_( useVertices_ ? cfg.getParameter<edm::InputTag>( "vertices" ) : edm::InputTag("NONE"))
00029   // now get the pset for each selector
00030 {
00031 
00032   useAnyMVA_ = false;
00033   forestLabel_ = "MVASelectorIter0";
00034   std::string type = "BDTG";
00035   useForestFromDB_ = true;
00036   dbFileName_ = "";
00037 
00038   forest_ = 0;
00039   gbrVals_ = 0;
00040 
00041   if(cfg.exists("useAnyMVA")) useAnyMVA_ = cfg.getParameter<bool>("useAnyMVA");
00042 
00043   if(useAnyMVA_){
00044     gbrVals_ = new float[11];
00045     if(cfg.exists("mvaType"))type = cfg.getParameter<std::string>("mvaType");
00046     if(cfg.exists("GBRForestLabel"))forestLabel_ = cfg.getParameter<std::string>("GBRForestLabel");
00047     if(cfg.exists("GBRForestFileName")){
00048       dbFileName_ = cfg.getParameter<std::string>("GBRForestFileName");
00049       useForestFromDB_ = false;
00050     }
00051 
00052      mvaType_ = type;
00053   }
00054   std::vector<edm::ParameterSet> trkSelectors( cfg.getParameter<std::vector< edm::ParameterSet> >("trackSelectors") );
00055   qualityToSet_.reserve(trkSelectors.size());
00056   vtxNumber_.reserve(trkSelectors.size());
00057   vertexCut_.reserve(trkSelectors.size());
00058   res_par_.reserve(trkSelectors.size());
00059   chi2n_par_.reserve(trkSelectors.size());
00060   chi2n_no1Dmod_par_.reserve(trkSelectors.size());
00061   d0_par1_.reserve(trkSelectors.size());
00062   dz_par1_.reserve(trkSelectors.size());
00063   d0_par2_.reserve(trkSelectors.size());
00064   dz_par2_.reserve(trkSelectors.size());
00065   applyAdaptedPVCuts_.reserve(trkSelectors.size());
00066   max_d0_.reserve(trkSelectors.size());
00067   max_z0_.reserve(trkSelectors.size());
00068   nSigmaZ_.reserve(trkSelectors.size());
00069   min_layers_.reserve(trkSelectors.size());
00070   min_3Dlayers_.reserve(trkSelectors.size());
00071   max_lostLayers_.reserve(trkSelectors.size());
00072   min_hits_bypass_.reserve(trkSelectors.size());
00073   applyAbsCutsIfNoPV_.reserve(trkSelectors.size());
00074   max_d0NoPV_.reserve(trkSelectors.size());
00075   max_z0NoPV_.reserve(trkSelectors.size());
00076   preFilter_.reserve(trkSelectors.size());
00077   max_relpterr_.reserve(trkSelectors.size());
00078   min_nhits_.reserve(trkSelectors.size());
00079   max_minMissHitOutOrIn_.reserve(trkSelectors.size());
00080   max_lostHitFraction_.reserve(trkSelectors.size());
00081   min_eta_.reserve(trkSelectors.size());
00082   max_eta_.reserve(trkSelectors.size());
00083   useMVA_.reserve(trkSelectors.size());
00084   //mvaReaders_.reserve(trkSelectors.size());
00085   min_MVA_.reserve(trkSelectors.size());
00086   //mvaType_.reserve(trkSelectors.size());
00087 
00088   produces<edm::ValueMap<float> >("MVAVals");
00089 
00090   for ( unsigned int i=0; i<trkSelectors.size(); i++) {
00091 
00092     qualityToSet_.push_back( TrackBase::undefQuality );
00093     // parameters for vertex selection
00094     vtxNumber_.push_back( useVertices_ ? trkSelectors[i].getParameter<int32_t>("vtxNumber") : 0 );
00095     vertexCut_.push_back( useVertices_ ? trkSelectors[i].getParameter<std::string>("vertexCut") : 0);
00096     //  parameters for adapted optimal cuts on chi2 and primary vertex compatibility
00097     res_par_.push_back(trkSelectors[i].getParameter< std::vector<double> >("res_par") );
00098     chi2n_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_par") );
00099     chi2n_no1Dmod_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_no1Dmod_par") );
00100     d0_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par1"));
00101     dz_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par1"));
00102     d0_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par2"));
00103     dz_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par2"));
00104     // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
00105     applyAdaptedPVCuts_.push_back(trkSelectors[i].getParameter<bool>("applyAdaptedPVCuts"));
00106     // Impact parameter absolute cuts.
00107     max_d0_.push_back(trkSelectors[i].getParameter<double>("max_d0"));
00108     max_z0_.push_back(trkSelectors[i].getParameter<double>("max_z0"));
00109     nSigmaZ_.push_back(trkSelectors[i].getParameter<double>("nSigmaZ"));
00110     // Cuts on numbers of layers with hits/3D hits/lost hits.
00111     min_layers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumberLayers") );
00112     min_3Dlayers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumber3DLayers") );
00113     max_lostLayers_.push_back(trkSelectors[i].getParameter<uint32_t>("maxNumberLostLayers"));
00114     min_hits_bypass_.push_back(trkSelectors[i].getParameter<uint32_t>("minHitsToBypassChecks"));
00115     // Flag to apply absolute cuts if no PV passes the selection
00116     applyAbsCutsIfNoPV_.push_back(trkSelectors[i].getParameter<bool>("applyAbsCutsIfNoPV"));
00117     keepAllTracks_.push_back( trkSelectors[i].getParameter<bool>("keepAllTracks")); 
00118     max_relpterr_.push_back(trkSelectors[i].getParameter<double>("max_relpterr"));
00119     min_nhits_.push_back(trkSelectors[i].getParameter<uint32_t>("min_nhits"));
00120     max_minMissHitOutOrIn_.push_back(
00121         trkSelectors[i].existsAs<int32_t>("max_minMissHitOutOrIn") ? 
00122         trkSelectors[i].getParameter<int32_t>("max_minMissHitOutOrIn") : 99);
00123     max_lostHitFraction_.push_back(
00124         trkSelectors[i].existsAs<double>("max_lostHitFraction") ?
00125         trkSelectors[i].getParameter<double>("max_lostHitFraction") : 1.0);
00126     min_eta_.push_back(trkSelectors[i].existsAs<double>("min_eta") ?
00127         trkSelectors[i].getParameter<double>("min_eta"):-9999);
00128     max_eta_.push_back(trkSelectors[i].existsAs<double>("max_eta") ?
00129         trkSelectors[i].getParameter<double>("max_eta"):9999);
00130   
00131     setQualityBit_.push_back( false );
00132     std::string qualityStr = trkSelectors[i].getParameter<std::string>("qualityBit");
00133     if (qualityStr != "") {
00134       setQualityBit_[i] = true;
00135       qualityToSet_[i]  = TrackBase::qualityByName(trkSelectors[i].getParameter<std::string>("qualityBit"));
00136     }
00137   
00138     if (setQualityBit_[i] && (qualityToSet_[i] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
00139     "You can't set the quality bit " << trkSelectors[i].getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
00140 
00141     if (applyAbsCutsIfNoPV_[i]) {
00142       max_d0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_d0NoPV"));
00143       max_z0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_z0NoPV"));
00144     }
00145     else{//dummy values
00146       max_d0NoPV_.push_back(0.);
00147       max_z0NoPV_.push_back(0.);
00148     }
00149   
00150     name_.push_back( trkSelectors[i].getParameter<std::string>("name") );
00151 
00152     preFilter_[i]=trkSelectors.size(); // no prefilter
00153 
00154     std::string pfName=trkSelectors[i].getParameter<std::string>("preFilterName");
00155     if (pfName!="") {
00156       bool foundPF=false;
00157       for ( unsigned int j=0; j<i; j++) 
00158         if (name_[j]==pfName ) {
00159           foundPF=true;
00160           preFilter_[i]=j;
00161         }
00162       if ( !foundPF)
00163         throw cms::Exception("Configuration") << "Invalid prefilter name in MultiTrackSelector " 
00164                                               << trkSelectors[i].getParameter<std::string>("preFilterName");
00165           
00166     }
00167 
00168     //    produces<std::vector<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
00169     produces<edm::ValueMap<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
00170     if(useAnyMVA_){
00171       bool thisMVA = false;
00172       if(trkSelectors[i].exists("useMVA"))thisMVA = trkSelectors[i].getParameter<bool>("useMVA");
00173       useMVA_.push_back(thisMVA);
00174       if(thisMVA){
00175         double minVal = -1;
00176         if(trkSelectors[i].exists("minMVA"))minVal = trkSelectors[i].getParameter<double>("minMVA");
00177         min_MVA_.push_back(minVal);
00178 
00179       }else{
00180         min_MVA_.push_back(-9999.0);
00181       }
00182     }else{
00183       min_MVA_.push_back(-9999.0);
00184     }
00185 
00186   }
00187 }
00188 
00189 MultiTrackSelector::~MultiTrackSelector() {
00190   if(gbrVals_)delete [] gbrVals_;
00191   if(!useForestFromDB_ && forest_)delete forest_;
00192 }
00193 
00194 void MultiTrackSelector::produce( edm::Event& evt, const edm::EventSetup& es ) 
00195 {
00196   using namespace std; 
00197   using namespace edm;
00198   using namespace reco;
00199 
00200   // Get tracks 
00201   Handle<TrackCollection> hSrcTrack;
00202   evt.getByLabel( src_, hSrcTrack );
00203   const TrackCollection& srcTracks(*hSrcTrack);
00204 
00205   // looking for the beam spot
00206   edm::Handle<reco::BeamSpot> hBsp;
00207   evt.getByLabel(beamspot_, hBsp);
00208   const reco::BeamSpot& vertexBeamSpot(*hBsp);
00209 
00210         
00211   // Select good primary vertices for use in subsequent track selection
00212   edm::Handle<reco::VertexCollection> hVtx;
00213   if (useVertices_) evt.getByLabel(vertices_, hVtx);
00214 
00215   unsigned int trkSize=srcTracks.size();
00216   std::vector<int> selTracksSave( qualityToSet_.size()*trkSize,0);
00217 
00218   processMVA(evt,es);
00219 
00220   for (unsigned int i=0; i<qualityToSet_.size(); i++) {  
00221     std::vector<int> selTracks(trkSize,0);
00222     auto_ptr<edm::ValueMap<int> > selTracksValueMap = auto_ptr<edm::ValueMap<int> >(new edm::ValueMap<int>);
00223     edm::ValueMap<int>::Filler filler(*selTracksValueMap);
00224 
00225     std::vector<Point> points;
00226     std::vector<float> vterr, vzerr;
00227     if (useVertices_) selectVertices(i,*hVtx, points, vterr, vzerr);
00228 
00229     // Loop over tracks
00230     size_t current = 0;
00231     for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
00232       const Track & trk = * it;
00233       // Check if this track passes cuts
00234 
00235       LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
00236 
00237       //already removed
00238       bool ok=true;
00239       if (preFilter_[i]<i && selTracksSave[preFilter_[i]*trkSize+current] < 0) {
00240         selTracks[current]=-1;
00241         ok=false;
00242         if ( !keepAllTracks_[i]) 
00243           continue;
00244       }
00245       else {
00246         double mvaVal = 0;
00247         if(useAnyMVA_)mvaVal = mvaVals_[current];
00248         ok = select(i,vertexBeamSpot, trk, points, vterr, vzerr,mvaVal);
00249         if (!ok) { 
00250           LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
00251           if (!keepAllTracks_[i]) { 
00252             selTracks[current]=-1;
00253             continue;
00254           }
00255         }
00256         else
00257           LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
00258       }
00259 
00260       if (preFilter_[i]<i ) {
00261         selTracks[current]=selTracksSave[preFilter_[i]*trkSize+current];
00262       }
00263       else {
00264         selTracks[current]=trk.qualityMask();
00265       }
00266       if ( ok && setQualityBit_[i]) {
00267         selTracks[current]= (selTracks[current] | (1<<qualityToSet_[i]));
00268         if (!points.empty()) {
00269           if (qualityToSet_[i]==TrackBase::loose) {
00270             selTracks[current]=(selTracks[current] | (1<<TrackBase::looseSetWithPV));
00271           }
00272           else if (qualityToSet_[i]==TrackBase::highPurity) {
00273             selTracks[current]=(selTracks[current] | (1<<TrackBase::highPuritySetWithPV));
00274           }
00275         }
00276       }
00277     }
00278     for ( unsigned int j=0; j< trkSize; j++ ) selTracksSave[j+i*trkSize]=selTracks[j];
00279     filler.insert(hSrcTrack, selTracks.begin(),selTracks.end());
00280     filler.fill();
00281 
00282     //    evt.put(selTracks,name_[i]);
00283     evt.put(selTracksValueMap,name_[i]);
00284   }
00285 }
00286 
00287 
00288  bool MultiTrackSelector::select(unsigned int tsNum, 
00289                                  const reco::BeamSpot &vertexBeamSpot, 
00290                                  const reco::Track &tk, 
00291                                  const std::vector<Point> &points,
00292                                  std::vector<float> &vterr,
00293                                  std::vector<float> &vzerr,
00294                                  double mvaVal) {
00295   // Decide if the given track passes selection cuts.
00296 
00297   using namespace std; 
00298   
00299   if(tk.found()>=min_hits_bypass_[tsNum]) return true;
00300   if ( tk.ndof() < 1E-5 ) return false;
00301 
00302   // Cuts on numbers of layers with hits/3D hits/lost hits.
00303   uint32_t nlayers     = tk.hitPattern().trackerLayersWithMeasurement();
00304   uint32_t nlayers3D   = tk.hitPattern().pixelLayersWithMeasurement() +
00305     tk.hitPattern().numberOfValidStripLayersWithMonoAndStereo();
00306   uint32_t nlayersLost = tk.hitPattern().trackerLayersWithoutMeasurement();
00307   LogDebug("TrackSelection") << "cuts on nlayers: " << nlayers << " " << nlayers3D << " " << nlayersLost << " vs " 
00308                              << min_layers_[tsNum] << " " << min_3Dlayers_[tsNum] << " " << max_lostLayers_[tsNum];
00309   if (nlayers < min_layers_[tsNum]) return false;
00310   if (nlayers3D < min_3Dlayers_[tsNum]) return false;
00311   if (nlayersLost > max_lostLayers_[tsNum]) return false;
00312   LogTrace("TrackSelection") << "cuts on nlayers passed";
00313 
00314   float chi2n =  tk.normalizedChi2();
00315   float chi2n_no1Dmod = chi2n;
00316 
00317   int count1dhits = 0;
00318   for (trackingRecHit_iterator ith = tk.recHitsBegin(), edh = tk.recHitsEnd(); ith != edh; ++ith) {
00319     const TrackingRecHit * hit = ith->get();
00320     if (hit->isValid()) {
00321       if (typeid(*hit) == typeid(SiStripRecHit1D)) ++count1dhits;
00322     }
00323   }
00324   if (count1dhits > 0) {
00325     float chi2 = tk.chi2();
00326     float ndof = tk.ndof();
00327     chi2n = (chi2+count1dhits)/float(ndof+count1dhits);
00328   }
00329   // For each 1D rechit, the chi^2 and ndof is increased by one.  This is a way of retaining approximately
00330   // the same normalized chi^2 distribution as with 2D rechits.
00331   if (chi2n > chi2n_par_[tsNum]*nlayers) return false;
00332 
00333   if (chi2n_no1Dmod > chi2n_no1Dmod_par_[tsNum]*nlayers) return false;
00334 
00335   // Get track parameters
00336   float pt = std::max(float(tk.pt()),0.000001f);
00337   float eta = tk.eta();
00338   if (eta<min_eta_[tsNum] || eta>max_eta_[tsNum]) return false;
00339 
00340   //cuts on relative error on pt and number of valid hits
00341   float relpterr = float(tk.ptError())/pt;
00342   uint32_t nhits = tk.numberOfValidHits();
00343   if(relpterr > max_relpterr_[tsNum]) return false;
00344   if(nhits < min_nhits_[tsNum]) return false;
00345 
00346   int lostIn = tk.trackerExpectedHitsInner().numberOfLostTrackerHits();
00347   int lostOut = tk.trackerExpectedHitsOuter().numberOfLostTrackerHits();
00348   int minLost = std::min(lostIn,lostOut);
00349   if (minLost > max_minMissHitOutOrIn_[tsNum]) return false;
00350   float lostMidFrac = tk.numberOfLostHits() / (tk.numberOfValidHits() + tk.numberOfLostHits());
00351   if (lostMidFrac > max_lostHitFraction_[tsNum]) return false;
00352 
00353 
00354 
00356   //Adding the MVA selection before vertex cuts//
00358 
00359   if(useAnyMVA_ && useMVA_[tsNum]){
00360     if(mvaVal < min_MVA_[tsNum])return false;
00361   }
00362 
00364   //End of MVA selection section//
00366 
00367   //other track parameters
00368   float d0 = -tk.dxy(vertexBeamSpot.position()), d0E =  tk.d0Error(),
00369     dz = tk.dz(vertexBeamSpot.position()), dzE =  tk.dzError();
00370 
00371   // parametrized d0 resolution for the track pt
00372   float nomd0E = sqrt(res_par_[tsNum][0]*res_par_[tsNum][0]+(res_par_[tsNum][1]/pt)*(res_par_[tsNum][1]/pt));
00373   // parametrized z0 resolution for the track pt and eta
00374   float nomdzE = nomd0E*(std::cosh(eta));
00375 
00376   float dzCut = min( pow(dz_par1_[tsNum][0]*nlayers,dz_par1_[tsNum][1])*nomdzE, 
00377                       pow(dz_par2_[tsNum][0]*nlayers,dz_par2_[tsNum][1])*dzE );
00378   float d0Cut = min( pow(d0_par1_[tsNum][0]*nlayers,d0_par1_[tsNum][1])*nomd0E, 
00379                       pow(d0_par2_[tsNum][0]*nlayers,d0_par2_[tsNum][1])*d0E );
00380 
00381 
00382   // ---- PrimaryVertex compatibility cut
00383   bool primaryVertexZCompatibility(false);   
00384   bool primaryVertexD0Compatibility(false);   
00385 
00386   if (points.empty()) { //If not primaryVertices are reconstructed, check just the compatibility with the BS
00387     //z0 within (n sigma + dzCut) of the beam spot z, if no good vertex is found
00388     if ( abs(dz) < hypot(vertexBeamSpot.sigmaZ()*nSigmaZ_[tsNum],dzCut) ) primaryVertexZCompatibility = true;  
00389     // d0 compatibility with beam line
00390     if (abs(d0) < d0Cut) primaryVertexD0Compatibility = true;     
00391   }
00392 
00393   int iv=0;
00394   for (std::vector<Point>::const_iterator point = points.begin(), end = points.end(); point != end; ++point) {
00395     LogTrace("TrackSelection") << "Test track w.r.t. vertex with z position " << point->z();
00396     if(primaryVertexZCompatibility && primaryVertexD0Compatibility) break;
00397     float dzPV = tk.dz(*point); //re-evaluate the dz with respect to the vertex position
00398     float d0PV = tk.dxy(*point); //re-evaluate the dxy with respect to the vertex position
00399     if(useVtxError_){
00400        float dzErrPV = sqrt(dzE*dzE+vzerr[iv]*vzerr[iv]); // include vertex error in z
00401        float d0ErrPV = sqrt(d0E*d0E+vterr[iv]*vterr[iv]); // include vertex error in xy
00402        iv++;
00403        if (abs(dzPV) < dz_par1_[tsNum][0]*pow(nlayers,dz_par1_[tsNum][1])*nomdzE &&
00404            abs(dzPV) < dz_par2_[tsNum][0]*pow(nlayers,dz_par2_[tsNum][1])*dzErrPV &&
00405            abs(dzPV) < max_z0_[tsNum])  primaryVertexZCompatibility = true;
00406        if (abs(d0PV) < d0_par1_[tsNum][0]*pow(nlayers,d0_par1_[tsNum][1])*nomd0E &&
00407            abs(d0PV) < d0_par2_[tsNum][0]*pow(nlayers,d0_par2_[tsNum][1])*d0ErrPV &&
00408            abs(d0PV) < max_d0_[tsNum]) primaryVertexD0Compatibility = true; 
00409     }else{
00410        if (abs(dzPV) < dzCut)  primaryVertexZCompatibility = true;
00411        if (abs(d0PV) < d0Cut) primaryVertexD0Compatibility = true;     
00412     }
00413     LogTrace("TrackSelection") << "distances " << dzPV << " " << d0PV << " vs " << dzCut << " " << d0Cut;
00414   }
00415 
00416   if (points.empty() && applyAbsCutsIfNoPV_[tsNum]) {
00417     if ( abs(dz) > max_z0NoPV_[tsNum] || abs(d0) > max_d0NoPV_[tsNum]) return false;
00418   }  else {
00419     // Absolute cuts on all tracks impact parameters with respect to beam-spot.
00420     // If BS is not compatible, verify if at least the reco-vertex is compatible (useful for incorrect BS settings)
00421     if (abs(d0) > max_d0_[tsNum] && !primaryVertexD0Compatibility) return false;
00422     LogTrace("TrackSelection") << "absolute cuts on d0 passed";
00423     if (abs(dz) > max_z0_[tsNum] && !primaryVertexZCompatibility) return false;
00424     LogTrace("TrackSelection") << "absolute cuts on dz passed";
00425   }
00426 
00427   LogTrace("TrackSelection") << "cuts on PV: apply adapted PV cuts? " << applyAdaptedPVCuts_[tsNum] 
00428                              << " d0 compatibility? " << primaryVertexD0Compatibility  
00429                              << " z compatibility? " << primaryVertexZCompatibility ;
00430 
00431   if (applyAdaptedPVCuts_[tsNum]) {
00432     return (primaryVertexD0Compatibility && primaryVertexZCompatibility);
00433   } else {
00434     return true;     
00435   }
00436 
00437 }
00438 
00439  void MultiTrackSelector::selectVertices(unsigned int tsNum, 
00440                                          const reco::VertexCollection &vtxs, 
00441                                          std::vector<Point> &points,
00442                                          std::vector<float> &vterr, 
00443                                          std::vector<float> &vzerr) {
00444   // Select good primary vertices
00445   using namespace reco;
00446   int32_t toTake = vtxNumber_[tsNum]; 
00447   for (VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
00448 
00449     LogDebug("SelectVertex") << " select vertex with z position " << it->z() << " " 
00450                              << it->chi2() << " " << it->ndof() << " " << TMath::Prob(it->chi2(), static_cast<int32_t>(it->ndof()));
00451     Vertex vtx = *it;
00452     bool pass = vertexCut_[tsNum]( vtx );
00453     if( pass ) { 
00454       points.push_back(it->position()); 
00455       vterr.push_back(sqrt(it->yError()*it->xError()));
00456       vzerr.push_back(it->zError());
00457       LogTrace("SelectVertex") << " SELECTED vertex with z position " << it->z();
00458       toTake--; if (toTake == 0) break;
00459     }
00460   }
00461 }
00462 
00463 void MultiTrackSelector::processMVA(edm::Event& evt, const edm::EventSetup& es)
00464 {
00465 
00466   using namespace std; 
00467   using namespace edm;
00468   using namespace reco;
00469 
00470   // Get tracks 
00471   Handle<TrackCollection> hSrcTrack;
00472   evt.getByLabel( src_, hSrcTrack );
00473   const TrackCollection& srcTracks(*hSrcTrack);
00474 
00475   auto_ptr<edm::ValueMap<float> >mvaValValueMap = auto_ptr<edm::ValueMap<float> >(new edm::ValueMap<float>);
00476   edm::ValueMap<float>::Filler mvaFiller(*mvaValValueMap);
00477 
00478   mvaVals_.clear();
00479 
00480   if(!useAnyMVA_){
00481     size_t current = 0;
00482     for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
00483       mvaVals_.push_back(-99.0);
00484     }
00485     mvaFiller.insert(hSrcTrack,mvaVals_.begin(),mvaVals_.end());
00486     mvaFiller.fill();
00487     evt.put(mvaValValueMap,"MVAVals");
00488     return;
00489   }
00490 
00491   if(!forest_){
00492     if(useForestFromDB_){
00493       edm::ESHandle<GBRForest> forestHandle;
00494       es.get<GBRWrapperRcd>().get(forestLabel_,forestHandle);
00495       forest_ = (GBRForest*)forestHandle.product();
00496     }else{
00497       TFile gbrfile(dbFileName_.c_str());
00498       forest_ = (GBRForest*)gbrfile.Get(forestLabel_.c_str());
00499     }
00500   }
00501     
00502 
00503 
00504   size_t current = 0;
00505   for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
00506     const Track & trk = * it;
00507     tmva_ndof_ = trk.ndof();
00508     tmva_nlayers_ = trk.hitPattern().trackerLayersWithMeasurement();
00509     tmva_nlayers3D_ = trk.hitPattern().pixelLayersWithMeasurement() + trk.hitPattern().numberOfValidStripLayersWithMonoAndStereo();
00510     tmva_nlayerslost_ = trk.hitPattern().trackerLayersWithoutMeasurement();
00511     float chi2n =  trk.normalizedChi2();
00512     float chi2n_no1Dmod = chi2n;
00513     
00514     int count1dhits = 0;
00515     for (trackingRecHit_iterator ith = trk.recHitsBegin(), edh = trk.recHitsEnd(); ith != edh; ++ith) {
00516       const TrackingRecHit * hit = ith->get();
00517       if (hit->isValid()) {
00518         if (typeid(*hit) == typeid(SiStripRecHit1D)) ++count1dhits;
00519       }
00520     }
00521     if (count1dhits > 0) {
00522       float chi2 = trk.chi2();
00523       float ndof = trk.ndof();
00524       chi2n = (chi2+count1dhits)/float(ndof+count1dhits);
00525     }
00526     tmva_chi2n_ = chi2n;
00527     tmva_chi2n_no1dmod_ = chi2n_no1Dmod;
00528     tmva_eta_ = trk.eta();
00529     tmva_relpterr_ = float(trk.ptError())/std::max(float(trk.pt()),0.000001f);
00530     tmva_nhits_ = trk.numberOfValidHits();
00531     int lostIn = trk.trackerExpectedHitsInner().numberOfLostTrackerHits();
00532     int lostOut = trk.trackerExpectedHitsOuter().numberOfLostTrackerHits();
00533     int minLost = std::min(lostIn,lostOut);      tmva_minlost_ = minLost;
00534     tmva_lostmidfrac_ = trk.numberOfLostHits() / (trk.numberOfValidHits() + trk.numberOfLostHits());
00535 
00536     gbrVals_[0] = tmva_lostmidfrac_;
00537     gbrVals_[1] = tmva_minlost_;
00538     gbrVals_[2] = tmva_nhits_;
00539     gbrVals_[3] = tmva_relpterr_;
00540     gbrVals_[4] = tmva_eta_;
00541     gbrVals_[5] = tmva_chi2n_no1dmod_;
00542     gbrVals_[6] = tmva_chi2n_;
00543     gbrVals_[7] = tmva_nlayerslost_;
00544     gbrVals_[8] = tmva_nlayers3D_;
00545     gbrVals_[9] = tmva_nlayers_;
00546     gbrVals_[10] = tmva_ndof_;
00547 
00548     double gbrVal = forest_->GetClassifier(gbrVals_);
00549     mvaVals_.push_back((float)gbrVal);
00550   }
00551   mvaFiller.insert(hSrcTrack,mvaVals_.begin(),mvaVals_.end());
00552   mvaFiller.fill();
00553   evt.put(mvaValValueMap,"MVAVals");
00554 
00555 }