Fork me on GitHub

source: git/modules/TrackCovariance.cc@ a3261d7

Last change on this file since a3261d7 was 3051ea17, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

added track curvature to candidate class and switched back to regular covariance matrix

  • Property mode set to 100644
File size: 4.9 KB
RevLine 
[ff9fb2d9]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2020 Universite catholique de Louvain (UCLouvain), 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 TrackCovariance
20 *
21 * Smears track parameters according to appropriate covariance matrix.
22 *
[1388e73]23 * \authors P. Demin - UCLouvain, Louvain-la-Neuve
24 * M. Selvaggi - CERN
[ff9fb2d9]25 *
26 */
27
[7d790c1]28//FIXME add reference to Bedeschi-code
[1388e73]29//FIXME make sure about units of P, X
30//FIXME fix pt > 200 GeV issue and angle > 6.41
[7d790c1]31
[ff9fb2d9]32#include "modules/TrackCovariance.h"
33
34#include "classes/DelphesClasses.h"
35
36#include "TrackCovariance/SolGeom.h"
37#include "TrackCovariance/SolGridCov.h"
38#include "TrackCovariance/ObsTrk.h"
39
40#include "TLorentzVector.h"
41#include "TMath.h"
42#include "TObjArray.h"
43
[1388e73]44#include <iostream>
45#include <sstream>
46
47using namespace std;
48
[ff9fb2d9]49//------------------------------------------------------------------------------
50
51TrackCovariance::TrackCovariance() :
52 fGeometry(0), fCovariance(0), fItInputArray(0)
53{
54 fGeometry = new SolGeom();
55 fCovariance = new SolGridCov();
56}
57
58//------------------------------------------------------------------------------
59
60TrackCovariance::~TrackCovariance()
61{
62 if(fGeometry) delete fGeometry;
63 if(fCovariance) delete fCovariance;
64}
65
66//------------------------------------------------------------------------------
67
68void TrackCovariance::Init()
69{
70 fBz = GetDouble("Bz", 0.0);
71 fGeometry->Read(GetString("DetectorGeometry", ""));
72
73 fCovariance->Calc(fGeometry);
74
75 // import input array
76
77 fInputArray = ImportArray(GetString("InputArray", "TrackMerger/tracks"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // create output array
81
82 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
83}
84
85//------------------------------------------------------------------------------
86
87void TrackCovariance::Finish()
88{
89 if(fItInputArray) delete fItInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void TrackCovariance::Process()
95{
96 Candidate *candidate, *mother;
[1388e73]97 Double_t mass, p, pt, q, ct;
[3051ea17]98 Double_t dd0, ddz, dphi, dct, dp, dpt, dC;
[1388e73]99
[ff9fb2d9]100
101 fItInputArray->Reset();
102 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
103 {
[1388e73]104 const TLorentzVector &candidatePosition = candidate->InitialPosition;
[ff9fb2d9]105 const TLorentzVector &candidateMomentum = candidate->Momentum;
106
[7b518f0]107 mass = candidateMomentum.M();
[ff9fb2d9]108
109 ObsTrk track(candidatePosition.Vect(), candidateMomentum.Vect(), candidate->Charge, fBz, fCovariance);
110
[46c8df8]111 mother = candidate;
[ff9fb2d9]112 candidate = static_cast<Candidate *>(candidate->Clone());
[1388e73]113
[7b518f0]114 candidate->Momentum.SetVectM(track.GetObsP(), mass);
[1388e73]115 candidate->InitialPosition.SetXYZT(track.GetObsX().X(),track.GetObsX().Y(),track.GetObsX().Z(),candidatePosition.T());
[c18dca6]116
[3051ea17]117 // save full covariance 5x5 matrix internally (D0, phi, Curvature, dz, ctg(theta))
118 candidate->TrackCovariance = track.GetCov();
[1388e73]119
120 pt = candidate->Momentum.Pt();
121 p = candidate->Momentum.P();
122 q = track.GetObsQ();
123 ct = track.GetObsPar()[4];
124
[46c8df8]125 candidate->Xd = track.GetObsX().X();
126 candidate->Yd = track.GetObsX().Y();
127 candidate->Zd = track.GetObsX().Z();
128
[1388e73]129 candidate->D0 = track.GetObsPar()[0];
[3051ea17]130 candidate->Phi = track.GetObsPar()[1];
131 candidate->C = track.GetObsPar()[2];
[1388e73]132 candidate->DZ = track.GetObsPar()[3];
133 candidate->CtgTheta = track.GetObsPar()[4];
[3051ea17]134 candidate->P = track.GetObsP().Mag();
[1388e73]135 candidate->PT = pt;
136 candidate->Charge = q;
137
138 dd0 = TMath::Sqrt(track.GetCov()(0, 0));
139 ddz = TMath::Sqrt(track.GetCov()(3, 3));
140 dphi = TMath::Sqrt(track.GetCov()(1, 1));
141 dct = TMath::Sqrt(track.GetCov()(4, 4));
142 dpt = 2 * TMath::Sqrt( track.GetCov()(2, 2))*pt*pt / (0.2998*fBz);
143 dp = TMath::Sqrt((1.+ct*ct)*dpt*dpt + 4*pt*pt*ct*ct*dct*dct/(1.+ct*ct)/(1.+ct*ct));
[3051ea17]144 dC = TMath::Sqrt(track.GetCov()(2, 2));
[1388e73]145
146 candidate->ErrorD0 = dd0;
147 candidate->ErrorDZ = ddz;
148 candidate->ErrorP = dp;
[3051ea17]149 candidate->ErrorC = dC;
[1388e73]150 candidate->ErrorCtgTheta = dct;
151 candidate->ErrorPhi = dphi;
152 candidate->ErrorPT = dpt;
153 //candidate->TrackResolution = dpt / pt;
154 candidate->TrackResolution = dp / p;
[ff9fb2d9]155
[c18dca6]156
[ff9fb2d9]157 candidate->AddCandidate(mother);
158
159 fOutputArray->Add(candidate);
160 }
161}
162
163//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.