Fork me on GitHub

source: git/modules/PhotonConversions.cc@ 989ee8f

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 989ee8f was c1ce3fe, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

adapt DelphesCylindricalFormula and PhotonConversions to ROOT 6.04

  • Property mode set to 100644
File size: 6.9 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 PhotonConversions
21 *
22 *
23 * Converts photons into e+ e- pairs according to mass ditribution in the detector.
24 *
25 * \author M. Selvaggi - UCL, Louvain-la-Neuve
26 *
27 */
28
29#include "modules/PhotonConversions.h"
30
31#include "classes/DelphesClasses.h"
32#include "classes/DelphesFactory.h"
33#include "classes/DelphesCylindricalFormula.h"
34
35#include "ExRootAnalysis/ExRootResult.h"
36#include "ExRootAnalysis/ExRootFilter.h"
37#include "ExRootAnalysis/ExRootClassifier.h"
38
39#include "TMath.h"
40#include "TF1.h"
41#include "TString.h"
42#include "TFormula.h"
43#include "TRandom3.h"
44#include "TObjArray.h"
45#include "TDatabasePDG.h"
46#include "TLorentzVector.h"
47#include "TVector3.h"
48
49#include <algorithm>
50#include <stdexcept>
51#include <iostream>
52#include <sstream>
53
54using namespace std;
55
56//------------------------------------------------------------------------------
57
58PhotonConversions::PhotonConversions() :
59 fItInputArray(0), fConversionMap(0), fDecayXsec(0)
60{
61 fDecayXsec = new TF1;
62 fConversionMap = new DelphesCylindricalFormula;
63}
64
65//------------------------------------------------------------------------------
66
67PhotonConversions::~PhotonConversions()
68{
69}
70
71//------------------------------------------------------------------------------
72
73void PhotonConversions::Init()
74{
75 fRadius = GetDouble("Radius", 1.0);
76 fRadius2 = fRadius*fRadius;
77
78 fHalfLength = GetDouble("HalfLength", 3.0);
79
80 fEtaMax = GetDouble("EtaMax", 5.0);
81 fEtaMin = GetDouble("EtaMin", 2.0);
82
83 fStep = GetDouble("Step", 0.1); // in meters
84
85 fConversionMap->Compile(GetString("ConversionMap", "0.0"));
86
87#if ROOT_VERSION_CODE < ROOT_VERSION(6,04,00)
88 fDecayXsec->Compile("1.0 - 4.0/3.0 * x * (1.0 - x)");
89#else
90 fDecayXsec->GetFormula()->Compile("1.0 - 4.0/3.0 * x * (1.0 - x)");
91#endif
92 fDecayXsec->SetRange(0.0, 1.0);
93
94 // import array with output from filter/classifier module
95
96 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
97 fItInputArray = fInputArray->MakeIterator();
98
99 // create output arrays
100
101 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
102}
103
104//------------------------------------------------------------------------------
105
106void PhotonConversions::Finish()
107{
108 if(fItInputArray) delete fItInputArray;
109 if(fDecayXsec) delete fDecayXsec;
110 if(fConversionMap) delete fConversionMap;
111}
112
113//------------------------------------------------------------------------------
114
115void PhotonConversions::Process()
116{
117 Candidate *candidate, *ep, *em;
118 TLorentzVector candidatePosition, candidateMomentum;
119 TVector3 pos_i;
120 Double_t px, py, pz, pt, pt2, e, eta, phi;
121 Double_t x, y, z, t;
122 Double_t x_t, y_t, z_t, r_t;
123 Double_t x_i, y_i, z_i, r_i, phi_i;
124 Double_t dt, t1, t2, t3, t4;
125 Double_t tmp, discr, discr2;
126 Int_t nsteps, i;
127 Double_t rate, p_conv, x1, x2;
128 Bool_t converted;
129
130 fItInputArray->Reset();
131 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
132 {
133
134 if(candidate->PID != 22)
135 {
136 fOutputArray->Add(candidate);
137 }
138 else
139 {
140 candidatePosition = candidate->Position;
141 candidateMomentum = candidate->Momentum;
142 x = candidatePosition.X()*1.0E-3;
143 y = candidatePosition.Y()*1.0E-3;
144 z = candidatePosition.Z()*1.0E-3;
145
146 // check that particle position is inside the cylinder
147 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength) continue;
148
149 px = candidateMomentum.Px();
150 py = candidateMomentum.Py();
151 pz = candidateMomentum.Pz();
152 pt = candidateMomentum.Pt();
153 pt2 = candidateMomentum.Perp2();
154 eta = candidateMomentum.Eta();
155 phi = candidateMomentum.Phi();
156 e = candidateMomentum.E();
157
158 if(eta < fEtaMin || eta > fEtaMax) continue;
159
160 // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
161 tmp = px*y - py*x;
162 discr2 = pt2*fRadius2 - tmp*tmp;
163
164 if(discr2 < 0.0)
165 {
166 // no solutions
167 continue;
168 }
169
170 tmp = px*x + py*y;
171 discr = TMath::Sqrt(discr2);
172 t1 = (-tmp + discr)/pt2;
173 t2 = (-tmp - discr)/pt2;
174 t = (t1 < 0.0) ? t2 : t1;
175
176 z_t = z + pz*t;
177 if(TMath::Abs(z_t) > fHalfLength)
178 {
179 t3 = (+fHalfLength - z) / pz;
180 t4 = (-fHalfLength - z) / pz;
181 t = (t3 < 0.0) ? t4 : t3;
182 }
183
184 // final position
185 x_t = x + px*t;
186 y_t = y + py*t;
187 z_t = z + pz*t;
188
189 r_t = TMath::Sqrt(x_t*x_t + y_t*y_t + z_t*z_t);
190
191
192 // here starts conversion code
193 nsteps = Int_t(r_t/fStep);
194
195 x_i = x;
196 y_i = y;
197 z_i = z;
198
199 dt = t/nsteps;
200
201 converted = false;
202
203 for(i = 0; i < nsteps; ++i)
204 {
205 x_i += px*dt;
206 y_i += py*dt;
207 z_i += pz*dt;
208 pos_i.SetXYZ(x_i,y_i,z_i);
209
210 // convert photon position into cylindrical coordinates, cylindrical r,phi,z !!
211
212 r_i = TMath::Sqrt(x_i*x_i + y_i*y_i);
213 phi_i = pos_i.Phi();
214
215 // read conversion rate/meter from card
216 rate = fConversionMap->Eval(r_i, phi_i, z_i);
217
218 // convert into conversion probability
219 p_conv = 1 - TMath::Exp(-7.0/9.0*fStep*rate);
220
221 // case conversion occurs
222 if(gRandom->Uniform() < p_conv)
223 {
224 converted = true;
225
226 // generate x1 and x2, the fraction of the photon energy taken resp. by e+ and e-
227 x1 = fDecayXsec->GetRandom();
228 x2 = 1 - x1;
229
230 ep = static_cast<Candidate*>(candidate->Clone());
231 em = static_cast<Candidate*>(candidate->Clone());
232
233 ep->Position.SetXYZT(x_i*1.0E3, y_i*1.0E3, z_i*1.0E3, candidatePosition.T() + nsteps*dt*e*1.0E3);
234 em->Position.SetXYZT(x_i*1.0E3, y_i*1.0E3, z_i*1.0E3, candidatePosition.T() + nsteps*dt*e*1.0E3);
235
236 ep->Momentum.SetPtEtaPhiE(x1*pt, eta, phi, x1*e);
237 em->Momentum.SetPtEtaPhiE(x2*pt, eta, phi, x2*e);
238
239 ep->PID = -11;
240 em->PID = 11;
241
242 ep->Charge = 1.0;
243 em->Charge = -1.0;
244
245 ep->IsFromConversion = 1;
246 em->IsFromConversion = 1;
247
248 fOutputArray->Add(em);
249 fOutputArray->Add(ep);
250
251 break;
252 }
253 }
254 if(!converted) fOutputArray->Add(candidate);
255 }
256 }
257}
258
259//------------------------------------------------------------------------------
260
Note: See TracBrowser for help on using the repository browser.