Fork me on GitHub

source: git/modules/Efficiency.cc@ 5a697dde

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5a697dde was 95aa610, checked in by Michele Selvaggi <michele.selvaggi@…>, 10 years ago

fixed phi, energy dependence in efficiency, smearing, identification, btagging, angular smearing modules

  • Property mode set to 100644
File size: 3.2 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19
20/** \class Efficiency
21 *
22 * Selects candidates from the InputArray according to the efficiency formula.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/Efficiency.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootResult.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootClassifier.h"
37
38#include "TMath.h"
39#include "TString.h"
40#include "TFormula.h"
41#include "TRandom3.h"
42#include "TObjArray.h"
43#include "TDatabasePDG.h"
44#include "TLorentzVector.h"
45
46#include <algorithm>
47#include <stdexcept>
48#include <iostream>
49#include <sstream>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55Efficiency::Efficiency() :
56 fFormula(0), fItInputArray(0)
57{
58 fFormula = new DelphesFormula;
59}
60
61//------------------------------------------------------------------------------
62
63Efficiency::~Efficiency()
64{
65 if(fFormula) delete fFormula;
66}
67
68//------------------------------------------------------------------------------
69
70void Efficiency::Init()
71{
72 // read efficiency formula
73
74 fFormula->Compile(GetString("EfficiencyFormula", "1.0"));
75
76 // import input array
77
78 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
79 fItInputArray = fInputArray->MakeIterator();
80
81 // create output array
82
83 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
84}
85
86//------------------------------------------------------------------------------
87
88void Efficiency::Finish()
89{
90 if(fItInputArray) delete fItInputArray;
91}
92
93//------------------------------------------------------------------------------
94
95void Efficiency::Process()
96{
97 Candidate *candidate;
98 Double_t pt, eta, phi, e;
99
100 fItInputArray->Reset();
101 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
102 {
103 const TLorentzVector &candidatePosition = candidate->Position;
104 const TLorentzVector &candidateMomentum = candidate->Momentum;
105 eta = candidatePosition.Eta();
106 phi = candidatePosition.Phi();
107 pt = candidateMomentum.Pt();
108 e = candidateMomentum.E();
109
110 // apply an efficency formula
111 if(gRandom->Uniform() > fFormula->Eval(pt, eta, phi, e)) continue;
112
113 fOutputArray->Add(candidate);
114 }
115}
116
117//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.