ATLAS Offline Software
Loading...
Searching...
No Matches
TrigInDetTrackFitParCnv_p2.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
7#include <cstring>
8using std::memset;
9
10
11bool TrigInDetTrackFitParCnv_p2 :: calculateSquareRoot(const TrigInDetTrackFitPar* p, TrigInDetTrackFitPar_p2* pP)
12{
13 int i{0},j{0},idx{0};
14 double a[5][5];
15 float L[5][5];
16
17 if((p->cov()==0)||(p->cov()->size()==0))
18 {
19 memset(&pP->m_cov[0],0,sizeof(pP->m_cov));
20 pP->m_cov[0]=-999.0;
21 pP->m_cov[2]=-999.0;
22 pP->m_cov[5]=-999.0;
23 pP->m_cov[9]=-999.0;
24 pP->m_cov[14]=-999.0;
25 return false;
26 }
27
28 for(i=0;i<5;i++)
29 for(j=i;j<5;j++) a[i][j]=a[j][i]=(*p->cov())[idx++];
30
32 {
33 idx=0;
34 for(i=0;i<5;i++) for(j=0;j<=i;j++) pP->m_cov[idx++]=L[i][j];
35 }
36 else
37 {
38 memset(&pP->m_cov[0],0,sizeof(pP->m_cov));
39 pP->m_cov[0] =(float)sqrt( (*p->cov())[0] > 0 ? (*p->cov())[0] : 0 );
40 pP->m_cov[2] =(float)sqrt( (*p->cov())[5] > 0 ? (*p->cov())[5] : 0 );
41 pP->m_cov[5] =(float)sqrt( (*p->cov())[9] > 0 ? (*p->cov())[9] : 0 );
42 pP->m_cov[9] =(float)sqrt( (*p->cov())[12] > 0 ? (*p->cov())[12] : 0 );
43 pP->m_cov[14]=(float)sqrt( (*p->cov())[14] > 0 ? (*p->cov())[14] : 0 );
44 }
45 return true;
46}
47
48std::unique_ptr<std::vector<double> >
49TrigInDetTrackFitParCnv_p2 :: restoreCovariance(const TrigInDetTrackFitPar_p2* pP)
50{
51
52 int i{0},j{0},k{0},idx{0};
53 float L[5][5], LT[5][5];
54 if(pP->m_cov[0]<0.0)
55 {
56 return nullptr;
57 }
58
59 memset(&L[0][0],0,sizeof(L));memset(&LT[0][0],0,sizeof(LT));
60
61 for(i=0;i<5;i++) for(j=0;j<=i;j++)
62 {
63 L[i][j]=LT[j][i]=pP->m_cov[idx++];
64 }
65
66 auto pV = std::make_unique<std::vector<double> >();
67
68 for(i=0;i<5;i++)
69
70 for(j=i;j<5;j++)
71 {
72 // Note: use extended (long double) precision here.
73 // That happens implicitly on x86 with optimization on; saying it
74 // explicitly ensures that we get the same results with and without
75 // optimization. (If this is a performance issue for platforms
76 // other than x86, one could change to double for those platforms.)
77 long double C=0.0;
78 for(k=0;k<5;k++)
79 C+=(long double)L[i][k]*LT[k][j];
80 pV->push_back(C);
81 }
82
83 return pV;
84}
85
86bool TrigInDetTrackFitParCnv_p2 :: CholeskyDecomposition(double a[5][5], float L[5][5])
87{
88
89 int i{0},j{0},k{0};
90 double sum{0};
91 float p[5];
92
93 for(i=0;i<5;i++)
94 {
95 for(j=i;j<5;j++)
96 {
97 sum=a[i][j];
98 for(k=i-1;k>=0;k--)
99 sum-=a[i][k]*a[j][k];
100 if(i==j)
101 {
102 if(sum<=0.0)
103 {
104 return false;
105 }
106 p[i]=sqrt(sum);L[i][i]=p[i];
107 }
108 else
109 {
110 a[j][i]=sum/p[i];
111 L[j][i]=a[j][i];
112 }
113 }
114 }
115 return true;
116}
117
118
119
120void TrigInDetTrackFitParCnv_p2 :: persToTrans( const TrigInDetTrackFitPar_p2 *persObj,
121 TrigInDetTrackFitPar *transObj,
122 MsgStream& log )
123{
124
125 log << MSG::DEBUG << "TrigInDetTrackFitParCnv_p2::persToTrans" << endmsg;
126
127 //restore param errors from the cov matrix
128 std::unique_ptr<std::vector<double> > cov (restoreCovariance(persObj));
129
130 double ea0 = -999;
131 double ephi0 = -999;
132 double ez0 = -999;
133 double eeta = -999;
134 double epT = -999;
135 if (cov) {
136 ea0 = sqrt((*cov)[0]);
137 ephi0 = sqrt((*cov)[5]);
138 ez0 = sqrt((*cov)[9]);
139 eeta = sqrt((*cov)[12]);
140 epT = sqrt((*cov)[14]);
141 }
142
143 *transObj = TrigInDetTrackFitPar (persObj->m_a0,
144 persObj->m_phi0,
145 persObj->m_z0,
146 persObj->m_eta,
147 persObj->m_pT,
148 ea0, ephi0, ez0, eeta, epT,
150 persObj->m_surfaceCoordinate,
151 cov.release());
152}
153
154void TrigInDetTrackFitParCnv_p2 :: transToPers( const TrigInDetTrackFitPar *transObj,
156 MsgStream& log )
157{
158
159 log << MSG::DEBUG << "TrigInDetTrackFitParCnv_p2::transToPers" << endmsg;
160
161 persObj->m_a0 = transObj->a0() ;
162 persObj->m_phi0 = transObj->phi0() ;
163 persObj->m_z0 = transObj->z0() ;
164 persObj->m_eta = transObj->eta() ;
165 persObj->m_pT = transObj->pT() ;
166
167 calculateSquareRoot(transObj,persObj);
168
169 persObj->m_surfaceType = transObj->surfaceType();
170 persObj->m_surfaceCoordinate = transObj->surfaceCoordinate();
171}
#define endmsg
static Double_t a
bool CholeskyDecomposition(double a[5][5], float L[5][5])
bool calculateSquareRoot(const TrigInDetTrackFitPar *, TrigInDetTrackFitPar_p2 *)
std::unique_ptr< std::vector< double > > restoreCovariance(const TrigInDetTrackFitPar_p2 *)
encapsulates LVL2 track parameters and covariance matrix The vector of track parameters consists of
void surfaceCoordinate(double c)
Setter: surface reference coordinate for non-perigee surfaces.
void eta(const double eta)
Setter: pseudorapidity.
void z0(const double z0)
Setter: longitudinal impact parameter.
void phi0(const double phi0)
Setter: azimuthal angle of the momentum.
void a0(const double a0)
Setter: transverse impact parameter.
void surfaceType(TrigSurfaceType s)
Setter: surface type PERIGEE=0, BARREL=1, ENDCAP=2.
void pT(const double pT)
Setter: transverse momentum.
struct color C