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