Fork me on GitHub

Ignore:
File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/PileUpMerger.cc

    r341014c r3e2bb2b  
    2929#include "classes/DelphesClasses.h"
    3030#include "classes/DelphesFactory.h"
     31#include "classes/DelphesTF2.h"
    3132#include "classes/DelphesPileUpReader.h"
    32 #include "classes/DelphesTF2.h"
    33 
     33
     34#include "ExRootAnalysis/ExRootResult.h"
     35#include "ExRootAnalysis/ExRootFilter.h"
    3436#include "ExRootAnalysis/ExRootClassifier.h"
    35 #include "ExRootAnalysis/ExRootFilter.h"
    36 #include "ExRootAnalysis/ExRootResult.h"
    37 
     37
     38#include "TMath.h"
     39#include "TString.h"
     40#include "TFormula.h"
     41#include "TRandom3.h"
     42#include "TObjArray.h"
    3843#include "TDatabasePDG.h"
    39 #include "TFormula.h"
    4044#include "TLorentzVector.h"
    41 #include "TMath.h"
    42 #include "TObjArray.h"
    43 #include "TRandom3.h"
    44 #include "TString.h"
    4545
    4646#include <algorithm>
     47#include <stdexcept>
    4748#include <iostream>
    4849#include <sstream>
    49 #include <stdexcept>
    5050
    5151using namespace std;
     
    5959}
    6060
     61
    6162//------------------------------------------------------------------------------
    6263
     
    7475  fPileUpDistribution = GetInt("PileUpDistribution", 0);
    7576
    76   fMeanPileUp = GetDouble("MeanPileUp", 10);
     77  fMeanPileUp  = GetDouble("MeanPileUp", 10);
    7778
    7879  fZVertexSpread = GetDouble("ZVertexSpread", 0.15);
     
    134135  dt0 = -1.0e6;
    135136
    136   dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
     137  dt *= c_light*1.0E3; // necessary in order to make t in mm/c
    137138  dz *= 1.0E3; // necessary in order to make z in mm
    138139
     
    149150  vertex = factory->NewCandidate();
    150151
    151   while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
     152  while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
    152153  {
    153154    vx += candidate->Position.X();
     
    158159
    159160    // take postion and time from first stable particle
    160     if(dz0 < -999999.0)
     161    if (dz0 < -999999.0)
    161162      dz0 = z;
    162     if(dt0 < -999999.0)
     163    if (dt0 < -999999.0)
    163164      dt0 = t;
    164165
     
    171172    fParticleOutputArray->Add(candidate);
    172173
    173     if(TMath::Abs(candidate->Charge) > 1.0E-9)
     174    if(TMath::Abs(candidate->Charge) >  1.0E-9)
    174175    {
    175176      nch++;
    176       sumpt2 += pt * pt;
     177      sumpt2 += pt*pt;
    177178      vertex->AddCandidate(candidate);
    178179    }
     
    197198  switch(fPileUpDistribution)
    198199  {
    199   case 0:
    200     numberOfEvents = gRandom->Poisson(fMeanPileUp);
    201     break;
    202   case 1:
    203     numberOfEvents = gRandom->Integer(2 * fMeanPileUp + 1);
    204     break;
    205   case 2:
    206     numberOfEvents = fMeanPileUp;
    207     break;
    208   default:
    209     numberOfEvents = gRandom->Poisson(fMeanPileUp);
    210     break;
     200    case 0:
     201      numberOfEvents = gRandom->Poisson(fMeanPileUp);
     202      break;
     203    case 1:
     204      numberOfEvents = gRandom->Integer(2*fMeanPileUp + 1);
     205      break;
     206    case 2:
     207      numberOfEvents = fMeanPileUp;
     208      break;
     209    default:
     210      numberOfEvents = gRandom->Poisson(fMeanPileUp);
     211      break;
    211212  }
    212213
    213214  allEntries = fReader->GetEntries();
     215
    214216
    215217  for(event = 0; event < numberOfEvents; ++event)
     
    217219    do
    218220    {
    219       entry = TMath::Nint(gRandom->Rndm() * allEntries);
    220     } while(entry >= allEntries);
     221      entry = TMath::Nint(gRandom->Rndm()*allEntries);
     222    }
     223    while(entry >= allEntries);
    221224
    222225    fReader->ReadEntry(entry);
    223226
    224     // --- Pile-up vertex smearing
     227   // --- Pile-up vertex smearing
    225228
    226229    fFunction->GetRandom2(dz, dt);
    227230
    228     dt *= c_light * 1.0E3; // necessary in order to make t in mm/c
     231    dt *= c_light*1.0E3; // necessary in order to make t in mm/c
    229232    dz *= 1.0E3; // necessary in order to make z in mm
    230233
     
    249252
    250253      pdgParticle = pdg->GetParticle(pid);
    251       candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge() / 3.0) : -999;
     254      candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
    252255      candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
    253256
     
    268271
    269272      ++numberOfParticles;
    270       if(TMath::Abs(candidate->Charge) > 1.0E-9)
     273      if(TMath::Abs(candidate->Charge) >  1.0E-9)
    271274      {
    272275        nch++;
    273         sumpt2 += pt * pt;
     276        sumpt2 += pt*pt;
    274277        vertex->AddCandidate(candidate);
    275278      }
     
    296299
    297300    fVertexOutputArray->Add(vertex);
     301
    298302  }
    299303}
Note: See TracChangeset for help on using the changeset viewer.