Fork me on GitHub

source: git/modules/IdentificationMap.cc@ d836f45

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

fix GPLv3 header

  • Property mode set to 100644
File size: 5.8 KB
RevLine 
[caba091]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[caba091]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.
[1fa50c2]9 *
[caba091]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.
[1fa50c2]14 *
[caba091]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
[e2a76ae]19
20/** \class IdentificationMap
21 *
[856efb6]22 * Converts particles with some PDG code into another particle, according to parametrized probability as function of pt eta
[e2a76ae]23 given by user.
24 *
25 * $Date: 2014-08-07 14:57:44 +0100 (Thu, 07 Aug 2014) $
26 * $Revision: 905 $
27 *
28 *
29 * \author M. Selvaggi - UCL, Louvain-la-Neuve
30 *
31 */
32
33#include "modules/IdentificationMap.h"
34
35#include "classes/DelphesClasses.h"
36#include "classes/DelphesFactory.h"
37#include "classes/DelphesFormula.h"
38
39#include "ExRootAnalysis/ExRootResult.h"
40#include "ExRootAnalysis/ExRootFilter.h"
41#include "ExRootAnalysis/ExRootClassifier.h"
42
43#include "TMath.h"
44#include "TString.h"
45#include "TFormula.h"
46#include "TRandom3.h"
47#include "TObjArray.h"
48#include "TDatabasePDG.h"
49#include "TLorentzVector.h"
50
51#include <algorithm>
52#include <stdexcept>
53#include <iostream>
54#include <sstream>
55
56using namespace std;
57
58//------------------------------------------------------------------------------
59
60IdentificationMap::IdentificationMap() :
61 fItInputArray(0)
62{
63}
64
65//------------------------------------------------------------------------------
66
67IdentificationMap::~IdentificationMap()
68{
69}
70
71//------------------------------------------------------------------------------
72
73void IdentificationMap::Init()
74{
75 // read efficiency formula
76
77
78 TMisIDMap::iterator itEfficiencyMap;
79 ExRootConfParam param;
80 DelphesFormula *formula;
81 Int_t i, size, pdg;
82
83 // read efficiency formulas
84 param = GetParam("EfficiencyFormula");
85 size = param.GetSize();
86
87 fEfficiencyMap.clear();
88 for(i = 0; i < size/3; ++i)
89 {
90 formula = new DelphesFormula;
91 formula->Compile(param[i*3 + 2].GetString());
92 pdg = param[i*3].GetInt();
93 fEfficiencyMap.insert(make_pair(pdg,make_pair(param[i*3 + 1].GetInt(),formula)));
94
[57e42c6]95 // cout<<param[i*3].GetInt()<<","<<param[i*3+1].GetInt()<<","<<param[i*3 + 2].GetString()<<endl;
[e2a76ae]96
97 }
98
99 // set default efficiency formula
100 itEfficiencyMap = fEfficiencyMap.find(0);
101 if(itEfficiencyMap == fEfficiencyMap.end())
102 {
103 formula = new DelphesFormula;
104 formula->Compile("1.0");
105
106 fEfficiencyMap.insert(make_pair(0,make_pair(0,formula)));
107 }
108
109 // import input array
110
111 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
112 fItInputArray = fInputArray->MakeIterator();
113
114 // create output array
115
116 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
117}
118
119//------------------------------------------------------------------------------
120
121void IdentificationMap::Finish()
122{
123 if(fItInputArray) delete fItInputArray;
124
125 TMisIDMap::iterator itEfficiencyMap;
126 DelphesFormula *formula;
127 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
128 {
129 formula = (itEfficiencyMap->second).second;
130 if(formula) delete formula;
131 }
132
133}
134
135//------------------------------------------------------------------------------
136
137void IdentificationMap::Process()
138{
139 Candidate *candidate;
140 Double_t pt, eta, phi;
141 TMisIDMap::iterator itEfficiencyMap;
142 pair <TMisIDMap::iterator, TMisIDMap::iterator> range;
143 DelphesFormula *formula;
144 Int_t pdgIn, pdgOut, charge;
145
146 Double_t P, Pi;
147
148 // cout<<"------------ New Event ------------"<<endl;
149
150 fItInputArray->Reset();
151 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
152 {
153 const TLorentzVector &candidatePosition = candidate->Position;
154 const TLorentzVector &candidateMomentum = candidate->Momentum;
155 eta = candidatePosition.Eta();
156 phi = candidatePosition.Phi();
157 pt = candidateMomentum.Pt();
158 pdgIn = candidate->PID;
159 charge = candidate->Charge;
160
161 // cout<<"------------ New Candidate ------------"<<endl;
162 // cout<<candidate->PID<<" "<<pt<<","<<eta<<","<<phi<<endl;
163
164 P = 1.0;
165
166 //first check that PID of this particle is specified in cfg, if not set look for PID=0
167
168 itEfficiencyMap = fEfficiencyMap.find(pdgIn);
169
170 range = fEfficiencyMap.equal_range(pdgIn);
171 if(range.first == range.second) range = fEfficiencyMap.equal_range(-pdgIn);
172 if(range.first == range.second) range = fEfficiencyMap.equal_range(0);
173
174 //loop over submap for this pid
175 for (TMisIDMap::iterator it=range.first; it!=range.second; ++it)
176 {
177
178 formula = (it->second).second;
179 pdgOut = (it->second).first;
180
181 Pi = formula->Eval(pt, eta);
182
183 // check that sum of probabilities does not exceed 1.
184 P = (P - Pi)/P;
185
[f355a2c]186 if( P < 0.0 ) continue;
[e2a76ae]187 else
188 {
189
[f355a2c]190 //randomly assign a PID to particle according to map
191 Double_t rndm = gRandom->Uniform();
192
193 if(rndm > P)
194 {
195 //change PID of particle
196 if(pdgOut != 0) candidate->PID = charge*pdgOut;
197 fOutputArray->Add(candidate);
198 break;
199 }
200 }
[e2a76ae]201
202 }
203
204 }
205
206}
207
208//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.