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