Fork me on GitHub

source: git/modules/ExampleModule.cc@ b4e1ab5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b4e1ab5 was cab38f6, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove svn tags and fix formatting

  • Property mode set to 100644
File size: 3.4 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 ExampleModule
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/ExampleModule.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
55ExampleModule::ExampleModule() :
56 fFormula(0), fItInputArray(0)
57{
58 fFormula = new DelphesFormula;
59}
60
61//------------------------------------------------------------------------------
62
63ExampleModule::~ExampleModule()
64{
65 if(fFormula) delete fFormula;
66}
67
68//------------------------------------------------------------------------------
69
70void ExampleModule::Init()
71{
72 // read parameters
73
74 fIntParam = GetInt("IntParam", 10);
75
76 fDoubleParam = GetDouble("DoubleParam", 1.0);
77
78 fFormula->Compile(GetString("EfficiencyFormula", "0.4"));
79
80 ExRootConfParam param = GetParam("ArrayParam");
81 Long_t i, size;
82 fArrayParam.clear();
83
84 size = param.GetSize();
85 for(i = 0; i < size; ++i)
86 {
87 fArrayParam.push_back(param[i].GetDouble());
88 }
89
90 // import input array(s)
91
92 fInputArray = ImportArray(GetString("InputArray", "FastJetFinder/jets"));
93 fItInputArray = fInputArray->MakeIterator();
94
95 // create output array(s)
96
97 fOutputArray = ExportArray(GetString("OutputArray", "jets"));
98
99}
100
101//------------------------------------------------------------------------------
102
103void ExampleModule::Finish()
104{
105 if(fItInputArray) delete fItInputArray;
106}
107
108//------------------------------------------------------------------------------
109
110void ExampleModule::Process()
111{
112 Candidate *candidate;
113 TLorentzVector candidatePosition, candidateMomentum;
114
115 // loop over all input candidates
116 fItInputArray->Reset();
117 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
118 {
119 candidatePosition = candidate->Position;
120 candidateMomentum = candidate->Momentum;
121
122 // apply an efficency formula
123 if(gRandom->Uniform() <= fFormula->Eval(candidateMomentum.Pt(), candidatePosition.Eta()))
124 {
125 fOutputArray->Add(candidate);
126 }
127 }
128}
129
130//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.