CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoTauTag/TauTagTools/src/Discriminants.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/TauTagTools/interface/Discriminants.h"
00002 
00003 namespace PFTauDiscriminants {
00004 using namespace std;
00005 
00006 typedef reco::Particle::LorentzVector LorentzVector;
00007 
00008 void
00009 DecayMode::doComputation(PFTauDiscriminantManager* input, std::vector<int>& result)
00010 {
00011    // convert to int for TTree
00012    result.push_back(static_cast<int>(input->getDecayMode()->getDecayMode()));
00013 }
00014 
00015 void
00016 OutlierNCharged::doComputation(PFTauDiscriminantManager* input, std::vector<int>& result)
00017 {
00018    size_t output = 0;
00019    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00020    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00021          iObject != theOutlierObjects.end();
00022          ++iObject)
00023    {
00024       const reco::Candidate* currentObject = *iObject;
00025       if (currentObject && currentObject->charge() != 0 )
00026          output++;
00027    }
00028    // convert to int for TTree
00029    result.push_back(static_cast<int>(output));
00030 }
00031 
00032 void
00033 OutlierN::doComputation(PFTauDiscriminantManager* input, std::vector<int>& result)
00034 {
00035    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00036    size_t output = theOutlierObjects.size();
00037    // convert to int for TTree
00038    result.push_back(static_cast<int>(output));
00039 }
00040 
00041 void
00042 Pt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00043 {
00044    result.push_back(input->getDecayMode()->pt());
00045 }
00046 
00047 void
00048 Eta::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00049 {
00050    result.push_back(std::abs(input->getDecayMode()->eta()));
00051 }
00052 
00053 void
00054 MainTrackPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00055 {
00056    const reco::Candidate* theMainTrack = input->mainTrack();
00057    if (theMainTrack)
00058       result.push_back(theMainTrack->pt());
00059 }
00060 
00061 void
00062 MainTrackAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00063 {
00064    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00065 
00066    const reco::Candidate* theMainTrack = input->mainTrack();
00067 
00068    DeltaR<math::XYZVector> myDRComputer;
00069 
00070    if (theMainTrack)
00071       result.push_back(myDRComputer(theMainTrack->momentum(), signalObjectsAxis));
00072 }
00073 
00074 void
00075 TrackPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00076 {
00077    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00078 
00079    const reco::Candidate* theMainTrack = input->mainTrack();
00080 
00081    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00082          iObject != theSignalObjects.end();
00083          ++iObject)
00084    {
00085       const reco::Candidate* currentObject = *iObject;
00086       if (currentObject->charge() && currentObject != theMainTrack)
00087          result.push_back(currentObject->pt());
00088    }
00089 }
00090 
00091 void
00092 PiZeroPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00093 {
00094    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00095 
00096    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00097          iObject != theSignalObjects.end();
00098          ++iObject)
00099    {
00100       const reco::Candidate* currentObject = *iObject;
00101       if (!currentObject->charge())
00102          result.push_back(currentObject->pt());
00103    }
00104 }
00105 
00106 void
00107 FilteredObjectPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00108 {
00109    const reco::PFTauDecayMode* theDecayMode = input->getDecayMode();
00110    if (!theDecayMode)
00111       return;
00112 
00113    const std::vector<const reco::Candidate*> theFilteredObjects = theDecayMode->filteredObjectCandidates();
00114 
00115    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theFilteredObjects.begin();
00116                                                       iObject != theFilteredObjects.end();
00117                                                     ++iObject)
00118    {
00119       const reco::Candidate* currentObject = *iObject;
00120       result.push_back(currentObject->pt());
00121    }
00122 }
00123 
00124 void
00125 GammaOccupancy::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00126 {
00127    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00128 
00129    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00130          iObject != theSignalObjects.end();
00131          ++iObject)
00132    {
00133       const reco::Candidate* currentObject = *iObject;
00134       if (!currentObject->charge())
00135          result.push_back(input->getLeafDaughters(currentObject).size());
00136    }
00137 }
00138 
00139 void
00140 GammaPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00141 {
00142    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00143 
00144    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00145          iObject != theSignalObjects.end();
00146          ++iObject)
00147    {
00148       const reco::Candidate* currentObject = *iObject;
00149       if (!currentObject->charge())
00150       {
00151          std::vector<const reco::Candidate*> daughters = input->getLeafDaughters(currentObject);
00152          for(std::vector<const reco::Candidate*>::const_iterator iDaughter  = daughters.begin();
00153                                                             iDaughter != daughters.end();
00154                                                           ++iDaughter)
00155          {
00156             result.push_back((*iDaughter)->pt());
00157          }
00158       }
00159    }
00160 }
00161 
00162 
00163 void
00164 TrackAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00165 {
00166    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00167 
00168    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00169 
00170    const reco::Candidate* theMainTrack = input->mainTrack();
00171 
00172    DeltaR<math::XYZVector> myDRComputer;
00173 
00174    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00175          iObject != theSignalObjects.end();
00176          ++iObject)
00177    {
00178       const reco::Candidate* currentObject = *iObject;
00179       if (currentObject->charge() && currentObject != theMainTrack)
00180          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00181    }
00182 }
00183 
00184 void
00185 PiZeroAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00186 {
00187    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00188 
00189    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00190 
00191    DeltaR<math::XYZVector> myDRComputer;
00192 
00193    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00194          iObject != theSignalObjects.end();
00195          ++iObject)
00196    {
00197       const reco::Candidate* currentObject = *iObject;
00198       if (!currentObject->charge())
00199          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00200    }
00201 }
00202 
00203 void
00204 Dalitz::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00205 {
00206    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00207    const reco::Candidate* theMainTrack = input->mainTrack();
00208    if (!theMainTrack)
00209       return;
00210    LorentzVector mainTrackFourVector = theMainTrack->p4();
00211 
00212    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00213          iObject != theSignalObjects.end();
00214          ++iObject)
00215    {
00216       const reco::Candidate* currentObject = *iObject;
00217       if (currentObject != theMainTrack)
00218       {
00219          LorentzVector combinedFourVector = mainTrackFourVector + currentObject->p4();
00220          result.push_back(combinedFourVector.mass());
00221       }
00222    }
00223 }
00224 
00225 // takes invariant mass of all objects in signal cone
00226 void
00227 InvariantMassOfSignal::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00228 {
00229    result.push_back(input->getDecayMode()->mass());
00230 }
00231 
00232 // takes invariant mass of all objects in signal cone + Filtered objects
00233 void
00234 InvariantMassOfSignalWithFiltered::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00235 {
00236    LorentzVector signalObjects = input->getDecayMode()->p4();
00237    signalObjects += input->getDecayMode()->filteredObjects().p4();
00238    result.push_back(signalObjects.M());
00239 }
00240 
00241 // returns vector of invariant masses of larger and larger subsets of all signal objects e.g. result[2] is
00242 // the invariant mass of the lead track with the next highest Pt object
00243 
00244 void
00245 InvariantMass::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00246 {
00247    const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00248    const reco::Candidate* theMainTrack = input->mainTrack();
00249    if (!theMainTrack)
00250       return;
00251    LorentzVector fourVectorSoFar = theMainTrack->p4();
00252 
00253    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00254          iObject != theSignalObjects.end();
00255          ++iObject)
00256    {
00257       const reco::Candidate* currentObject = *iObject;
00258       if (currentObject != theMainTrack)
00259       {
00260          fourVectorSoFar += currentObject->p4();
00261          result.push_back(fourVectorSoFar.mass());
00262       }
00263    }
00264 }
00265 
00266 void
00267 OutlierPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00268 {
00269    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00270    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00271          iObject != theOutlierObjects.end();
00272          ++iObject)
00273    {
00274       const reco::Candidate* currentObject = *iObject;
00275       if (currentObject)
00276          result.push_back(currentObject->pt());
00277    }
00278 }
00279 
00280 void
00281 OutlierSumPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00282 {
00283    LorentzVector totalFourVector;
00284    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00285    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00286          iObject != theOutlierObjects.end();
00287          ++iObject)
00288    {
00289       const reco::Candidate* currentObject = *iObject;
00290       if (currentObject)
00291          totalFourVector += currentObject->p4();
00292    }
00293    result.push_back(totalFourVector.pt());
00294 }
00295 
00296 void
00297 OutlierMass::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00298 {
00299    LorentzVector totalFourVector;
00300    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00301    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00302          iObject != theOutlierObjects.end();
00303          ++iObject)
00304    {
00305       const reco::Candidate* currentObject = *iObject;
00306       if (currentObject)
00307          totalFourVector += currentObject->p4();
00308    }
00309    result.push_back(totalFourVector.M());
00310 }
00311 
00312 void
00313 OutlierAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00314 {
00315    const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00316    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00317    DeltaR<math::XYZVector> myDRComputer;
00318    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00319          iObject != theoutlierObjects.end();
00320          ++iObject)
00321    {
00322       const reco::Candidate* currentObject = *iObject;
00323       if (currentObject)
00324          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00325    }
00326 }
00327 
00328 void
00329 ChargedOutlierPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00330 {
00331    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00332    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00333          iObject != theOutlierObjects.end();
00334          ++iObject)
00335    {
00336       const reco::Candidate* currentObject = *iObject;
00337       if (currentObject && currentObject->charge())
00338          result.push_back(currentObject->pt());
00339    }
00340 }
00341 
00342 void
00343 ChargedOutlierSumPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00344 {
00345    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00346    double output = 0.0;
00347    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00348          iObject != theOutlierObjects.end();
00349          ++iObject)
00350    {
00351       const reco::Candidate* currentObject = *iObject;
00352       if (currentObject && currentObject->charge())
00353          output += currentObject->pt();
00354    }
00355    result.push_back(output);
00356 }
00357 
00358 void
00359 ChargedOutlierAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00360 {
00361    const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00362    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00363    DeltaR<math::XYZVector> myDRComputer;
00364    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00365          iObject != theoutlierObjects.end();
00366          ++iObject)
00367    {
00368       const reco::Candidate* currentObject = *iObject;
00369       if (currentObject && currentObject->charge())
00370          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00371    }
00372 }
00373 
00374 void
00375 NeutralOutlierPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00376 {
00377    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00378    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00379          iObject != theOutlierObjects.end();
00380          ++iObject)
00381    {
00382       const reco::Candidate* currentObject = *iObject;
00383       if (currentObject && !currentObject->charge())
00384          result.push_back(currentObject->pt());
00385    }
00386 }
00387 
00388 void
00389 NeutralOutlierSumPt::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00390 {
00391    const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00392    double output = 0.0;
00393    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00394          iObject != theOutlierObjects.end();
00395          ++iObject)
00396    {
00397       const reco::Candidate* currentObject = *iObject;
00398       if (currentObject && !currentObject->charge())
00399          output += currentObject->pt();
00400    }
00401    result.push_back(output);
00402 }
00403 
00404 void
00405 NeutralOutlierAngle::doComputation(PFTauDiscriminantManager* input, std::vector<double>& result)
00406 {
00407    const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00408    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00409    DeltaR<math::XYZVector> myDRComputer;
00410    for(std::vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00411          iObject != theoutlierObjects.end();
00412          ++iObject)
00413    {
00414       const reco::Candidate* currentObject = *iObject;
00415       if (currentObject && !currentObject->charge())
00416          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00417    }
00418 }
00419 
00420 
00421 }
00422