Fork me on GitHub

Changeset 341014c in git for modules/VertexSorter.cc


Ignore:
Timestamp:
Feb 12, 2019, 9:29:17 PM (6 years ago)
Author:
Pavel Demin <pavel-demin@…>
Branches:
ImprovedOutputFile, Timing, llp, master
Children:
6455202
Parents:
45e58be
Message:

apply .clang-format to all .h, .cc and .cpp files

File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/VertexSorter.cc

    r45e58be r341014c  
    1616#include "classes/DelphesPileUpReader.h"
    1717
     18#include "ExRootAnalysis/ExRootClassifier.h"
     19#include "ExRootAnalysis/ExRootFilter.h"
    1820#include "ExRootAnalysis/ExRootResult.h"
    19 #include "ExRootAnalysis/ExRootFilter.h"
    20 #include "ExRootAnalysis/ExRootClassifier.h"
    21 
     21
     22#include "TDatabasePDG.h"
     23#include "TFormula.h"
     24#include "TLorentzVector.h"
    2225#include "TMath.h"
     26#include "TMatrixT.h"
     27#include "TObjArray.h"
     28#include "TRandom3.h"
    2329#include "TString.h"
    24 #include "TFormula.h"
    25 #include "TRandom3.h"
    26 #include "TObjArray.h"
    27 #include "TDatabasePDG.h"
    28 #include "TLorentzVector.h"
    29 #include "TMatrixT.h"
    3030#include "TVector3.h"
    3131
     32#include <algorithm>
     33#include <iostream>
     34#include <map>
     35#include <stdexcept>
     36#include <string>
    3237#include <utility>
    33 #include <algorithm>
    34 #include <stdexcept>
    35 #include <iostream>
    3638#include <vector>
    37 #include <map>
    38 #include <string>
    3939
    4040using namespace std;
    4141
    42 static const Double_t mm  = 1.;
    43 static const Double_t m = 1000.*mm;
    44 static const Double_t ns  = 1.;
    45 static const Double_t s = 1.e+9 *ns;
    46 static const Double_t c_light   = 2.99792458e+8 * m/s;
     42static const Double_t mm = 1.;
     43static const Double_t m = 1000. * mm;
     44static const Double_t ns = 1.;
     45static const Double_t s = 1.e+9 * ns;
     46static const Double_t c_light = 2.99792458e+8 * m / s;
    4747
    4848//------------------------------------------------------------------------------
     
    6868  fItTrackInputArray = fTrackInputArray->MakeIterator();
    6969
    70   if (string (GetString("JetInputArray", "")) != "")
    71     {
    72       fJetInputArray = ImportArray(GetString("JetInputArray", ""));
    73       fItJetInputArray = fJetInputArray->MakeIterator();
    74     }
     70  if(string(GetString("JetInputArray", "")) != "")
     71  {
     72    fJetInputArray = ImportArray(GetString("JetInputArray", ""));
     73    fItJetInputArray = fJetInputArray->MakeIterator();
     74  }
    7575
    7676  // import beamspot
     
    8282  {
    8383    fBeamSpotInputArray = 0;
    84   } 
    85  
     84  }
     85
    8686  fOutputArray = ExportArray(GetString("OutputArray", "clusters"));
    8787
     
    9999//------------------------------------------------------------------------------
    100100
    101 static Bool_t secondDescending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
     101static Bool_t secondDescending(pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
    102102{
    103103  return (pair0.second > pair1.second);
    104104}
    105105
    106 static Bool_t secondAscending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
     106static Bool_t secondAscending(pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
    107107{
    108108  return (pair0.second < pair1.second);
     
    118118  map<Int_t, Double_t> clusterIDToSumPT2;
    119119  map<Int_t, Double_t>::const_iterator itClusterIDToSumPT2;
    120   vector<pair<Int_t, Double_t> > sortedClusterIDs;
    121   vector<pair<Int_t, Double_t> >::const_iterator itSortedClusterIDs;
    122 
    123   for (Int_t iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
    124     {
    125       const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
    126       clusterIDToIndex[cluster.ClusterIndex] = iCluster;
    127       clusterIDToSumPT2[cluster.ClusterIndex] = 0.0;
    128     }
     120  vector<pair<Int_t, Double_t>> sortedClusterIDs;
     121  vector<pair<Int_t, Double_t>>::const_iterator itSortedClusterIDs;
     122
     123  for(Int_t iCluster = 0; iCluster < fInputArray->GetEntries(); iCluster++)
     124  {
     125    const Candidate &cluster = *((Candidate *)fInputArray->At(iCluster));
     126    clusterIDToIndex[cluster.ClusterIndex] = iCluster;
     127    clusterIDToSumPT2[cluster.ClusterIndex] = 0.0;
     128  }
    129129
    130130  if(fMethod == "BTV")
    131     {
    132       if (!fJetInputArray)
    133         {
    134           cout << "BTV PV sorting selected, but no jet collection given!" << endl;
    135           throw 0;
    136         }
    137 
    138       fItTrackInputArray->Reset();
    139       while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
    140         {
    141           if (candidate->Momentum.Pt () < 1.0)
    142             continue;
    143           if (candidate->ClusterIndex < 0)
    144             continue;
    145           TLorentzVector p (candidate->Momentum.Px (), candidate->Momentum.Py (), candidate->Momentum.Pz (), candidate->Momentum.E ());
    146           Bool_t isInJet = false;
    147 
    148           fItJetInputArray->Reset();
    149           while((jetCandidate = static_cast<Candidate*>(fItJetInputArray->Next())))
    150             {
    151               if (jetCandidate->Momentum.Pt () < 30.0)
    152                 continue;
    153               TLorentzVector q (jetCandidate->Momentum.Px (), jetCandidate->Momentum.Py (), jetCandidate->Momentum.Pz (), jetCandidate->Momentum.E ());
    154 
    155               if (p.DeltaR (q) > 0.4)
    156                 continue;
    157               isInJet = true;
    158               break;
    159             }
    160           if (!isInJet)
    161             continue;
    162 
    163           clusterIDToSumPT2.at (candidate->ClusterIndex) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
    164         }
    165 
    166       for (itClusterIDToSumPT2 = clusterIDToSumPT2.begin(); itClusterIDToSumPT2 != clusterIDToSumPT2.end(); ++itClusterIDToSumPT2)
    167         sortedClusterIDs.push_back (make_pair (itClusterIDToSumPT2->first, itClusterIDToSumPT2->second));
    168       sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
    169     }
    170   else if (fMethod == "GenClosest")
    171     {
    172       if (!fBeamSpotInputArray)
    173         {
    174           cout << "GenClosest PV sorting selected, but no beamspot collection given!" << endl;
    175           throw 0;
    176         }
    177       if (!fBeamSpotInputArray->GetSize ())
    178         {
    179           cout << "Beamspot collection is empty!" << endl;
    180           throw 0;
    181         }
    182 
    183       beamSpotCandidate = (Candidate *) fBeamSpotInputArray->At (0);
    184       for (Int_t iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
    185         {
    186           const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
    187           sortedClusterIDs.push_back (make_pair (cluster.ClusterIndex, fabs (cluster.Position.Z () - beamSpotCandidate->Position.Z ())));
    188         }
    189       sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondAscending);
    190     }
    191   else if (fMethod == "GenBest")
    192     {
    193       fItTrackInputArray->Reset();
    194       while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
    195         {
    196           if (candidate->IsPU)
    197             continue;
    198           for (itClusterIDToIndex = clusterIDToIndex.begin(); itClusterIDToIndex != clusterIDToIndex.end(); ++itClusterIDToIndex)
    199             {
    200               if (candidate->ClusterIndex != itClusterIDToIndex->first)
    201                 continue;
    202               clusterIDToSumPT2.at (itClusterIDToIndex->first) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
    203             }
    204         }
    205 
    206       for (itClusterIDToSumPT2 = clusterIDToSumPT2.begin(); itClusterIDToSumPT2 != clusterIDToSumPT2.end(); ++itClusterIDToSumPT2)
    207         sortedClusterIDs.push_back (make_pair (itClusterIDToSumPT2->first, itClusterIDToSumPT2->second));
    208       sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
    209     }
     131  {
     132    if(!fJetInputArray)
     133    {
     134      cout << "BTV PV sorting selected, but no jet collection given!" << endl;
     135      throw 0;
     136    }
     137
     138    fItTrackInputArray->Reset();
     139    while((candidate = static_cast<Candidate *>(fItTrackInputArray->Next())))
     140    {
     141      if(candidate->Momentum.Pt() < 1.0)
     142        continue;
     143      if(candidate->ClusterIndex < 0)
     144        continue;
     145      TLorentzVector p(candidate->Momentum.Px(), candidate->Momentum.Py(), candidate->Momentum.Pz(), candidate->Momentum.E());
     146      Bool_t isInJet = false;
     147
     148      fItJetInputArray->Reset();
     149      while((jetCandidate = static_cast<Candidate *>(fItJetInputArray->Next())))
     150      {
     151        if(jetCandidate->Momentum.Pt() < 30.0)
     152          continue;
     153        TLorentzVector q(jetCandidate->Momentum.Px(), jetCandidate->Momentum.Py(), jetCandidate->Momentum.Pz(), jetCandidate->Momentum.E());
     154
     155        if(p.DeltaR(q) > 0.4)
     156          continue;
     157        isInJet = true;
     158        break;
     159      }
     160      if(!isInJet)
     161        continue;
     162
     163      clusterIDToSumPT2.at(candidate->ClusterIndex) += candidate->Momentum.Pt() * candidate->Momentum.Pt();
     164    }
     165
     166    for(itClusterIDToSumPT2 = clusterIDToSumPT2.begin(); itClusterIDToSumPT2 != clusterIDToSumPT2.end(); ++itClusterIDToSumPT2)
     167      sortedClusterIDs.push_back(make_pair(itClusterIDToSumPT2->first, itClusterIDToSumPT2->second));
     168    sort(sortedClusterIDs.begin(), sortedClusterIDs.end(), secondDescending);
     169  }
     170  else if(fMethod == "GenClosest")
     171  {
     172    if(!fBeamSpotInputArray)
     173    {
     174      cout << "GenClosest PV sorting selected, but no beamspot collection given!" << endl;
     175      throw 0;
     176    }
     177    if(!fBeamSpotInputArray->GetSize())
     178    {
     179      cout << "Beamspot collection is empty!" << endl;
     180      throw 0;
     181    }
     182
     183    beamSpotCandidate = (Candidate *)fBeamSpotInputArray->At(0);
     184    for(Int_t iCluster = 0; iCluster < fInputArray->GetEntries(); iCluster++)
     185    {
     186      const Candidate &cluster = *((Candidate *)fInputArray->At(iCluster));
     187      sortedClusterIDs.push_back(make_pair(cluster.ClusterIndex, fabs(cluster.Position.Z() - beamSpotCandidate->Position.Z())));
     188    }
     189    sort(sortedClusterIDs.begin(), sortedClusterIDs.end(), secondAscending);
     190  }
     191  else if(fMethod == "GenBest")
     192  {
     193    fItTrackInputArray->Reset();
     194    while((candidate = static_cast<Candidate *>(fItTrackInputArray->Next())))
     195    {
     196      if(candidate->IsPU)
     197        continue;
     198      for(itClusterIDToIndex = clusterIDToIndex.begin(); itClusterIDToIndex != clusterIDToIndex.end(); ++itClusterIDToIndex)
     199      {
     200        if(candidate->ClusterIndex != itClusterIDToIndex->first)
     201          continue;
     202        clusterIDToSumPT2.at(itClusterIDToIndex->first) += candidate->Momentum.Pt() * candidate->Momentum.Pt();
     203      }
     204    }
     205
     206    for(itClusterIDToSumPT2 = clusterIDToSumPT2.begin(); itClusterIDToSumPT2 != clusterIDToSumPT2.end(); ++itClusterIDToSumPT2)
     207      sortedClusterIDs.push_back(make_pair(itClusterIDToSumPT2->first, itClusterIDToSumPT2->second));
     208    sort(sortedClusterIDs.begin(), sortedClusterIDs.end(), secondDescending);
     209  }
    210210  else
    211     {
    212       cout << "\"" << fMethod << "\" is not a valid sorting method!" << endl;
    213       cout << "Valid methods are:" << endl;
    214       cout << "  BTV" << endl;
    215       cout << "  GenClosest" << endl;
    216       cout << "  GenBest" << endl;
    217       throw 0;
    218     }
    219   for (itSortedClusterIDs = sortedClusterIDs.begin(); itSortedClusterIDs != sortedClusterIDs.end(); ++itSortedClusterIDs)
    220     {
    221       Candidate *cluster = (Candidate *) fInputArray->At (clusterIDToIndex.at (itSortedClusterIDs->first));
    222       if (fMethod == "BTV")
    223         cluster->BTVSumPT2 = itSortedClusterIDs->second;
    224       else if (fMethod == "GenClosest")
    225         cluster->GenDeltaZ = itSortedClusterIDs->second;
    226       else if (fMethod == "GenBest")
    227         cluster->GenSumPT2 = itSortedClusterIDs->second;
    228       fOutputArray->Add (cluster);
    229     }
    230 }
    231 
    232 //------------------------------------------------------------------------------
     211  {
     212    cout << "\"" << fMethod << "\" is not a valid sorting method!" << endl;
     213    cout << "Valid methods are:" << endl;
     214    cout << "  BTV" << endl;
     215    cout << "  GenClosest" << endl;
     216    cout << "  GenBest" << endl;
     217    throw 0;
     218  }
     219  for(itSortedClusterIDs = sortedClusterIDs.begin(); itSortedClusterIDs != sortedClusterIDs.end(); ++itSortedClusterIDs)
     220  {
     221    Candidate *cluster = (Candidate *)fInputArray->At(clusterIDToIndex.at(itSortedClusterIDs->first));
     222    if(fMethod == "BTV")
     223      cluster->BTVSumPT2 = itSortedClusterIDs->second;
     224    else if(fMethod == "GenClosest")
     225      cluster->GenDeltaZ = itSortedClusterIDs->second;
     226    else if(fMethod == "GenBest")
     227      cluster->GenSumPT2 = itSortedClusterIDs->second;
     228    fOutputArray->Add(cluster);
     229  }
     230}
     231
     232//------------------------------------------------------------------------------
Note: See TracChangeset for help on using the changeset viewer.