CMS 3D CMS Logo

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, 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, vector<int>& result)
00017 {
00018    size_t output = 0;
00019    const vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00020    for(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 Pt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00034 {
00035    result.push_back(input->getDecayMode()->pt());
00036 }
00037 
00038 void
00039 Eta::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00040 {
00041    result.push_back(input->getDecayMode()->eta());
00042 }
00043 
00044 void
00045 MainTrackPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00046 {
00047    const reco::Candidate* theMainTrack = input->mainTrack();
00048    if (theMainTrack)
00049       result.push_back(theMainTrack->pt());
00050 }
00051 
00052 void
00053 MainTrackAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00054 {
00055    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00056 
00057    const reco::Candidate* theMainTrack = input->mainTrack();
00058 
00059    DeltaR<math::XYZVector> myDRComputer;
00060 
00061    if (theMainTrack)
00062       result.push_back(myDRComputer(theMainTrack->momentum(), signalObjectsAxis));
00063 }
00064 
00065 void
00066 TrackPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00067 {
00068    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00069 
00070    const reco::Candidate* theMainTrack = input->mainTrack();
00071 
00072    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00073          iObject != theSignalObjects.end();
00074          ++iObject)
00075    {
00076       const reco::Candidate* currentObject = *iObject;
00077       if (currentObject->charge() && currentObject != theMainTrack)
00078          result.push_back(currentObject->pt());
00079    }
00080 }
00081 
00082 void
00083 PiZeroPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00084 {
00085    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00086 
00087    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00088          iObject != theSignalObjects.end();
00089          ++iObject)
00090    {
00091       const reco::Candidate* currentObject = *iObject;
00092       if (!currentObject->charge())
00093          result.push_back(currentObject->pt());
00094    }
00095 }
00096 
00097 void
00098 TrackAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00099 {
00100    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00101 
00102    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00103 
00104    const reco::Candidate* theMainTrack = input->mainTrack();
00105 
00106    DeltaR<math::XYZVector> myDRComputer;
00107 
00108    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00109          iObject != theSignalObjects.end();
00110          ++iObject)
00111    {
00112       const reco::Candidate* currentObject = *iObject;
00113       if (currentObject->charge() && currentObject != theMainTrack)
00114          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00115    }
00116 }
00117 
00118 void
00119 PiZeroAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00120 {
00121    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00122 
00123    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00124 
00125    DeltaR<math::XYZVector> myDRComputer;
00126 
00127    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00128          iObject != theSignalObjects.end();
00129          ++iObject)
00130    {
00131       const reco::Candidate* currentObject = *iObject;
00132       if (!currentObject->charge())
00133          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00134    }
00135 }
00136 
00137 void
00138 Dalitz::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00139 {
00140    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00141    const reco::Candidate* theMainTrack = input->mainTrack();
00142    if (!theMainTrack)
00143       return;
00144    LorentzVector mainTrackFourVector = theMainTrack->p4();
00145 
00146    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00147          iObject != theSignalObjects.end();
00148          ++iObject)
00149    {
00150       const reco::Candidate* currentObject = *iObject;
00151       if (currentObject != theMainTrack)
00152       {
00153          LorentzVector combinedFourVector = mainTrackFourVector + currentObject->p4();
00154          result.push_back(combinedFourVector.mass());
00155       }
00156    }
00157 }
00158 
00159 // takes invariant mass of all objects in signal cone
00160 void
00161 InvariantMassOfSignal::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00162 {
00163    result.push_back(input->getDecayMode()->mass());
00164 }
00165 
00166 // returns vector of invariant masses of larger and larger subsets of all signal objects e.g. result[2] is
00167 // the invariant mass of the lead track with the next highest Pt object
00168 
00169 void
00170 InvariantMass::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00171 {
00172    const vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
00173    const reco::Candidate* theMainTrack = input->mainTrack();
00174    if (!theMainTrack)
00175       return;
00176    LorentzVector fourVectorSoFar = theMainTrack->p4();
00177 
00178    for(vector<const reco::Candidate*>::const_iterator iObject  = theSignalObjects.begin();
00179          iObject != theSignalObjects.end();
00180          ++iObject)
00181    {
00182       const reco::Candidate* currentObject = *iObject;
00183       if (currentObject != theMainTrack)
00184       {
00185          fourVectorSoFar += currentObject->p4();
00186          result.push_back(fourVectorSoFar.mass());
00187       }
00188    }
00189 }
00190 
00191 void
00192 OutlierPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00193 {
00194    const vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00195    for(vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00196          iObject != theOutlierObjects.end();
00197          ++iObject)
00198    {
00199       const reco::Candidate* currentObject = *iObject;
00200       if (currentObject)
00201          result.push_back(currentObject->pt());
00202    }
00203 }
00204 
00205 void
00206 OutlierAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00207 {
00208    const vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00209    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00210    DeltaR<math::XYZVector> myDRComputer;
00211    for(vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00212          iObject != theoutlierObjects.end();
00213          ++iObject)
00214    {
00215       const reco::Candidate* currentObject = *iObject;
00216       if (currentObject)
00217          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00218    }
00219 }
00220 
00221 void
00222 ChargedOutlierPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00223 {
00224    const vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00225    for(vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00226          iObject != theOutlierObjects.end();
00227          ++iObject)
00228    {
00229       const reco::Candidate* currentObject = *iObject;
00230       if (currentObject && currentObject->charge())
00231          result.push_back(currentObject->pt());
00232    }
00233 }
00234 
00235 void
00236 ChargedOutlierAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00237 {
00238    const vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00239    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00240    DeltaR<math::XYZVector> myDRComputer;
00241    for(vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00242          iObject != theoutlierObjects.end();
00243          ++iObject)
00244    {
00245       const reco::Candidate* currentObject = *iObject;
00246       if (currentObject && currentObject->charge())
00247          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00248    }
00249 }
00250 
00251 void
00252 NeutralOutlierPt::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00253 {
00254    const vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
00255    for(vector<const reco::Candidate*>::const_iterator iObject  = theOutlierObjects.begin();
00256          iObject != theOutlierObjects.end();
00257          ++iObject)
00258    {
00259       const reco::Candidate* currentObject = *iObject;
00260       if (currentObject && !currentObject->charge())
00261          result.push_back(currentObject->pt());
00262    }
00263 }
00264 
00265 void
00266 NeutralOutlierAngle::doComputation(PFTauDiscriminantManager* input, vector<double>& result)
00267 {
00268    const vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
00269    math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
00270    DeltaR<math::XYZVector> myDRComputer;
00271    for(vector<const reco::Candidate*>::const_iterator iObject  = theoutlierObjects.begin();
00272          iObject != theoutlierObjects.end();
00273          ++iObject)
00274    {
00275       const reco::Candidate* currentObject = *iObject;
00276       if (currentObject && !currentObject->charge())
00277          result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
00278    }
00279 }
00280 
00281 
00282 }
00283 
00284 
00285 
00286 
00287 

Generated on Tue Jun 9 17:45:10 2009 for CMSSW by  doxygen 1.5.4