[d1e6379] | 1 | /** \class TrackSmearing
|
---|
| 2 | *
|
---|
| 3 | * Performs d0, dZ, p, Theta, Phi smearing of tracks.
|
---|
| 4 | *
|
---|
[5658083] | 5 | * \authors A. Hart, M. Selvaggi
|
---|
[d1e6379] | 6 | *
|
---|
| 7 | */
|
---|
| 8 |
|
---|
| 9 | #include "modules/TrackSmearing.h"
|
---|
| 10 |
|
---|
| 11 | #include "classes/DelphesClasses.h"
|
---|
| 12 | #include "classes/DelphesFactory.h"
|
---|
| 13 | #include "classes/DelphesFormula.h"
|
---|
| 14 |
|
---|
| 15 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
[341014c] | 16 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 17 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
[d1e6379] | 18 |
|
---|
| 19 | #include "TDatabasePDG.h"
|
---|
| 20 | #include "TFile.h"
|
---|
[341014c] | 21 | #include "TFormula.h"
|
---|
| 22 | #include "TLorentzVector.h"
|
---|
| 23 | #include "TMath.h"
|
---|
| 24 | #include "TObjArray.h"
|
---|
[d1e6379] | 25 | #include "TProfile2D.h"
|
---|
[341014c] | 26 | #include "TRandom3.h"
|
---|
| 27 | #include "TString.h"
|
---|
[d1e6379] | 28 |
|
---|
[341014c] | 29 | #include <algorithm>
|
---|
[d1e6379] | 30 | #include <iostream>
|
---|
| 31 | #include <sstream>
|
---|
[341014c] | 32 | #include <stdexcept>
|
---|
[d1e6379] | 33 |
|
---|
| 34 | using namespace std;
|
---|
| 35 |
|
---|
| 36 | //------------------------------------------------------------------------------
|
---|
| 37 |
|
---|
| 38 | TrackSmearing::TrackSmearing() :
|
---|
| 39 | fD0Formula(0), fDZFormula(0), fPFormula(0), fCtgThetaFormula(0), fPhiFormula(0), fItInputArray(0)
|
---|
| 40 | {
|
---|
| 41 | fD0Formula = new DelphesFormula;
|
---|
| 42 | fDZFormula = new DelphesFormula;
|
---|
| 43 | fPFormula = new DelphesFormula;
|
---|
| 44 | fCtgThetaFormula = new DelphesFormula;
|
---|
| 45 | fPhiFormula = new DelphesFormula;
|
---|
| 46 | }
|
---|
| 47 |
|
---|
| 48 | //------------------------------------------------------------------------------
|
---|
| 49 |
|
---|
| 50 | TrackSmearing::~TrackSmearing()
|
---|
| 51 | {
|
---|
| 52 | if(fD0Formula) delete fD0Formula;
|
---|
| 53 | if(fDZFormula) delete fDZFormula;
|
---|
| 54 | if(fPFormula) delete fPFormula;
|
---|
| 55 | if(fCtgThetaFormula) delete fCtgThetaFormula;
|
---|
| 56 | if(fPhiFormula) delete fPhiFormula;
|
---|
| 57 | }
|
---|
| 58 |
|
---|
| 59 | //------------------------------------------------------------------------------
|
---|
| 60 |
|
---|
| 61 | void TrackSmearing::Init()
|
---|
| 62 | {
|
---|
[cde9f31] | 63 | fBz = GetDouble("Bz", 0.0);
|
---|
| 64 |
|
---|
[d1e6379] | 65 | // read resolution formula
|
---|
| 66 |
|
---|
| 67 | // !!! IF WE WANT TO KEEP ROOT INPUT !!!
|
---|
[341014c] | 68 | if(string(GetString("D0ResolutionFormula", "0.0")) != "0.0")
|
---|
| 69 | {
|
---|
| 70 | fD0Formula->Compile(GetString("D0ResolutionFormula", "0.0"));
|
---|
| 71 | fUseD0Formula = true;
|
---|
| 72 | }
|
---|
[d1e6379] | 73 | else
|
---|
[341014c] | 74 | {
|
---|
| 75 | fD0ResolutionFile = GetString("D0ResolutionFile", "errors.root");
|
---|
| 76 | fD0ResolutionHist = GetString("D0ResolutionHist", "d0");
|
---|
| 77 | fUseD0Formula = false;
|
---|
| 78 | }
|
---|
| 79 | if(string(GetString("DZResolutionFormula", "0.0")) != "0.0")
|
---|
| 80 | {
|
---|
| 81 | fDZFormula->Compile(GetString("DZResolutionFormula", "0.0"));
|
---|
| 82 | fUseDZFormula = true;
|
---|
| 83 | }
|
---|
[d1e6379] | 84 | else
|
---|
[341014c] | 85 | {
|
---|
| 86 | fDZResolutionFile = GetString("DZResolutionFile", "errors.root");
|
---|
| 87 | fDZResolutionHist = GetString("DZResolutionHist", "dz");
|
---|
| 88 | fUseDZFormula = false;
|
---|
| 89 | }
|
---|
| 90 | if(string(GetString("PResolutionFormula", "0.0")) != "0.0")
|
---|
| 91 | {
|
---|
| 92 | fPFormula->Compile(GetString("PResolutionFormula", "0.0"));
|
---|
| 93 | fUsePFormula = true;
|
---|
| 94 | }
|
---|
[d1e6379] | 95 | else
|
---|
[341014c] | 96 | {
|
---|
| 97 | fPResolutionFile = GetString("PResolutionFile", "errors.root");
|
---|
| 98 | fPResolutionHist = GetString("PResolutionHist", "p");
|
---|
| 99 | fUsePFormula = false;
|
---|
| 100 | }
|
---|
| 101 | if(string(GetString("CtgThetaResolutionFormula", "0.0")) != "0.0")
|
---|
| 102 | {
|
---|
| 103 | fCtgThetaFormula->Compile(GetString("CtgThetaResolutionFormula", "0.0"));
|
---|
| 104 | fUseCtgThetaFormula = true;
|
---|
| 105 | }
|
---|
[d1e6379] | 106 | else
|
---|
[341014c] | 107 | {
|
---|
| 108 | fCtgThetaResolutionFile = GetString("CtgThetaResolutionFile", "errors.root");
|
---|
| 109 | fCtgThetaResolutionHist = GetString("CtgThetaResolutionHist", "ctgTheta");
|
---|
| 110 | fUseCtgThetaFormula = false;
|
---|
| 111 | }
|
---|
| 112 | if(string(GetString("PhiResolutionFormula", "0.0")) != "0.0")
|
---|
| 113 | {
|
---|
| 114 | fPhiFormula->Compile(GetString("PhiResolutionFormula", "0.0"));
|
---|
| 115 | fUsePhiFormula = true;
|
---|
| 116 | }
|
---|
[d1e6379] | 117 | else
|
---|
[341014c] | 118 | {
|
---|
| 119 | fPhiResolutionFile = GetString("PhiResolutionFile", "errors.root");
|
---|
| 120 | fPhiResolutionHist = GetString("PhiResolutionHist", "phi");
|
---|
| 121 | fUsePhiFormula = false;
|
---|
| 122 | }
|
---|
[d1e6379] | 123 |
|
---|
| 124 | fApplyToPileUp = GetBool("ApplyToPileUp", true);
|
---|
| 125 |
|
---|
| 126 | // import input array
|
---|
| 127 |
|
---|
| 128 | fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
|
---|
| 129 | fItInputArray = fInputArray->MakeIterator();
|
---|
[341014c] | 130 |
|
---|
[bff2e33] | 131 | // import beamspot
|
---|
| 132 | try
|
---|
| 133 | {
|
---|
| 134 | fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
|
---|
| 135 | }
|
---|
| 136 | catch(runtime_error &e)
|
---|
| 137 | {
|
---|
| 138 | fBeamSpotInputArray = 0;
|
---|
[341014c] | 139 | }
|
---|
| 140 |
|
---|
[d1e6379] | 141 | // create output array
|
---|
| 142 |
|
---|
| 143 | fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
|
---|
| 144 | }
|
---|
| 145 |
|
---|
| 146 | //------------------------------------------------------------------------------
|
---|
| 147 |
|
---|
| 148 | void TrackSmearing::Finish()
|
---|
| 149 | {
|
---|
| 150 | if(fItInputArray) delete fItInputArray;
|
---|
| 151 | }
|
---|
| 152 |
|
---|
| 153 | //------------------------------------------------------------------------------
|
---|
| 154 |
|
---|
| 155 | void TrackSmearing::Process()
|
---|
| 156 | {
|
---|
| 157 | Int_t iCandidate = 0;
|
---|
| 158 | TLorentzVector beamSpotPosition;
|
---|
| 159 | Candidate *candidate, *mother;
|
---|
| 160 | Double_t pt, eta, d0, d0Error, trueD0, dz, dzError, trueDZ, p, pError, trueP, ctgTheta, ctgThetaError, trueCtgTheta, phi, phiError, truePhi;
|
---|
| 161 | Double_t x, y, z, t, px, py, pz, theta;
|
---|
[cde9f31] | 162 | Double_t q, r;
|
---|
| 163 | Double_t x_c, y_c, r_c, phi_0;
|
---|
| 164 | Double_t rcu, rc2, xd, yd, zd;
|
---|
| 165 | const Double_t c_light = 2.99792458E8;
|
---|
[d1e6379] | 166 | TProfile2D *d0ErrorHist = NULL,
|
---|
| 167 | *dzErrorHist = NULL,
|
---|
| 168 | *pErrorHist = NULL,
|
---|
| 169 | *ctgThetaErrorHist = NULL,
|
---|
| 170 | *phiErrorHist = NULL;
|
---|
| 171 |
|
---|
[341014c] | 172 | if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
|
---|
[d1e6379] | 173 | beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 174 | else
|
---|
| 175 | {
|
---|
[341014c] | 176 | Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
|
---|
[d1e6379] | 177 | beamSpotPosition = beamSpotCandidate.Position;
|
---|
| 178 | }
|
---|
| 179 |
|
---|
[341014c] | 180 | if(!fUseD0Formula)
|
---|
[d1e6379] | 181 | {
|
---|
[341014c] | 182 | TFile *fin = TFile::Open(fD0ResolutionFile.c_str());
|
---|
| 183 | d0ErrorHist = (TProfile2D *)fin->Get(fD0ResolutionHist.c_str());
|
---|
| 184 | d0ErrorHist->SetDirectory(0);
|
---|
| 185 | fin->Close();
|
---|
[d1e6379] | 186 | }
|
---|
[341014c] | 187 | if(!fUseDZFormula)
|
---|
[d1e6379] | 188 | {
|
---|
[341014c] | 189 | TFile *fin = TFile::Open(fDZResolutionFile.c_str());
|
---|
| 190 | dzErrorHist = (TProfile2D *)fin->Get(fDZResolutionHist.c_str());
|
---|
| 191 | dzErrorHist->SetDirectory(0);
|
---|
| 192 | fin->Close();
|
---|
[d1e6379] | 193 | }
|
---|
[341014c] | 194 | if(!fUsePFormula)
|
---|
[d1e6379] | 195 | {
|
---|
[341014c] | 196 | TFile *fin = TFile::Open(fPResolutionFile.c_str());
|
---|
| 197 | pErrorHist = (TProfile2D *)fin->Get(fPResolutionHist.c_str());
|
---|
| 198 | pErrorHist->SetDirectory(0);
|
---|
| 199 | fin->Close();
|
---|
[d1e6379] | 200 | }
|
---|
[341014c] | 201 | if(!fUseCtgThetaFormula)
|
---|
[d1e6379] | 202 | {
|
---|
[341014c] | 203 | TFile *fin = TFile::Open(fCtgThetaResolutionFile.c_str());
|
---|
| 204 | ctgThetaErrorHist = (TProfile2D *)fin->Get(fCtgThetaResolutionHist.c_str());
|
---|
| 205 | ctgThetaErrorHist->SetDirectory(0);
|
---|
| 206 | fin->Close();
|
---|
[d1e6379] | 207 | }
|
---|
[341014c] | 208 | if(!fUsePhiFormula)
|
---|
[d1e6379] | 209 | {
|
---|
[341014c] | 210 | TFile *fin = TFile::Open(fPhiResolutionFile.c_str());
|
---|
| 211 | phiErrorHist = (TProfile2D *)fin->Get(fPhiResolutionHist.c_str());
|
---|
| 212 | phiErrorHist->SetDirectory(0);
|
---|
| 213 | fin->Close();
|
---|
[d1e6379] | 214 | }
|
---|
| 215 |
|
---|
| 216 | fItInputArray->Reset();
|
---|
[341014c] | 217 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d1e6379] | 218 | {
|
---|
[341014c] | 219 |
|
---|
[d1e6379] | 220 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 221 | const TLorentzVector &position = candidate->InitialPosition;
|
---|
[341014c] | 222 |
|
---|
[d1e6379] | 223 | pt = momentum.Pt();
|
---|
| 224 | eta = momentum.Eta();
|
---|
| 225 |
|
---|
| 226 | d0 = trueD0 = candidate->D0;
|
---|
| 227 | dz = trueDZ = candidate->DZ;
|
---|
[341014c] | 228 |
|
---|
[d1e6379] | 229 | p = trueP = candidate->P;
|
---|
| 230 | ctgTheta = trueCtgTheta = candidate->CtgTheta;
|
---|
| 231 | phi = truePhi = candidate->Phi;
|
---|
| 232 |
|
---|
[341014c] | 233 | if(fUseD0Formula)
|
---|
[d1e6379] | 234 | d0Error = fD0Formula->Eval(pt, eta);
|
---|
| 235 | else
|
---|
| 236 | {
|
---|
[341014c] | 237 | Int_t xbin, ybin;
|
---|
| 238 |
|
---|
| 239 | xbin = pt < d0ErrorHist->GetXaxis()->GetXmax() ? d0ErrorHist->GetXaxis()->FindBin(pt) : d0ErrorHist->GetXaxis()->GetBinCenter(d0ErrorHist->GetXaxis()->GetNbins());
|
---|
| 240 | ybin = d0ErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
|
---|
| 241 | d0Error = d0ErrorHist->GetBinContent(xbin, ybin);
|
---|
| 242 | if(!d0Error)
|
---|
| 243 | d0Error = -1.0;
|
---|
[d1e6379] | 244 | }
|
---|
[341014c] | 245 | if(d0Error < 0.0)
|
---|
[d1e6379] | 246 | continue;
|
---|
| 247 |
|
---|
[341014c] | 248 | if(fUseDZFormula)
|
---|
[d1e6379] | 249 | dzError = fDZFormula->Eval(pt, eta);
|
---|
| 250 | else
|
---|
| 251 | {
|
---|
[341014c] | 252 | Int_t xbin, ybin;
|
---|
| 253 |
|
---|
| 254 | xbin = pt < dzErrorHist->GetXaxis()->GetXmax() ? dzErrorHist->GetXaxis()->FindBin(pt) : dzErrorHist->GetXaxis()->GetBinCenter(dzErrorHist->GetXaxis()->GetNbins());
|
---|
| 255 | ybin = dzErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
|
---|
| 256 | dzError = dzErrorHist->GetBinContent(xbin, ybin);
|
---|
| 257 | if(!dzError)
|
---|
| 258 | dzError = -1.0;
|
---|
[d1e6379] | 259 | }
|
---|
[341014c] | 260 | if(dzError < 0.0)
|
---|
[d1e6379] | 261 | continue;
|
---|
| 262 |
|
---|
[341014c] | 263 | if(fUsePFormula)
|
---|
[d1e6379] | 264 | pError = fPFormula->Eval(pt, eta) * p;
|
---|
| 265 | else
|
---|
| 266 | {
|
---|
[341014c] | 267 | Int_t xbin, ybin;
|
---|
| 268 |
|
---|
| 269 | xbin = pt < pErrorHist->GetXaxis()->GetXmax() ? pErrorHist->GetXaxis()->FindBin(pt) : pErrorHist->GetXaxis()->GetBinCenter(pErrorHist->GetXaxis()->GetNbins());
|
---|
| 270 | ybin = pErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
|
---|
| 271 | pError = pErrorHist->GetBinContent(xbin, ybin) * p;
|
---|
| 272 | if(!pError)
|
---|
| 273 | pError = -1.0;
|
---|
[d1e6379] | 274 | }
|
---|
[341014c] | 275 | if(pError < 0.0)
|
---|
[d1e6379] | 276 | continue;
|
---|
| 277 |
|
---|
[341014c] | 278 | if(fUseCtgThetaFormula)
|
---|
[d1e6379] | 279 | ctgThetaError = fCtgThetaFormula->Eval(pt, eta);
|
---|
| 280 | else
|
---|
| 281 | {
|
---|
[341014c] | 282 | Int_t xbin, ybin;
|
---|
| 283 |
|
---|
| 284 | xbin = pt < ctgThetaErrorHist->GetXaxis()->GetXmax() ? ctgThetaErrorHist->GetXaxis()->FindBin(pt) : ctgThetaErrorHist->GetXaxis()->GetBinCenter(ctgThetaErrorHist->GetXaxis()->GetNbins());
|
---|
| 285 | ybin = ctgThetaErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
|
---|
| 286 | ctgThetaError = ctgThetaErrorHist->GetBinContent(xbin, ybin);
|
---|
| 287 | if(!ctgThetaError)
|
---|
| 288 | ctgThetaError = -1.0;
|
---|
[d1e6379] | 289 | }
|
---|
[341014c] | 290 | if(ctgThetaError < 0.0)
|
---|
[d1e6379] | 291 | continue;
|
---|
| 292 |
|
---|
[341014c] | 293 | if(fUsePhiFormula)
|
---|
[d1e6379] | 294 | phiError = fPhiFormula->Eval(pt, eta);
|
---|
| 295 | else
|
---|
| 296 | {
|
---|
[341014c] | 297 | Int_t xbin, ybin;
|
---|
| 298 |
|
---|
| 299 | xbin = pt < phiErrorHist->GetXaxis()->GetXmax() ? phiErrorHist->GetXaxis()->FindBin(pt) : phiErrorHist->GetXaxis()->GetBinCenter(phiErrorHist->GetXaxis()->GetNbins());
|
---|
| 300 | ybin = phiErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
|
---|
| 301 | phiError = phiErrorHist->GetBinContent(xbin, ybin);
|
---|
| 302 | if(!phiError)
|
---|
| 303 | phiError = -1.0;
|
---|
[d1e6379] | 304 | }
|
---|
[341014c] | 305 | if(phiError < 0.0)
|
---|
[d1e6379] | 306 | continue;
|
---|
| 307 |
|
---|
[341014c] | 308 | if(fApplyToPileUp || !candidate->IsPU)
|
---|
[d1e6379] | 309 | {
|
---|
[341014c] | 310 | d0 = gRandom->Gaus(d0, d0Error);
|
---|
| 311 | dz = gRandom->Gaus(dz, dzError);
|
---|
| 312 | p = gRandom->Gaus(p, pError);
|
---|
| 313 | ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
|
---|
| 314 | phi = gRandom->Gaus(phi, phiError);
|
---|
[d1e6379] | 315 | }
|
---|
| 316 |
|
---|
| 317 | if(p < 0.0) continue;
|
---|
[341014c] | 318 | while(phi > TMath::Pi()) phi -= TMath::TwoPi();
|
---|
| 319 | while(phi <= -TMath::Pi()) phi += TMath::TwoPi();
|
---|
[d1e6379] | 320 |
|
---|
| 321 | mother = candidate;
|
---|
[341014c] | 322 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d1e6379] | 323 | candidate->D0 = d0;
|
---|
| 324 | candidate->DZ = dz;
|
---|
| 325 | candidate->P = p;
|
---|
| 326 | candidate->CtgTheta = ctgTheta;
|
---|
| 327 | candidate->Phi = phi;
|
---|
| 328 |
|
---|
[341014c] | 329 | theta = TMath::ACos(ctgTheta / TMath::Sqrt(1.0 + ctgTheta * ctgTheta));
|
---|
| 330 | candidate->Momentum.SetPx(p * TMath::Cos(phi) * TMath::Sin(theta));
|
---|
| 331 | candidate->Momentum.SetPy(p * TMath::Sin(phi) * TMath::Sin(theta));
|
---|
| 332 | candidate->Momentum.SetPz(p * TMath::Cos(theta));
|
---|
| 333 | candidate->Momentum.SetE(candidate->Momentum.Pt() * TMath::CosH(eta));
|
---|
| 334 | candidate->PT = candidate->Momentum.Pt();
|
---|
| 335 |
|
---|
| 336 | x = position.X();
|
---|
| 337 | y = position.Y();
|
---|
| 338 | z = position.Z();
|
---|
| 339 | t = position.T();
|
---|
| 340 | px = candidate->Momentum.Px();
|
---|
| 341 | py = candidate->Momentum.Py();
|
---|
| 342 | pz = candidate->Momentum.Pz();
|
---|
| 343 | pt = candidate->Momentum.Pt();
|
---|
| 344 |
|
---|
[6e91f77] | 345 | // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
|
---|
[341014c] | 346 |
|
---|
| 347 | candidate->InitialPosition.SetX(x + ((px * y - py * x + d0 * pt) / (py - px)));
|
---|
| 348 | candidate->InitialPosition.SetY(y + ((px * y - py * x + d0 * pt) / (py - px)));
|
---|
| 349 | x = candidate->InitialPosition.X();
|
---|
| 350 | y = candidate->InitialPosition.Y();
|
---|
| 351 | candidate->InitialPosition.SetZ(z + ((pz * (px * (x - beamSpotPosition.X()) + py * (y - beamSpotPosition.Y())) + pt * pt * (dz - z)) / (pt * pt)));
|
---|
| 352 | z = candidate->InitialPosition.Z();
|
---|
[cde9f31] | 353 |
|
---|
[9fccf4f] | 354 | candidate->InitialPosition.SetT(t);
|
---|
[cde9f31] | 355 |
|
---|
| 356 | // update closest approach
|
---|
| 357 | x *= 1.0E-3;
|
---|
| 358 | y *= 1.0E-3;
|
---|
| 359 | z *= 1.0E-3;
|
---|
| 360 |
|
---|
| 361 | q = candidate->Charge;
|
---|
| 362 |
|
---|
[341014c] | 363 | r = pt / (q * fBz) * 1.0E9 / c_light; // in [m]
|
---|
[cde9f31] | 364 | phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
|
---|
| 365 |
|
---|
| 366 | // 2. helix axis coordinates
|
---|
[341014c] | 367 | x_c = x + r * TMath::Sin(phi_0);
|
---|
| 368 | y_c = y - r * TMath::Cos(phi_0);
|
---|
[cde9f31] | 369 | r_c = TMath::Hypot(x_c, y_c);
|
---|
| 370 |
|
---|
| 371 | rcu = TMath::Abs(r);
|
---|
[341014c] | 372 | rc2 = r_c * r_c;
|
---|
[cde9f31] | 373 |
|
---|
| 374 | // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
|
---|
[341014c] | 375 | xd = x_c * x_c * x_c - x_c * rcu * r_c + x_c * y_c * y_c;
|
---|
[cde9f31] | 376 | xd = (rc2 > 0.0) ? xd / rc2 : -999;
|
---|
[341014c] | 377 | yd = y_c * (-rcu * r_c + rc2);
|
---|
[cde9f31] | 378 | yd = (rc2 > 0.0) ? yd / rc2 : -999;
|
---|
[341014c] | 379 | zd = z + (TMath::Sqrt(xd * xd + yd * yd) - TMath::Sqrt(x * x + y * y)) * pz / pt;
|
---|
[cde9f31] | 380 |
|
---|
[341014c] | 381 | candidate->Xd = xd * 1.0E3;
|
---|
| 382 | candidate->Yd = yd * 1.0E3;
|
---|
| 383 | candidate->Zd = zd * 1.0E3;
|
---|
[cde9f31] | 384 |
|
---|
[341014c] | 385 | if(fApplyToPileUp || !candidate->IsPU)
|
---|
[d1e6379] | 386 | {
|
---|
[341014c] | 387 | candidate->ErrorD0 = d0Error;
|
---|
| 388 | candidate->ErrorDZ = dzError;
|
---|
| 389 | candidate->ErrorP = pError;
|
---|
| 390 | candidate->ErrorCtgTheta = ctgThetaError;
|
---|
| 391 | candidate->ErrorPhi = phiError;
|
---|
| 392 | candidate->ErrorPT = ptError(p, ctgTheta, pError, ctgThetaError);
|
---|
| 393 | candidate->TrackResolution = pError / p;
|
---|
[d1e6379] | 394 | }
|
---|
[341014c] | 395 |
|
---|
[d1e6379] | 396 | candidate->AddCandidate(mother);
|
---|
| 397 | fOutputArray->Add(candidate);
|
---|
| 398 |
|
---|
| 399 | iCandidate++;
|
---|
| 400 | }
|
---|
| 401 | }
|
---|
| 402 |
|
---|
[341014c] | 403 | Double_t TrackSmearing::ptError(const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
|
---|
[d1e6379] | 404 | {
|
---|
| 405 | Double_t a, b;
|
---|
| 406 | a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
|
---|
| 407 | b = (dP * dP) / (ctgTheta * ctgTheta + 1);
|
---|
[341014c] | 408 | return sqrt(a + b);
|
---|
[d1e6379] | 409 | }
|
---|
| 410 |
|
---|
| 411 | //------------------------------------------------------------------------------
|
---|