[d1e6379] | 1 | /** \class TrackSmearing
|
---|
| 2 | *
|
---|
| 3 | * Performs d0, dZ, p, Theta, Phi smearing of tracks.
|
---|
| 4 | *
|
---|
| 5 | *
|
---|
| 6 | *
|
---|
| 7 | * \author A. Hart, M. Selvaggi
|
---|
| 8 | *
|
---|
| 9 | */
|
---|
| 10 |
|
---|
| 11 | #include "modules/TrackSmearing.h"
|
---|
| 12 |
|
---|
| 13 | #include "classes/DelphesClasses.h"
|
---|
| 14 | #include "classes/DelphesFactory.h"
|
---|
| 15 | #include "classes/DelphesFormula.h"
|
---|
| 16 |
|
---|
| 17 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 18 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 19 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 20 |
|
---|
| 21 | #include "TMath.h"
|
---|
| 22 | #include "TString.h"
|
---|
| 23 | #include "TFormula.h"
|
---|
| 24 | #include "TRandom3.h"
|
---|
| 25 | #include "TObjArray.h"
|
---|
| 26 | #include "TDatabasePDG.h"
|
---|
| 27 | #include "TLorentzVector.h"
|
---|
| 28 | #include "TFile.h"
|
---|
| 29 | #include "TProfile2D.h"
|
---|
| 30 |
|
---|
| 31 | #include <algorithm>
|
---|
| 32 | #include <stdexcept>
|
---|
| 33 | #include <iostream>
|
---|
| 34 | #include <sstream>
|
---|
| 35 |
|
---|
| 36 | using namespace std;
|
---|
| 37 |
|
---|
| 38 | //------------------------------------------------------------------------------
|
---|
| 39 |
|
---|
| 40 | TrackSmearing::TrackSmearing() :
|
---|
| 41 | fD0Formula(0), fDZFormula(0), fPFormula(0), fCtgThetaFormula(0), fPhiFormula(0), fItInputArray(0)
|
---|
| 42 | {
|
---|
| 43 | fD0Formula = new DelphesFormula;
|
---|
| 44 | fDZFormula = new DelphesFormula;
|
---|
| 45 | fPFormula = new DelphesFormula;
|
---|
| 46 | fCtgThetaFormula = new DelphesFormula;
|
---|
| 47 | fPhiFormula = new DelphesFormula;
|
---|
| 48 | }
|
---|
| 49 |
|
---|
| 50 | //------------------------------------------------------------------------------
|
---|
| 51 |
|
---|
| 52 | TrackSmearing::~TrackSmearing()
|
---|
| 53 | {
|
---|
| 54 | if(fD0Formula) delete fD0Formula;
|
---|
| 55 | if(fDZFormula) delete fDZFormula;
|
---|
| 56 | if(fPFormula) delete fPFormula;
|
---|
| 57 | if(fCtgThetaFormula) delete fCtgThetaFormula;
|
---|
| 58 | if(fPhiFormula) delete fPhiFormula;
|
---|
| 59 | }
|
---|
| 60 |
|
---|
| 61 | //------------------------------------------------------------------------------
|
---|
| 62 |
|
---|
| 63 | void TrackSmearing::Init()
|
---|
| 64 | {
|
---|
| 65 | // read resolution formula
|
---|
| 66 |
|
---|
| 67 | // !!! IF WE WANT TO KEEP ROOT INPUT !!!
|
---|
| 68 | if (string (GetString("D0ResolutionFormula", "0.0")) != "0.0")
|
---|
| 69 | {
|
---|
| 70 | fD0Formula->Compile(GetString("D0ResolutionFormula", "0.0"));
|
---|
| 71 | fUseD0Formula = true;
|
---|
| 72 | }
|
---|
| 73 | else
|
---|
| 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 | }
|
---|
| 84 | else
|
---|
| 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 | }
|
---|
| 95 | else
|
---|
| 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 | }
|
---|
| 106 | else
|
---|
| 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 | }
|
---|
| 117 | else
|
---|
| 118 | {
|
---|
| 119 | fPhiResolutionFile = GetString("PhiResolutionFile", "errors.root");
|
---|
| 120 | fPhiResolutionHist = GetString("PhiResolutionHist", "phi");
|
---|
| 121 | fUsePhiFormula = false;
|
---|
| 122 | }
|
---|
| 123 |
|
---|
| 124 | fApplyToPileUp = GetBool("ApplyToPileUp", true);
|
---|
| 125 |
|
---|
| 126 | // import input array
|
---|
| 127 |
|
---|
| 128 | fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
|
---|
| 129 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 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;
|
---|
| 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;
|
---|
| 162 | TProfile2D *d0ErrorHist = NULL,
|
---|
| 163 | *dzErrorHist = NULL,
|
---|
| 164 | *pErrorHist = NULL,
|
---|
| 165 | *ctgThetaErrorHist = NULL,
|
---|
| 166 | *phiErrorHist = NULL;
|
---|
| 167 |
|
---|
[bff2e33] | 168 | //cout<<fBeamSpotInputArray->GetSize ()<<endl;
|
---|
| 169 | if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
|
---|
[d1e6379] | 170 | beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 171 | else
|
---|
| 172 | {
|
---|
| 173 | Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At (0));
|
---|
| 174 | beamSpotPosition = beamSpotCandidate.Position;
|
---|
| 175 | }
|
---|
| 176 |
|
---|
| 177 |
|
---|
| 178 | if (!fUseD0Formula)
|
---|
| 179 | {
|
---|
| 180 | TFile *fin = TFile::Open (fD0ResolutionFile.c_str ());
|
---|
| 181 | d0ErrorHist = (TProfile2D *) fin->Get (fD0ResolutionHist.c_str ());
|
---|
| 182 | d0ErrorHist->SetDirectory (0);
|
---|
| 183 | fin->Close ();
|
---|
| 184 | }
|
---|
| 185 | if (!fUseDZFormula)
|
---|
| 186 | {
|
---|
| 187 | TFile *fin = TFile::Open (fDZResolutionFile.c_str ());
|
---|
| 188 | dzErrorHist = (TProfile2D *) fin->Get (fDZResolutionHist.c_str ());
|
---|
| 189 | dzErrorHist->SetDirectory (0);
|
---|
| 190 | fin->Close ();
|
---|
| 191 | }
|
---|
| 192 | if (!fUsePFormula)
|
---|
| 193 | {
|
---|
| 194 | TFile *fin = TFile::Open (fPResolutionFile.c_str ());
|
---|
| 195 | pErrorHist = (TProfile2D *) fin->Get (fPResolutionHist.c_str ());
|
---|
| 196 | pErrorHist->SetDirectory (0);
|
---|
| 197 | fin->Close ();
|
---|
| 198 | }
|
---|
| 199 | if (!fUseCtgThetaFormula)
|
---|
| 200 | {
|
---|
| 201 | TFile *fin = TFile::Open (fCtgThetaResolutionFile.c_str ());
|
---|
| 202 | ctgThetaErrorHist = (TProfile2D *) fin->Get (fCtgThetaResolutionHist.c_str ());
|
---|
| 203 | ctgThetaErrorHist->SetDirectory (0);
|
---|
| 204 | fin->Close ();
|
---|
| 205 | }
|
---|
| 206 | if (!fUsePhiFormula)
|
---|
| 207 | {
|
---|
| 208 | TFile *fin = TFile::Open (fPhiResolutionFile.c_str ());
|
---|
| 209 | phiErrorHist = (TProfile2D *) fin->Get (fPhiResolutionHist.c_str ());
|
---|
| 210 | phiErrorHist->SetDirectory (0);
|
---|
| 211 | fin->Close ();
|
---|
| 212 | }
|
---|
| 213 |
|
---|
| 214 | fItInputArray->Reset();
|
---|
| 215 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
| 216 | {
|
---|
| 217 |
|
---|
| 218 | const TLorentzVector &momentum = candidate->Momentum;
|
---|
| 219 | const TLorentzVector &position = candidate->InitialPosition;
|
---|
| 220 |
|
---|
| 221 | pt = momentum.Pt();
|
---|
| 222 | eta = momentum.Eta();
|
---|
| 223 |
|
---|
| 224 | d0 = trueD0 = candidate->D0;
|
---|
| 225 | dz = trueDZ = candidate->DZ;
|
---|
| 226 | p = trueP = candidate->P;
|
---|
| 227 | ctgTheta = trueCtgTheta = candidate->CtgTheta;
|
---|
| 228 | phi = truePhi = candidate->Phi;
|
---|
| 229 |
|
---|
| 230 | if (fUseD0Formula)
|
---|
| 231 | d0Error = fD0Formula->Eval(pt, eta);
|
---|
| 232 | else
|
---|
| 233 | {
|
---|
| 234 | Int_t xbin, ybin;
|
---|
| 235 |
|
---|
| 236 | xbin = pt < d0ErrorHist->GetXaxis ()->GetXmax () ? d0ErrorHist->GetXaxis ()->FindBin (pt)
|
---|
| 237 | : d0ErrorHist->GetXaxis ()->GetBinCenter (d0ErrorHist->GetXaxis ()->GetNbins ());
|
---|
| 238 | ybin = d0ErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
|
---|
| 239 | d0Error = d0ErrorHist->GetBinContent (xbin, ybin);
|
---|
| 240 | if (!d0Error)
|
---|
| 241 | d0Error = -1.0;
|
---|
| 242 | }
|
---|
| 243 | if (d0Error < 0.0)
|
---|
| 244 | continue;
|
---|
| 245 |
|
---|
| 246 | if (fUseDZFormula)
|
---|
| 247 | dzError = fDZFormula->Eval(pt, eta);
|
---|
| 248 | else
|
---|
| 249 | {
|
---|
| 250 | Int_t xbin, ybin;
|
---|
| 251 |
|
---|
| 252 | xbin = pt < dzErrorHist->GetXaxis ()->GetXmax () ? dzErrorHist->GetXaxis ()->FindBin (pt)
|
---|
| 253 | : dzErrorHist->GetXaxis ()->GetBinCenter (dzErrorHist->GetXaxis ()->GetNbins ());
|
---|
| 254 | ybin = dzErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
|
---|
| 255 | dzError = dzErrorHist->GetBinContent (xbin, ybin);
|
---|
| 256 | if (!dzError)
|
---|
| 257 | dzError = -1.0;
|
---|
| 258 | }
|
---|
| 259 | if (dzError < 0.0)
|
---|
| 260 | continue;
|
---|
| 261 |
|
---|
| 262 | if (fUsePFormula)
|
---|
| 263 | pError = fPFormula->Eval(pt, eta) * p;
|
---|
| 264 | else
|
---|
| 265 | {
|
---|
| 266 | Int_t xbin, ybin;
|
---|
| 267 |
|
---|
| 268 | xbin = pt < pErrorHist->GetXaxis ()->GetXmax () ? pErrorHist->GetXaxis ()->FindBin (pt)
|
---|
| 269 | : 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;
|
---|
| 274 | }
|
---|
| 275 | if (pError < 0.0)
|
---|
| 276 | continue;
|
---|
| 277 |
|
---|
| 278 | if (fUseCtgThetaFormula)
|
---|
| 279 | ctgThetaError = fCtgThetaFormula->Eval(pt, eta);
|
---|
| 280 | else
|
---|
| 281 | {
|
---|
| 282 | Int_t xbin, ybin;
|
---|
| 283 |
|
---|
| 284 | xbin = pt < ctgThetaErrorHist->GetXaxis ()->GetXmax () ? ctgThetaErrorHist->GetXaxis ()->FindBin (pt)
|
---|
| 285 | : ctgThetaErrorHist->GetXaxis ()->GetBinCenter (ctgThetaErrorHist->GetXaxis ()->GetNbins ());
|
---|
| 286 | ybin = ctgThetaErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
|
---|
| 287 | ctgThetaError = ctgThetaErrorHist->GetBinContent (xbin, ybin);
|
---|
| 288 | if (!ctgThetaError)
|
---|
| 289 | ctgThetaError = -1.0;
|
---|
| 290 | }
|
---|
| 291 | if (ctgThetaError < 0.0)
|
---|
| 292 | continue;
|
---|
| 293 |
|
---|
| 294 | if (fUsePhiFormula)
|
---|
| 295 | phiError = fPhiFormula->Eval(pt, eta);
|
---|
| 296 | else
|
---|
| 297 | {
|
---|
| 298 | Int_t xbin, ybin;
|
---|
| 299 |
|
---|
| 300 | xbin = pt < phiErrorHist->GetXaxis ()->GetXmax () ? phiErrorHist->GetXaxis ()->FindBin (pt)
|
---|
| 301 | : phiErrorHist->GetXaxis ()->GetBinCenter (phiErrorHist->GetXaxis ()->GetNbins ());
|
---|
| 302 | ybin = phiErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
|
---|
| 303 | phiError = phiErrorHist->GetBinContent (xbin, ybin);
|
---|
| 304 | if (!phiError)
|
---|
| 305 | phiError = -1.0;
|
---|
| 306 | }
|
---|
| 307 | if (phiError < 0.0)
|
---|
| 308 | continue;
|
---|
| 309 |
|
---|
| 310 | if (fApplyToPileUp || !candidate->IsPU)
|
---|
| 311 | {
|
---|
| 312 | d0 = gRandom->Gaus(d0, d0Error);
|
---|
| 313 | dz = gRandom->Gaus(dz, dzError);
|
---|
| 314 | p = gRandom->Gaus(p, pError);
|
---|
| 315 | ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
|
---|
| 316 | phi = gRandom->Gaus(phi, phiError);
|
---|
| 317 | }
|
---|
| 318 |
|
---|
| 319 | if(p < 0.0) continue;
|
---|
| 320 | while (phi > TMath::Pi ()) phi -= TMath::TwoPi ();
|
---|
| 321 | while (phi <= -TMath::Pi ()) phi += TMath::TwoPi ();
|
---|
| 322 |
|
---|
| 323 | mother = candidate;
|
---|
| 324 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
| 325 | candidate->D0 = d0;
|
---|
| 326 | candidate->DZ = dz;
|
---|
| 327 | candidate->P = p;
|
---|
| 328 | candidate->CtgTheta = ctgTheta;
|
---|
| 329 | candidate->Phi = phi;
|
---|
| 330 |
|
---|
| 331 | theta = TMath::ACos(ctgTheta / TMath::Sqrt (1.0 + ctgTheta * ctgTheta));
|
---|
| 332 | candidate->Momentum.SetPx (p * TMath::Cos (phi) * TMath::Sin (theta));
|
---|
| 333 | candidate->Momentum.SetPy (p * TMath::Sin (phi) * TMath::Sin (theta));
|
---|
| 334 | candidate->Momentum.SetPz (p * TMath::Cos (theta));
|
---|
| 335 | candidate->Momentum.SetE (candidate->Momentum.Pt () * TMath::CosH (eta));
|
---|
| 336 | candidate->PT = candidate->Momentum.Pt ();
|
---|
| 337 |
|
---|
| 338 | x = position.X ();
|
---|
| 339 | y = position.Y ();
|
---|
| 340 | z = position.Z ();
|
---|
| 341 | t = position.T ();
|
---|
| 342 | px = candidate->Momentum.Px ();
|
---|
| 343 | py = candidate->Momentum.Py ();
|
---|
| 344 | pz = candidate->Momentum.Pz ();
|
---|
| 345 | pt = candidate->Momentum.Pt ();
|
---|
| 346 |
|
---|
[6e91f77] | 347 | // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
|
---|
[d1e6379] | 348 |
|
---|
| 349 | candidate->InitialPosition.SetX (x + ((px * y - py * x + d0 * pt) / (py - px)));
|
---|
| 350 | candidate->InitialPosition.SetY (y + ((px * y - py * x + d0 * pt) / (py - px)));
|
---|
| 351 | x = candidate->InitialPosition.X ();
|
---|
| 352 | y = candidate->InitialPosition.Y ();
|
---|
| 353 | candidate->InitialPosition.SetZ (z + ((pz * (px * (x - beamSpotPosition.X ()) + py * (y - beamSpotPosition.Y ())) + pt * pt * (dz - z)) / (pt * pt)));
|
---|
| 354 | candidate->InitialPosition.SetT (t);
|
---|
| 355 |
|
---|
| 356 | if (fApplyToPileUp || !candidate->IsPU)
|
---|
| 357 | {
|
---|
| 358 | candidate->ErrorD0 = d0Error;
|
---|
| 359 | candidate->ErrorDZ = dzError;
|
---|
| 360 | candidate->ErrorP = pError;
|
---|
| 361 | candidate->ErrorCtgTheta = ctgThetaError;
|
---|
| 362 | candidate->ErrorPhi = phiError;
|
---|
| 363 | candidate->ErrorPT = ptError (p, ctgTheta, pError, ctgThetaError);
|
---|
| 364 | candidate->TrackResolution = pError;
|
---|
| 365 | }
|
---|
| 366 |
|
---|
| 367 | candidate->AddCandidate(mother);
|
---|
| 368 | fOutputArray->Add(candidate);
|
---|
| 369 |
|
---|
| 370 | iCandidate++;
|
---|
| 371 | }
|
---|
| 372 | }
|
---|
| 373 |
|
---|
| 374 | Double_t TrackSmearing::ptError (const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
|
---|
| 375 | {
|
---|
| 376 | Double_t a, b;
|
---|
| 377 | a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
|
---|
| 378 | b = (dP * dP) / (ctgTheta * ctgTheta + 1);
|
---|
| 379 | return sqrt (a + b);
|
---|
| 380 | }
|
---|
| 381 |
|
---|
| 382 | //------------------------------------------------------------------------------
|
---|