[b443089] | 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 | *
|
---|
[b443089] | 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 | *
|
---|
[b443089] | 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 | *
|
---|
[b443089] | 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 |
|
---|
[d7d2da3] | 19 | /** \class ParticlePropagator
|
---|
| 20 | *
|
---|
| 21 | * Propagates charged and neutral particles
|
---|
[d41ba4a] | 22 | * from a given vertex to a cylinder defined by its radius,
|
---|
[d7d2da3] | 23 | * its half-length, centered at (0,0,0) and with its axis
|
---|
| 24 | * oriented along the z-axis.
|
---|
| 25 | *
|
---|
| 26 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 27 | *
|
---|
| 28 | */
|
---|
| 29 |
|
---|
| 30 | #include "modules/ParticlePropagator.h"
|
---|
| 31 |
|
---|
| 32 | #include "classes/DelphesClasses.h"
|
---|
| 33 | #include "classes/DelphesFactory.h"
|
---|
| 34 | #include "classes/DelphesFormula.h"
|
---|
| 35 |
|
---|
| 36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
[341014c] | 37 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 38 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
[d7d2da3] | 39 |
|
---|
| 40 | #include "TDatabasePDG.h"
|
---|
[341014c] | 41 | #include "TFormula.h"
|
---|
[d7d2da3] | 42 | #include "TLorentzVector.h"
|
---|
[341014c] | 43 | #include "TMath.h"
|
---|
| 44 | #include "TObjArray.h"
|
---|
| 45 | #include "TRandom3.h"
|
---|
| 46 | #include "TString.h"
|
---|
[d7d2da3] | 47 |
|
---|
[d41ba4a] | 48 | #include <algorithm>
|
---|
[d7d2da3] | 49 | #include <iostream>
|
---|
| 50 | #include <sstream>
|
---|
[341014c] | 51 | #include <stdexcept>
|
---|
[d7d2da3] | 52 |
|
---|
| 53 | using namespace std;
|
---|
| 54 |
|
---|
| 55 | //------------------------------------------------------------------------------
|
---|
| 56 |
|
---|
| 57 | ParticlePropagator::ParticlePropagator() :
|
---|
| 58 | fItInputArray(0)
|
---|
| 59 | {
|
---|
| 60 | }
|
---|
| 61 |
|
---|
| 62 | //------------------------------------------------------------------------------
|
---|
| 63 |
|
---|
| 64 | ParticlePropagator::~ParticlePropagator()
|
---|
| 65 | {
|
---|
| 66 | }
|
---|
| 67 |
|
---|
| 68 | //------------------------------------------------------------------------------
|
---|
| 69 |
|
---|
| 70 | void ParticlePropagator::Init()
|
---|
| 71 | {
|
---|
| 72 | fRadius = GetDouble("Radius", 1.0);
|
---|
[341014c] | 73 | fRadius2 = fRadius * fRadius;
|
---|
[d7d2da3] | 74 | fHalfLength = GetDouble("HalfLength", 3.0);
|
---|
| 75 | fBz = GetDouble("Bz", 0.0);
|
---|
| 76 | if(fRadius < 1.0E-2)
|
---|
[d41ba4a] | 77 | {
|
---|
[d7d2da3] | 78 | cout << "ERROR: magnetic field radius is too low\n";
|
---|
| 79 | return;
|
---|
| 80 | }
|
---|
| 81 | if(fHalfLength < 1.0E-2)
|
---|
| 82 | {
|
---|
| 83 | cout << "ERROR: magnetic field length is too low\n";
|
---|
| 84 | return;
|
---|
| 85 | }
|
---|
| 86 |
|
---|
[9330b7b] | 87 | fRadiusMax = GetDouble("RadiusMax", fRadius);
|
---|
| 88 | fHalfLengthMax = GetDouble("HalfLengthMax", fHalfLength);
|
---|
| 89 |
|
---|
[d7d2da3] | 90 | // import array with output from filter/classifier module
|
---|
| 91 |
|
---|
| 92 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
| 93 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 94 |
|
---|
[187fc41] | 95 | // import beamspot
|
---|
[bff2e33] | 96 | try
|
---|
| 97 | {
|
---|
| 98 | fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
|
---|
| 99 | }
|
---|
| 100 | catch(runtime_error &e)
|
---|
| 101 | {
|
---|
| 102 | fBeamSpotInputArray = 0;
|
---|
[e55f5b0] | 103 | }
|
---|
[d7d2da3] | 104 | // create output arrays
|
---|
| 105 |
|
---|
| 106 | fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
|
---|
[751cb9c] | 107 | fNeutralOutputArray = ExportArray(GetString("NeutralOutputArray", "neutralParticles"));
|
---|
[d7d2da3] | 108 | fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
|
---|
| 109 | fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
|
---|
| 110 | fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
|
---|
| 111 | }
|
---|
| 112 |
|
---|
| 113 | //------------------------------------------------------------------------------
|
---|
| 114 |
|
---|
| 115 | void ParticlePropagator::Finish()
|
---|
| 116 | {
|
---|
| 117 | if(fItInputArray) delete fItInputArray;
|
---|
| 118 | }
|
---|
| 119 |
|
---|
| 120 | //------------------------------------------------------------------------------
|
---|
| 121 |
|
---|
| 122 | void ParticlePropagator::Process()
|
---|
| 123 | {
|
---|
[60e1de6] | 124 | Candidate *candidate, *mother, *particle;
|
---|
| 125 | TLorentzVector particlePosition, particleMomentum, beamSpotPosition;
|
---|
[d7d2da3] | 126 | Double_t px, py, pz, pt, pt2, e, q;
|
---|
| 127 | Double_t x, y, z, t, r, phi;
|
---|
| 128 | Double_t x_c, y_c, r_c, phi_c, phi_0;
|
---|
| 129 | Double_t x_t, y_t, z_t, r_t;
|
---|
| 130 | Double_t t1, t2, t3, t4, t5, t6;
|
---|
| 131 | Double_t t_z, t_r, t_ra, t_rb;
|
---|
| 132 | Double_t tmp, discr, discr2;
|
---|
| 133 | Double_t delta, gammam, omega, asinrho;
|
---|
[187fc41] | 134 | Double_t rcu, rc2, xd, yd, zd;
|
---|
| 135 | Double_t l, d0, dz, p, ctgTheta, phip, etap, alpha;
|
---|
| 136 | Double_t bsx, bsy, bsz;
|
---|
[e55f5b0] | 137 |
|
---|
[d7d2da3] | 138 | const Double_t c_light = 2.99792458E8;
|
---|
[e55f5b0] | 139 |
|
---|
[341014c] | 140 | if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
|
---|
[187fc41] | 141 | beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 142 | else
|
---|
| 143 | {
|
---|
[341014c] | 144 | Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
|
---|
[187fc41] | 145 | beamSpotPosition = beamSpotCandidate.Position;
|
---|
| 146 | }
|
---|
[e55f5b0] | 147 |
|
---|
[d7d2da3] | 148 | fItInputArray->Reset();
|
---|
[341014c] | 149 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d7d2da3] | 150 | {
|
---|
[60e1de6] | 151 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 152 | {
|
---|
| 153 | particle = candidate;
|
---|
| 154 | }
|
---|
| 155 | else
|
---|
| 156 | {
|
---|
[341014c] | 157 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
[60e1de6] | 158 | }
|
---|
| 159 |
|
---|
| 160 | particlePosition = particle->Position;
|
---|
| 161 | particleMomentum = particle->Momentum;
|
---|
[341014c] | 162 | x = particlePosition.X() * 1.0E-3;
|
---|
| 163 | y = particlePosition.Y() * 1.0E-3;
|
---|
| 164 | z = particlePosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 165 |
|
---|
[341014c] | 166 | bsx = beamSpotPosition.X() * 1.0E-3;
|
---|
| 167 | bsy = beamSpotPosition.Y() * 1.0E-3;
|
---|
| 168 | bsz = beamSpotPosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 169 |
|
---|
[60e1de6] | 170 | q = particle->Charge;
|
---|
[d7d2da3] | 171 |
|
---|
| 172 | // check that particle position is inside the cylinder
|
---|
[9330b7b] | 173 | if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
|
---|
[d7d2da3] | 174 | {
|
---|
| 175 | continue;
|
---|
| 176 | }
|
---|
| 177 |
|
---|
[60e1de6] | 178 | px = particleMomentum.Px();
|
---|
| 179 | py = particleMomentum.Py();
|
---|
| 180 | pz = particleMomentum.Pz();
|
---|
| 181 | pt = particleMomentum.Pt();
|
---|
| 182 | pt2 = particleMomentum.Perp2();
|
---|
| 183 | e = particleMomentum.E();
|
---|
[d7d2da3] | 184 |
|
---|
| 185 | if(pt2 < 1.0E-9)
|
---|
| 186 | {
|
---|
| 187 | continue;
|
---|
| 188 | }
|
---|
| 189 |
|
---|
[9330b7b] | 190 | if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
|
---|
| 191 | {
|
---|
| 192 | mother = candidate;
|
---|
[341014c] | 193 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[9330b7b] | 194 |
|
---|
[60e1de6] | 195 | candidate->InitialPosition = particlePosition;
|
---|
| 196 | candidate->Position = particlePosition;
|
---|
[9330b7b] | 197 | candidate->L = 0.0;
|
---|
| 198 |
|
---|
[60e1de6] | 199 | candidate->Momentum = particleMomentum;
|
---|
[9330b7b] | 200 | candidate->AddCandidate(mother);
|
---|
| 201 |
|
---|
| 202 | fOutputArray->Add(candidate);
|
---|
| 203 | }
|
---|
| 204 | else if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
|
---|
[d7d2da3] | 205 | {
|
---|
[fcdb8bc] | 206 | // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
|
---|
[341014c] | 207 | tmp = px * y - py * x;
|
---|
| 208 | discr2 = pt2 * fRadius2 - tmp * tmp;
|
---|
[d41ba4a] | 209 |
|
---|
[b594101] | 210 | if(discr2 < 0.0)
|
---|
[d7d2da3] | 211 | {
|
---|
| 212 | // no solutions
|
---|
| 213 | continue;
|
---|
| 214 | }
|
---|
| 215 |
|
---|
[341014c] | 216 | tmp = px * x + py * y;
|
---|
[d7d2da3] | 217 | discr = TMath::Sqrt(discr2);
|
---|
[341014c] | 218 | t1 = (-tmp + discr) / pt2;
|
---|
| 219 | t2 = (-tmp - discr) / pt2;
|
---|
[b594101] | 220 | t = (t1 < 0.0) ? t2 : t1;
|
---|
[d7d2da3] | 221 |
|
---|
[341014c] | 222 | z_t = z + pz * t;
|
---|
[d7d2da3] | 223 | if(TMath::Abs(z_t) > fHalfLength)
|
---|
| 224 | {
|
---|
| 225 | t3 = (+fHalfLength - z) / pz;
|
---|
| 226 | t4 = (-fHalfLength - z) / pz;
|
---|
[b594101] | 227 | t = (t3 < 0.0) ? t4 : t3;
|
---|
[d7d2da3] | 228 | }
|
---|
| 229 |
|
---|
[341014c] | 230 | x_t = x + px * t;
|
---|
| 231 | y_t = y + py * t;
|
---|
| 232 | z_t = z + pz * t;
|
---|
[e55f5b0] | 233 |
|
---|
[341014c] | 234 | l = TMath::Sqrt((x_t - x) * (x_t - x) + (y_t - y) * (y_t - y) + (z_t - z) * (z_t - z));
|
---|
[d7d2da3] | 235 |
|
---|
| 236 | mother = candidate;
|
---|
[341014c] | 237 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d7d2da3] | 238 |
|
---|
[60e1de6] | 239 | candidate->InitialPosition = particlePosition;
|
---|
[341014c] | 240 | candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * e * 1.0E3);
|
---|
| 241 | candidate->L = l * 1.0E3;
|
---|
[e55f5b0] | 242 |
|
---|
[60e1de6] | 243 | candidate->Momentum = particleMomentum;
|
---|
[d7d2da3] | 244 | candidate->AddCandidate(mother);
|
---|
[d41ba4a] | 245 |
|
---|
[d7d2da3] | 246 | fOutputArray->Add(candidate);
|
---|
[d41ba4a] | 247 | if(TMath::Abs(q) > 1.0E-9)
|
---|
[d7d2da3] | 248 | {
|
---|
| 249 | switch(TMath::Abs(candidate->PID))
|
---|
| 250 | {
|
---|
[341014c] | 251 | case 11:
|
---|
| 252 | fElectronOutputArray->Add(candidate);
|
---|
| 253 | break;
|
---|
| 254 | case 13:
|
---|
| 255 | fMuonOutputArray->Add(candidate);
|
---|
| 256 | break;
|
---|
| 257 | default:
|
---|
| 258 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 259 | }
|
---|
| 260 | }
|
---|
[751cb9c] | 261 | else
|
---|
| 262 | {
|
---|
[341014c] | 263 | fNeutralOutputArray->Add(candidate);
|
---|
[751cb9c] | 264 | }
|
---|
[d7d2da3] | 265 | }
|
---|
| 266 | else
|
---|
| 267 | {
|
---|
| 268 |
|
---|
[b594101] | 269 | // 1. initial transverse momentum p_{T0}: Part->pt
|
---|
| 270 | // initial transverse momentum direction phi_0 = -atan(p_X0/p_Y0)
|
---|
| 271 | // relativistic gamma: gamma = E/mc^2; gammam = gamma * m
|
---|
| 272 | // gyration frequency omega = q/(gamma m) fBz
|
---|
| 273 | // helix radius r = p_{T0} / (omega gamma m)
|
---|
[d7d2da3] | 274 |
|
---|
[341014c] | 275 | gammam = e * 1.0E9 / (c_light * c_light); // gammam in [eV/c^2]
|
---|
| 276 | omega = q * fBz / (gammam); // omega is here in [89875518/s]
|
---|
| 277 | r = pt / (q * fBz) * 1.0E9 / c_light; // in [m]
|
---|
[d7d2da3] | 278 |
|
---|
[b594101] | 279 | phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
|
---|
[d7d2da3] | 280 |
|
---|
| 281 | // 2. helix axis coordinates
|
---|
[341014c] | 282 | x_c = x + r * TMath::Sin(phi_0);
|
---|
| 283 | y_c = y - r * TMath::Cos(phi_0);
|
---|
[d7d2da3] | 284 | r_c = TMath::Hypot(x_c, y_c);
|
---|
| 285 | phi_c = TMath::ATan2(y_c, x_c);
|
---|
| 286 | phi = phi_c;
|
---|
| 287 | if(x_c < 0.0) phi += TMath::Pi();
|
---|
| 288 |
|
---|
[a0431dc] | 289 | rcu = TMath::Abs(r);
|
---|
[341014c] | 290 | rc2 = r_c * r_c;
|
---|
[b594101] | 291 |
|
---|
[a0431dc] | 292 | // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
|
---|
[341014c] | 293 | xd = x_c * x_c * x_c - x_c * rcu * r_c + x_c * y_c * y_c;
|
---|
[b594101] | 294 | xd = (rc2 > 0.0) ? xd / rc2 : -999;
|
---|
[341014c] | 295 | yd = y_c * (-rcu * r_c + rc2);
|
---|
[b594101] | 296 | yd = (rc2 > 0.0) ? yd / rc2 : -999;
|
---|
[341014c] | 297 | zd = z + (TMath::Sqrt(xd * xd + yd * yd) - TMath::Sqrt(x * x + y * y)) * pz / pt;
|
---|
[e55f5b0] | 298 |
|
---|
[187fc41] | 299 | // use perigee momentum rather than original particle
|
---|
[e55f5b0] | 300 | // momentum, since the orignal particle momentum isn't known
|
---|
| 301 |
|
---|
[9330b7b] | 302 | px = TMath::Sign(1.0, r) * pt * (-y_c / r_c);
|
---|
| 303 | py = TMath::Sign(1.0, r) * pt * (x_c / r_c);
|
---|
[60e1de6] | 304 | etap = particleMomentum.Eta();
|
---|
[187fc41] | 305 | phip = TMath::ATan2(py, px);
|
---|
[e55f5b0] | 306 |
|
---|
[60e1de6] | 307 | particleMomentum.SetPtEtaPhiE(pt, etap, phip, particleMomentum.E());
|
---|
[e55f5b0] | 308 |
|
---|
[187fc41] | 309 | // calculate additional track parameters (correct for beamspot position)
|
---|
[e55f5b0] | 310 |
|
---|
[341014c] | 311 | d0 = ((x - bsx) * py - (y - bsy) * px) / pt;
|
---|
| 312 | dz = z - ((x - bsx) * px + (y - bsy) * py) / pt * (pz / pt);
|
---|
| 313 | p = particleMomentum.P();
|
---|
| 314 | ctgTheta = 1.0 / TMath::Tan(particleMomentum.Theta());
|
---|
[e55f5b0] | 315 |
|
---|
[d7d2da3] | 316 | // 3. time evaluation t = TMath::Min(t_r, t_z)
|
---|
| 317 | // t_r : time to exit from the sides
|
---|
| 318 | // t_z : time to exit from the front or the back
|
---|
| 319 | t_r = 0.0; // in [ns]
|
---|
| 320 | int sign_pz = (pz > 0.0) ? 1 : -1;
|
---|
[341014c] | 321 | if(pz == 0.0)
|
---|
| 322 | t_z = 1.0E99;
|
---|
| 323 | else
|
---|
| 324 | t_z = gammam / (pz * 1.0E9 / c_light) * (-z + fHalfLength * sign_pz);
|
---|
[d7d2da3] | 325 |
|
---|
[341014c] | 326 | if(r_c + TMath::Abs(r) < fRadius)
|
---|
[d7d2da3] | 327 | {
|
---|
| 328 | // helix does not cross the cylinder sides
|
---|
| 329 | t = t_z;
|
---|
| 330 | }
|
---|
| 331 | else
|
---|
| 332 | {
|
---|
[341014c] | 333 | asinrho = TMath::ASin((fRadius * fRadius - r_c * r_c - r * r) / (2 * TMath::Abs(r) * r_c));
|
---|
[d7d2da3] | 334 | delta = phi_0 - phi;
|
---|
[341014c] | 335 | if(delta < -TMath::Pi()) delta += 2 * TMath::Pi();
|
---|
| 336 | if(delta > TMath::Pi()) delta -= 2 * TMath::Pi();
|
---|
[d7d2da3] | 337 | t1 = (delta + asinrho) / omega;
|
---|
| 338 | t2 = (delta + TMath::Pi() - asinrho) / omega;
|
---|
| 339 | t3 = (delta + TMath::Pi() + asinrho) / omega;
|
---|
| 340 | t4 = (delta - asinrho) / omega;
|
---|
| 341 | t5 = (delta - TMath::Pi() - asinrho) / omega;
|
---|
| 342 | t6 = (delta - TMath::Pi() + asinrho) / omega;
|
---|
| 343 |
|
---|
[b594101] | 344 | if(t1 < 0.0) t1 = 1.0E99;
|
---|
| 345 | if(t2 < 0.0) t2 = 1.0E99;
|
---|
| 346 | if(t3 < 0.0) t3 = 1.0E99;
|
---|
| 347 | if(t4 < 0.0) t4 = 1.0E99;
|
---|
| 348 | if(t5 < 0.0) t5 = 1.0E99;
|
---|
| 349 | if(t6 < 0.0) t6 = 1.0E99;
|
---|
[d7d2da3] | 350 |
|
---|
| 351 | t_ra = TMath::Min(t1, TMath::Min(t2, t3));
|
---|
| 352 | t_rb = TMath::Min(t4, TMath::Min(t5, t6));
|
---|
| 353 | t_r = TMath::Min(t_ra, t_rb);
|
---|
[d41ba4a] | 354 | t = TMath::Min(t_r, t_z);
|
---|
[d7d2da3] | 355 | }
|
---|
| 356 |
|
---|
| 357 | // 4. position in terms of x(t), y(t), z(t)
|
---|
| 358 | x_t = x_c + r * TMath::Sin(omega * t - phi_0);
|
---|
| 359 | y_t = y_c + r * TMath::Cos(omega * t - phi_0);
|
---|
[341014c] | 360 | z_t = z + pz * 1.0E9 / c_light / gammam * t;
|
---|
[d7d2da3] | 361 | r_t = TMath::Hypot(x_t, y_t);
|
---|
| 362 |
|
---|
[187fc41] | 363 | // compute path length for an helix
|
---|
[e55f5b0] | 364 |
|
---|
[341014c] | 365 | alpha = pz * 1.0E9 / c_light / gammam;
|
---|
| 366 | l = t * TMath::Sqrt(alpha * alpha + r * r * omega * omega);
|
---|
[e55f5b0] | 367 |
|
---|
[d7d2da3] | 368 | if(r_t > 0.0)
|
---|
| 369 | {
|
---|
[e55f5b0] | 370 |
|
---|
[acd0621] | 371 | // store these variables before cloning
|
---|
[5b51d33] | 372 | if(particle == candidate)
|
---|
| 373 | {
|
---|
[341014c] | 374 | particle->D0 = d0 * 1.0E3;
|
---|
| 375 | particle->DZ = dz * 1.0E3;
|
---|
[5b51d33] | 376 | particle->P = p;
|
---|
| 377 | particle->PT = pt;
|
---|
| 378 | particle->CtgTheta = ctgTheta;
|
---|
| 379 | particle->Phi = phip;
|
---|
| 380 | }
|
---|
[e55f5b0] | 381 |
|
---|
[d7d2da3] | 382 | mother = candidate;
|
---|
[341014c] | 383 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d7d2da3] | 384 |
|
---|
[60e1de6] | 385 | candidate->InitialPosition = particlePosition;
|
---|
[341014c] | 386 | candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * c_light * 1.0E3);
|
---|
[d7d2da3] | 387 |
|
---|
[60e1de6] | 388 | candidate->Momentum = particleMomentum;
|
---|
[e55f5b0] | 389 |
|
---|
[341014c] | 390 | candidate->L = l * 1.0E3;
|
---|
[e55f5b0] | 391 |
|
---|
[341014c] | 392 | candidate->Xd = xd * 1.0E3;
|
---|
| 393 | candidate->Yd = yd * 1.0E3;
|
---|
| 394 | candidate->Zd = zd * 1.0E3;
|
---|
[b594101] | 395 |
|
---|
| 396 | candidate->AddCandidate(mother);
|
---|
[d7d2da3] | 397 |
|
---|
| 398 | fOutputArray->Add(candidate);
|
---|
| 399 | switch(TMath::Abs(candidate->PID))
|
---|
| 400 | {
|
---|
[341014c] | 401 | case 11:
|
---|
| 402 | fElectronOutputArray->Add(candidate);
|
---|
| 403 | break;
|
---|
| 404 | case 13:
|
---|
| 405 | fMuonOutputArray->Add(candidate);
|
---|
| 406 | break;
|
---|
| 407 | default:
|
---|
| 408 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 409 | }
|
---|
| 410 | }
|
---|
| 411 | }
|
---|
| 412 | }
|
---|
| 413 | }
|
---|
| 414 |
|
---|
| 415 | //------------------------------------------------------------------------------
|
---|