ATLAS Offline Software
Loading...
Searching...
No Matches
VKalVrtFitSvc.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Header include
11//-------------------------------------------------
12// Other stuff
13#include "GaudiKernel/IChronoStatSvc.h"
14//
15#include <algorithm>
16#include <cmath>
17
18 namespace Trk {
19
20//__________________________________________________________________________
21//&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
22// Interface
23//
24
25 StatusCode TrkVKalVrtFitter::VKalVrtFit(const std::vector<const Perigee*>& InpPerigee,
27 TLorentzVector& Momentum,
28 long int& Charge,
30 dvect& Chi2PerTrk,
31 std::vector< std::vector<double> >& TrkAtVrt,
32 double& Chi2,
33 IVKalState& istate,
34 bool ifCovV0 /*= false*/) const
35{
36 assert(dynamic_cast<State*> (&istate)!=nullptr);
37 State& state = static_cast<State&> (istate);
38//
39//------ extract information about selected tracks
40//
41
42 int ntrk=0;
43 StatusCode sc = CvtPerigee(InpPerigee, ntrk, state);
44 if(sc.isFailure())return StatusCode::FAILURE;
45
46 int ierr = VKalVrtFit3( ntrk, Vertex, Momentum, Charge, ErrorMatrix,
47 Chi2PerTrk, TrkAtVrt,Chi2, state, ifCovV0 ) ;
48 if (ierr) return StatusCode::FAILURE;
49 return StatusCode::SUCCESS;
50}
51
52
53
54StatusCode TrkVKalVrtFitter::VKalVrtFit(const std::vector<const xAOD::TrackParticle*> & InpTrkC,
55 const std::vector<const xAOD::NeutralParticle*> & InpTrkN,
57 TLorentzVector& Momentum,
58 long int& Charge,
60 dvect& Chi2PerTrk,
61 std::vector< std::vector<double> >& TrkAtVrt,
62 double& Chi2,
63 IVKalState& istate,
64 bool ifCovV0 /*= false*/) const
65{
66 assert(dynamic_cast<State*> (&istate)!=nullptr);
67 State& state = static_cast<State&> (istate);
68
69//
70//------ extract information about selected tracks
71//
72 int ntrk=0;
73
74 // The tmpInputC will be not owning just holding plain ptr
75 // the ownership is handled via the TParamOwner
76 // and is unique_ptr so we do not leak
77 std::vector<const TrackParameters*> tmpInputC(0);
78 std::vector<std::unique_ptr<const TrackParameters>> TParamOwner(0);
79 StatusCode sc;
80 double closestHitR=1.e6; //VK needed for FirstMeasuredPointLimit if this hit itself is absent
81 if(m_firstMeasuredPoint){ //First measured point strategy
82 //------
83 if(!InpTrkC.empty()){
84 if( m_InDetExtrapolator == nullptr ){
85 if(msgLvl(MSG::WARNING))msg()<< "No InDet extrapolator given."<<
86 "Can't use FirstMeasuredPoint with xAOD::TrackParticle!!!" << endmsg;
87 return StatusCode::FAILURE;
88 }
89 std::vector<const xAOD::TrackParticle*>::const_iterator i_ntrk;
90 if(msgLvl(MSG::DEBUG))msg()<< "Start FirstMeasuredPoint handling"<<'\n';
91 unsigned int indexFMP;
92 for (i_ntrk = InpTrkC.begin(); i_ntrk < InpTrkC.end(); ++i_ntrk) {
93 if ((*i_ntrk)->indexOfParameterAtPosition(indexFMP, xAOD::FirstMeasurement)){
94 if(msgLvl(MSG::DEBUG))msg()<< "FirstMeasuredPoint on track is discovered. Use it."<<'\n';
95 // create parameters
96 TParamOwner.emplace_back(std::make_unique<CurvilinearParameters>(
97 (*i_ntrk)->curvilinearParameters(indexFMP)));
98 //For the last one we created, push also a not owning / view ptr to tmpInputC
99 tmpInputC.push_back((TParamOwner.back()).get());
100 }else{
101 if(msgLvl(MSG::DEBUG)){
102 msg()<< "FirstMeasuredPoint on track is absent."<<
103 "Try extrapolation from Perigee to FisrtMeasuredPoint radius"<<endmsg;
104 }
105
106 TParamOwner.emplace_back(m_fitPropagator->myxAODFstPntOnTrk((*i_ntrk)));
107 //For the last one we created, push also a not owning / view ptr to tmpInputC
108 tmpInputC.push_back((TParamOwner.back()).get());
109 if(tmpInputC[tmpInputC.size()-1]==nullptr){
110 //Extrapolation failure
111 if(msgLvl(MSG::WARNING)){
112 msg()<< "InDetExtrapolator can't etrapolate xAOD::TrackParticle Perigee "<<
113 "to FirstMeasuredPoint radius! Stop vertex fit!" << endmsg;
114 }
115 return StatusCode::FAILURE;
116 }
117 }
118 if( (*i_ntrk)->radiusOfFirstHit() < closestHitR ) {
119 closestHitR=(*i_ntrk)->radiusOfFirstHit();
120 }
121 }
122 sc=CvtTrackParameters(tmpInputC,ntrk,state);
123 if(sc.isFailure()){
124 return StatusCode::FAILURE;
125 }
126 }
127 }else{
128 if(!InpTrkC.empty()) {
129 sc=CvtTrackParticle(InpTrkC,ntrk,state);
130 }
131 }
132 if(sc.isFailure())return StatusCode::FAILURE;
133 if(!InpTrkN.empty()){sc=CvtNeutralParticle(InpTrkN,ntrk,state); if(sc.isFailure())return StatusCode::FAILURE;}
134 //--
135 int ierr = VKalVrtFit3( ntrk, Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, ifCovV0 ) ;
136 if (ierr) return StatusCode::FAILURE;
137 //
138 //-- Check vertex position with respect to first measured hit and refit with plane constraint if needed
139 state.m_planeCnstNDOF = 0;
141 Amg::Vector3D cnstRefPoint(0.,0.,0.);
142 if(closestHitR==1.e6){ // Not found previously
143 const xAOD::TrackParticle * trkFMP=nullptr;
144 for(auto &trka : InpTrkC){
145 double hitR=trka->radiusOfFirstHit();
146 if(closestHitR>hitR){
147 closestHitR=hitR;
148 trkFMP=trka;
149 }
150 }
151 if(closestHitR<1.e6){
152 auto perFMP=m_fitPropagator->myxAODFstPntOnTrk(trkFMP); //FMP is calculated by extrapolation to radiusOfFirstHit
153 if(perFMP) cnstRefPoint=perFMP->position();
154 }
155 }
156 Amg::Vector3D unitMom=Amg::Vector3D(Momentum.Px()/Momentum.P(),Momentum.Py()/Momentum.P(),Momentum.Pz()/Momentum.P());
157 //----------- Use as reference either hit(state.m_globalFirstHit) or its radius(closestHitR) if hit is absent
158 if(state.m_globalFirstHit)cnstRefPoint=state.m_globalFirstHit->position();
159 //------------
160 if(Vertex.perp()>closestHitR && cnstRefPoint.perp()>0.){
161 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<<"Vertex behind first measured point is detected. Constraint is applied!"<<endmsg;
162 state.m_planeCnstNDOF = 1; // Additional NDOF due to plane constraint
163 double D= unitMom.x()*(cnstRefPoint.x()-state.m_refFrameX)
164 +unitMom.y()*(cnstRefPoint.y()-state.m_refFrameY)
165 +unitMom.z()*(cnstRefPoint.z()-state.m_refFrameZ);
166 state.m_parPlaneCnst[0]=unitMom.x();
167 state.m_parPlaneCnst[1]=unitMom.y();
168 state.m_parPlaneCnst[2]=unitMom.z();
169 state.m_parPlaneCnst[3]=D;
170 state.m_cnstRadius=std::sqrt(std::pow(cnstRefPoint.x()-state.m_refFrameX,2.)+std::pow(cnstRefPoint.y()-state.m_refFrameY,2.));
171 state.m_cnstRadiusRef[0]=-state.m_refFrameX;
172 state.m_cnstRadiusRef[1]=-state.m_refFrameY;
173 std::vector<double> saveApproxV(3,0.); state.m_ApproximateVertex.swap(saveApproxV);
174 state.m_ApproximateVertex[0]=cnstRefPoint.x();
175 state.m_ApproximateVertex[1]=cnstRefPoint.y();
176 state.m_ApproximateVertex[2]=cnstRefPoint.z();
177 ierr = VKalVrtFit3( ntrk, Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, ifCovV0 );
178 state.m_vkalFitControl.setUsePlaneCnst(0.,0.,0.,0.);
179 if (ierr) { // refit without plane cnst
180 ierr = VKalVrtFit3(ntrk,Vertex,Momentum,Charge,ErrorMatrix,Chi2PerTrk,TrkAtVrt,Chi2, state, ifCovV0); // if fit with it failed
181 state.m_planeCnstNDOF = 0;
182 }
183 state.m_ApproximateVertex.swap(saveApproxV);
184 }
185 }
186 //--
187 if (ierr) return StatusCode::FAILURE;
188 return StatusCode::SUCCESS;
189}
190
191
192StatusCode TrkVKalVrtFitter::VKalVrtFit(const std::vector<const TrackParameters*> & InpTrkC,
193 const std::vector<const NeutralParameters*> & InpTrkN,
195 TLorentzVector& Momentum,
196 long int& Charge,
198 dvect& Chi2PerTrk,
199 std::vector< std::vector<double> >& TrkAtVrt,
200 double& Chi2,
201 IVKalState& istate,
202 bool ifCovV0 /*= false*/) const
203{
204 assert(dynamic_cast<State*> (&istate)!=nullptr);
205 State& state = static_cast<State&> (istate);
206
207//
208//------ extract information about selected tracks
209//
210 int ntrk=0;
211 StatusCode sc;
212 if(!InpTrkC.empty()){
213 sc=CvtTrackParameters(InpTrkC,ntrk,state);
214 if(sc.isFailure())return StatusCode::FAILURE;
215 }
216 if(!InpTrkN.empty()){
217 sc=CvtNeutralParameters(InpTrkN,ntrk,state);
218 if(sc.isFailure())return StatusCode::FAILURE;
219 }
220
221 if(state.m_ApproximateVertex.empty() && state.m_globalFirstHit){ //Initial guess if absent
222 state.m_ApproximateVertex.reserve(3);
223 state.m_ApproximateVertex.push_back(state.m_globalFirstHit->position().x());
224 state.m_ApproximateVertex.push_back(state.m_globalFirstHit->position().y());
225 state.m_ApproximateVertex.push_back(state.m_globalFirstHit->position().z());
226 }
227 int ierr = VKalVrtFit3( ntrk, Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt,Chi2, state, ifCovV0 ) ;
228 if (ierr) return StatusCode::FAILURE;
229//
230//-- Check vertex position with respect to first measured hit and refit with plane constraint if needed
231 state.m_planeCnstNDOF = 0;
232 if(state.m_globalFirstHit && m_firstMeasuredPointLimit && !ierr){
233 if(Vertex.perp()>state.m_globalFirstHit->position().perp()){
234 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<<"Vertex behind first measured point is detected. Constraint is applied!"<<endmsg;
235 state.m_planeCnstNDOF = 1; // Additional NDOF due to plane constraint
236 double pp[3]={Momentum.Px()/Momentum.Rho(),Momentum.Py()/Momentum.Rho(),Momentum.Pz()/Momentum.Rho()};
237 double D= pp[0]*(state.m_globalFirstHit->position().x()-state.m_refFrameX)
238 +pp[1]*(state.m_globalFirstHit->position().y()-state.m_refFrameY)
239 +pp[2]*(state.m_globalFirstHit->position().z()-state.m_refFrameZ);
240 state.m_vkalFitControl.setUsePlaneCnst( pp[0], pp[1], pp[2], D);
241 std::vector<double> saveApproxV(3,0.); state.m_ApproximateVertex.swap(saveApproxV);
242 state.m_ApproximateVertex[0]=state.m_globalFirstHit->position().x();
243 state.m_ApproximateVertex[1]=state.m_globalFirstHit->position().y();
244 state.m_ApproximateVertex[2]=state.m_globalFirstHit->position().z();
245 ierr = VKalVrtFit3( ntrk, Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt,Chi2, state, ifCovV0 ) ;
246 state.m_vkalFitControl.setUsePlaneCnst(0.,0.,0.,0.);
247 if (ierr) { // refit without plane cnst
248 ierr = VKalVrtFit3(ntrk,Vertex,Momentum,Charge,ErrorMatrix,Chi2PerTrk,TrkAtVrt,Chi2, state, ifCovV0 ) ; // if fit with it failed
249 state.m_planeCnstNDOF = 0;
250 }
251 state.m_ApproximateVertex.swap(saveApproxV);
252 }
253 }
254 if (ierr) return StatusCode::FAILURE;
255 return StatusCode::SUCCESS;
256}
257
258
259
260
261
262
263//--------------------------------------------------------------------------------------------------
264// Main code
265//
268 TLorentzVector& Momentum,
269 long int& Charge,
271 dvect& Chi2PerTrk,
272 std::vector< std::vector<double> >& TrkAtVrt,
273 double& Chi2,
274 State& state,
275 bool ifCovV0) const
276{
277//
278//------ Variables and arrays needed for fitting kernel
279//
280 int ierr,i;
281 double xyz0[3],covf[21],chi2f=-10.;
282 double ptot[4]={0.};
283 double xyzfit[3]={0.};
284//
285//--- Set field value at (0.,0.,0.) - some safety
286//
287 double Bx,By,Bz;
288 state.m_fitField.getMagFld(-state.m_refFrameX,-state.m_refFrameY,-state.m_refFrameZ,Bx,By,Bz);
289//
290//------ Fit option setting
291//
292 VKalVrtConfigureFitterCore(ntrk, state);
293//
294//------ Fit itself
295//
296 state.m_FitStatus=0;
297 state.m_vkalFitControl.renewFullCovariance(nullptr); //
300 if(state.m_ApproximateVertex.size()==3 && fabs(state.m_ApproximateVertex[2])<m_IDsizeZ &&
302 {
303 xyz0[0]=(double)state.m_ApproximateVertex[0] - state.m_refFrameX;
304 xyz0[1]=(double)state.m_ApproximateVertex[1] - state.m_refFrameY;
305 xyz0[2]=(double)state.m_ApproximateVertex[2] - state.m_refFrameZ;
306 } else {
307 xyz0[0]=xyz0[1]=xyz0[2]=0.;
308 }
309 double par0[NTrMaxVFit][3]; //used only for fit preparation
310 Trk::cfpest( ntrk, xyz0, state.m_ich, state.m_apar, par0);
311
312 Chi2PerTrk.resize (ntrk);
313 ierr=Trk::CFit( &state.m_vkalFitControl, ifCovV0, ntrk, state.m_ich, xyz0, par0, state.m_apar, state.m_awgt,
314 xyzfit, state.m_parfs, ptot, covf, chi2f,
315 Chi2PerTrk.data());
316
317 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG) << "VKalVrt fit status="<<ierr<<" Chi2="<<chi2f<<endmsg;
318
319 Chi2 = 100000000.;
320 if(ierr){
321 return ierr;
322 }
323 if(ptot[0]*ptot[0]+ptot[1]*ptot[1] == 0.) return -5; // Bad (divergent) fit
324//
325// Postfit operation. Creation of array for different error calculations and full error matrix copy
326//
327 state.m_FitStatus=ntrk;
328 if(ifCovV0 && state.m_vkalFitControl.getFullCovariance()){ //If full fit error matrix is returned by VKalVrtCORE
329 int SymCovMtxSize=(3*ntrk+3)*(3*ntrk+4)/2;
330 state.m_ErrMtx.assign (state.m_vkalFitControl.getFullCovariance(),
331 state.m_vkalFitControl.getFullCovariance()+SymCovMtxSize);
333 ErrorMatrix.clear(); ErrorMatrix.reserve(21); ErrorMatrix.assign(covf,covf+21);
334 } else {
335 ErrorMatrix.clear(); ErrorMatrix.reserve(6); ErrorMatrix.assign(covf,covf+6);
336 }
337//---------------------------------------------------------------------------
338 Momentum.SetPxPyPzE( ptot[0], ptot[1], ptot[2], ptot[3] );
339 Chi2 = (double) chi2f;
340
341 Vertex[0]= xyzfit[0] + state.m_refFrameX;
342 Vertex[1]= xyzfit[1] + state.m_refFrameY;
343 Vertex[2]= xyzfit[2] + state.m_refFrameZ;
344
345 double sizeR = state.m_allowUltraDisplaced ? m_MSsizeR : m_IDsizeR;
346 double sizeZ = state.m_allowUltraDisplaced ? m_MSsizeZ : m_IDsizeZ;
347 if (Vertex.perp() > sizeR || std::abs(Vertex.z()) > sizeZ) return -5; // Solution outside acceptable volume due to divergence
348
349 state.m_save_xyzfit[0]=xyzfit[0]; // saving of vertex position
350 state.m_save_xyzfit[1]=xyzfit[1]; // for full error matrix
351 state.m_save_xyzfit[2]=xyzfit[2];
352//
353// ------ Magnetic field in fitted vertex
354//
355 double fx,fy,fz;
356 state.m_fitField.getMagFld(xyzfit[0] ,xyzfit[1] ,xyzfit[2] ,fx,fy,fz);
357
358 Charge=0; for(i=0; i<ntrk; i++){Charge+=state.m_ich[i];};
359 Charge=-Charge; //VK 30.11.2009 Change sign acoording to ATLAS
360
361
362 TrkAtVrt.clear(); TrkAtVrt.reserve(ntrk);
363 for(i=0; i<ntrk; i++){
364 std::vector<double> TrkPar(3);
365 double effectiveBMAG=state.m_fitField.getEffField(fx, fy, fz, state.m_parfs[i][1], state.m_parfs[i][0]);
366 if(std::abs(effectiveBMAG) < 0.01) effectiveBMAG=0.01; //safety
367 VKalToTrkTrack(effectiveBMAG,(double)state.m_parfs[i][0],(double)state.m_parfs[i][1],(double)state.m_parfs[i][2],
368 TrkPar[0],TrkPar[1],TrkPar[2]);
369 TrkPar[2] = -TrkPar[2]; // Change of sign needed
370 TrkAtVrt.push_back( TrkPar );
371 }
372 return 0;
373 }
374
375
376
377// Converts Vertex, Mom, CovVrtMom in GLOBAL SYSTEM into perigee
378// Works correctly only in ID (solenoidal field)!
379//
380
382 const TLorentzVector& Momentum,
383 const dvect& CovVrtMom,
384 const long int& Charge,
385 dvect& Perigee,
386 dvect& CovPerigee,
387 IVKalState& istate) const
388 {
389 assert(dynamic_cast<State*> (&istate)!=nullptr);
390 State& state = static_cast<State&> (istate);
391 int i,j,ij;
392 double Vrt[3],PMom[4],Cov0[21],Per[5],CovPer[15];
393
394 for(i=0; i<3; i++) Vrt[i]=Vertex[i];
395 for(i=0; i<3; i++) PMom[i]=Momentum[i];
396 for(ij=i=0; i<6; i++){
397 for(j=0; j<=i; j++){
398 Cov0[ij]=CovVrtMom[ij];
399 ij++;
400 }
401 }
402 state.m_refFrameX=state.m_refFrameY=state.m_refFrameZ=0.; //VK Work in ATLAS ref frame ONLY!!!
403 long int vkCharge=-Charge; //VK 30.11.2009 Change sign according to ATLAS
404//
405// ------ Magnetic field in vertex (solenoidal field, ID only)
406//
407 double fx,fy,BMAG_CUR;
408 state.m_fitField.getMagFld(Vrt[0], Vrt[1], Vrt[2] ,fx,fy,BMAG_CUR);
409 if(fabs(BMAG_CUR) < 0.01) BMAG_CUR=0.01; // Safety
410
411 Trk::xyztrp( vkCharge, Vrt, PMom, Cov0, BMAG_CUR, Per, CovPer );
412
413 Perigee.clear();
414 CovPerigee.clear();
415
416 for(i=0; i<5; i++) Perigee.push_back((double)Per[i]);
417 for(i=0; i<15; i++) CovPerigee.push_back((double)CovPer[i]);
418
419 return StatusCode::SUCCESS;
420 }
421
422
423 void TrkVKalVrtFitter::VKalToTrkTrack( double effectiveBMAG, double vp1, double vp2, double vp3,
424 double& tp1, double& tp2, double& tp3) const
425//tp - ATLAS parameters, vp - VKalVrt parameters//
426 { tp1= vp2; //phi angle
427 tp2= vp1; //theta angle
428 tp3= vp3 * std::sin( vp1 ) /(m_CNVMAG*effectiveBMAG);
429 constexpr double pi = M_PI;
430 // -pi < phi < pi range
431 while ( tp1 > pi) tp1 -= 2.*pi;
432 while ( tp1 <-pi) tp1 += 2.*pi;
433 // 0 < Theta < pi range
434 while ( tp2 > pi) tp2 -= 2.*pi;
435 while ( tp2 <-pi) tp2 += 2.*pi;
436 if ( tp2 < 0.) {
437 tp2 = fabs(tp2); tp1 += pi;
438 while ( tp1 > pi) tp1 -= 2.*pi;
439 }
440
441 }
442
443
444/* Returns a complete error matrix after fit
445 useMom=0 def (V,Perigee1,Perigee2....PerigeeNTrk)
446 useMom=1 (V,PxPyPz1,PxPyPz2....PxPyPzNTrk)
447*/
448
449
450 StatusCode
451 TrkVKalVrtFitter::VKalGetFullCov( long int NTrk, dvect& CovVrtTrk,
452 IVKalState& istate,
453 bool useMom) const
454 {
455 assert(dynamic_cast<State*> (&istate)!=nullptr);
456 State& state = static_cast<State&> (istate);
457 if(!state.m_FitStatus) return StatusCode::FAILURE;
458 if(NTrk<1) return StatusCode::FAILURE;
459 if(NTrk>NTrMaxVFit) return StatusCode::FAILURE;
460 if(state.m_ErrMtx.empty()) return StatusCode::FAILURE; //Now error matrix is taken from CORE in VKalVrtFit3.
461//
462// ------ Magnetic field access
463//
464 double fx,fy,fz;
465 state.m_fitField.getMagFld(state.m_save_xyzfit[0],state.m_save_xyzfit[1],state.m_save_xyzfit[2],fx,fy,fz);
466//
467// ------ Base code
468//
469 int i,j,ik,jk,ip,iTrk;
470 int DIM=3*NTrk+3; //Current size of full covariance matrix
471 std::vector<std::vector<double> > Deriv (DIM);
472 for (std::vector<double>& v : Deriv) v.resize (DIM);
473 std::vector<double> CovMtxOld(DIM*DIM);
474
475
476 CovVrtTrk.resize(DIM*(DIM+1)/2);
477
478 ip=0;
479 for( i=0; i<DIM;i++) {
480 for( j=0; j<=i; j++) {
481 CovMtxOld[i*DIM+j]=CovMtxOld[j*DIM+i]=state.m_ErrMtx[ip++];
482 }
483 }
484
485 //delete [] ErrMtx;
486
487 for(i=0;i<DIM;i++){ for(j=0;j<DIM;j++) {Deriv[i][j]=0.;}}
488 Deriv[0][0]= 1.;
489 Deriv[1][1]= 1.;
490 Deriv[2][2]= 1.;
491
492 int iSt=0;
493 double Theta,invR,Phi;
494 for( iTrk=0; iTrk<NTrk; iTrk++){
495 Theta=state.m_parfs[iTrk][0];
496 Phi =state.m_parfs[iTrk][1];
497 invR =state.m_parfs[iTrk][2];
498 double effectiveBMAG=state.m_fitField.getEffField(fx, fy, fz, Phi, Theta);
499 if(std::abs(effectiveBMAG) < 0.01) effectiveBMAG = 0.01;
500
501 /*-----------*/
502 /* dNew/dOld */
503 iSt = 3 + iTrk*3;
504 if( !useMom ){
505 Deriv[iSt ][iSt+1] = 1; // Phi <-> Theta
506 Deriv[iSt+1][iSt ] = 1; // Phi <-> Theta
507 Deriv[iSt+2][iSt ] = -(cos(Theta)/(m_CNVMAG*effectiveBMAG)) * invR ; // d1/p / dTheta
508 Deriv[iSt+2][iSt+2] = -(sin(Theta)/(m_CNVMAG*effectiveBMAG)) ; // d1/p / d1/R
509 }else{
510 double pt=std::abs(m_CNVMAG*effectiveBMAG/invR);
511 double px=pt*cos(Phi);
512 double py=pt*sin(Phi);
513 double pz=pt/tan(Theta);
514 Deriv[iSt ][iSt ]= 0; //dPx/dTheta
515 Deriv[iSt ][iSt+1]= -py; //dPx/dPhi
516 Deriv[iSt ][iSt+2]= -px/invR; //dPx/dinvR
517
518 Deriv[iSt+1][iSt ]= 0; //dPy/dTheta
519 Deriv[iSt+1][iSt+1]= px; //dPy/dPhi
520 Deriv[iSt+1][iSt+2]= -py/invR; //dPy/dinvR
521
522 Deriv[iSt+2][iSt ]= -pt/sin(Theta)/sin(Theta); //dPz/dTheta
523 Deriv[iSt+2][iSt+1]= 0; //dPz/dPhi
524 Deriv[iSt+2][iSt+2]= -pz/invR; //dPz/dinvR
525 }
526 }
527//---------- Only upper half if filled and saved
528 int ipnt=0;
529 double tmp, tmpTmp;
530 for(i=0;i<DIM;i++){
531 for(j=0;j<=i;j++){
532 tmp=0.;
533 for(ik=0;ik<DIM;ik++){
534 if(Deriv[i][ik] == 0.) continue;
535 tmpTmp=0;
536 for(jk=DIM-1;jk>=0;jk--){
537 if(Deriv[j][jk] == 0.) continue;
538 tmpTmp += CovMtxOld[ik*DIM+jk]*Deriv[j][jk];
539 }
540 tmp += Deriv[i][ik]*tmpTmp;
541 }
542 CovVrtTrk[ipnt++]=tmp;
543 }}
544
545 return StatusCode::SUCCESS;
546
547 }
548
549
550
551
552
553 StatusCode TrkVKalVrtFitter::VKalGetMassError( double& dM, double& MassError,
554 const IVKalState& istate) const
555 {
556 assert(dynamic_cast<const State*> (&istate)!=nullptr);
557 const State& state = static_cast<const State&> (istate);
558 if(!state.m_FitStatus) return StatusCode::FAILURE;
559 dM = state.m_vkalFitControl.getVertexMass();
560 MassError = state.m_vkalFitControl.getVrtMassError();
561 return StatusCode::SUCCESS;
562 }
563
564
566 const IVKalState& istate) const
567 {
568 assert(dynamic_cast<const State*> (&istate)!=nullptr);
569 const State& state = static_cast<const State&> (istate);
570 if(!state.m_FitStatus) return StatusCode::FAILURE; // no fit made
571 trkWeights.clear();
572
573 int NTRK=state.m_FitStatus;
574
575 for (int i=0; i<NTRK; i++) trkWeights.push_back(state.m_vkalFitControl.vk_forcft.robres[i]);
576
577 return StatusCode::SUCCESS;
578 }
579
580
582 {
583 if(!state.m_FitStatus) return 0;
584 int NDOF=2*state.m_FitStatus-3;
585 if(state.m_usePointingCnst) { NDOF+=2; }
586 else if(state.m_useZPointingCnst) { NDOF+=1; }
587 if( state.m_usePassNear || state.m_usePassWithTrkErr ) { NDOF+= 2; }
588
589 if( state.m_massForConstraint>0. ) { NDOF+=1; }
590 if( !state.m_partMassCnst.empty() ) { NDOF+= state.m_partMassCnst.size(); }
591 if( state.m_useAprioriVertex ) { NDOF+= 3; }
592 if( state.m_usePhiCnst ) { NDOF+=1; }
593 if( state.m_useThetaCnst ) { NDOF+=1; }
594 return NDOF;
595 }
596
597}
#define M_PI
#define endmsg
static Double_t sc
#define pi
const Amg::Vector3D & position() const
Access method for the position.
std::vector< double > m_partMassCnst
std::vector< double > m_ApproximateVertex
double m_apar[NTrMaxVFit][5]
double m_awgt[NTrMaxVFit][15]
double m_parfs[NTrMaxVFit][3]
const TrackParameters * m_globalFirstHit
std::vector< double > m_ErrMtx
StatusCode CvtTrackParameters(const std::vector< const TrackParameters * > &InpTrk, int &ntrk, State &state) const
void VKalVrtConfigureFitterCore(int NTRK, State &state) const
const IExtrapolator * m_InDetExtrapolator
Pointer to Extrapolator AlgTool.
Gaudi::Property< double > m_MSsizeZ
virtual StatusCode VKalVrtCvtTool(const Amg::Vector3D &Vertex, const TLorentzVector &Momentum, const dvect &CovVrtMom, const long int &Charge, dvect &Perigee, dvect &CovPerigee, IVKalState &istate) const override final
StatusCode CvtPerigee(const std::vector< const Perigee * > &list, int &ntrk, State &state) const
virtual StatusCode VKalGetMassError(double &Mass, double &MassError, const IVKalState &istate) const override final
Gaudi::Property< double > m_IDsizeZ
int VKalVrtFit3(int ntrk, Amg::Vector3D &Vertex, TLorentzVector &Momentum, long int &Charge, dvect &ErrorMatrix, dvect &Chi2PerTrk, std::vector< std::vector< double > > &TrkAtVrt, double &Chi2, State &state, bool ifCovV0) const
virtual StatusCode VKalGetFullCov(long int, dvect &CovMtx, IVKalState &istate, bool=false) const override final
virtual StatusCode VKalVrtFit(const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::NeutralParticle * > &, Amg::Vector3D &Vertex, TLorentzVector &Momentum, long int &Charge, dvect &ErrorMatrix, dvect &Chi2PerTrk, std::vector< std::vector< double > > &TrkAtVrt, double &Chi2, IVKalState &istate, bool ifCovV0=false) const override final
Gaudi::Property< bool > m_firstMeasuredPoint
virtual StatusCode VKalGetTrkWeights(dvect &Weights, const IVKalState &istate) const override final
static int VKalGetNDOF(const State &state)
Gaudi::Property< bool > m_firstMeasuredRadiusLimit
StatusCode CvtTrackParticle(std::span< const xAOD::TrackParticle *const > list, int &ntrk, State &state) const
VKalExtPropagator * m_fitPropagator
Gaudi::Property< double > m_IDsizeR
Gaudi::Property< bool > m_firstMeasuredPointLimit
StatusCode CvtNeutralParticle(const std::vector< const xAOD::NeutralParticle * > &list, int &ntrk, State &state) const
Gaudi::Property< double > m_MSsizeR
void VKalToTrkTrack(double curBMAG, double vp1, double vp2, double vp3, double &tp1, double &tp2, double &tp3) const
StatusCode CvtNeutralParameters(const std::vector< const NeutralParameters * > &InpTrk, int &ntrk, State &state) const
virtual void getMagFld(const double, const double, const double, double &, double &, double &) override
void renewFullCovariance(double *)
double getVertexMass() const
void setVertexMass(double mass)
double getVrtMassError() const
void setUsePlaneCnst(double a, double b, double c, double d)
void setVrtMassError(double error)
const double * getFullCovariance() const
This class is a simplest representation of a vertex candidate.
double getEffField(double bx, double by, double bz, double phi, double theta)
Definition VKalVrtBMag.h:41
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the ATLAS eigen extensions are properly loaded.
void xyztrp(const long int ich, double *vrt0, double *pv0, double *covi, double BMAG, double *paro, double *errt)
Definition XYZtrp.cxx:16
int CFit(VKalVrtControl *FitCONTROL, int ifCovV0, int NTRK, long int *ich, double xyz0[3], double(*par0)[3], double(*inp_Trk5)[5], double(*inp_CovTrk5)[15], double xyzfit[3], double(*parfs)[3], double ptot[4], double covf[21], double &chi2, double *chi2tr)
Definition CFit.cxx:435
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
void cfpest(int ntrk, double *xyz, long int *ich, double(*parst)[5], double(*parf)[3])
Definition cfPEst.cxx:10
@ x
Definition ParamDefs.h:55
@ z
global position (cartesian)
Definition ParamDefs.h:57
@ pz
global momentum (cartesian)
Definition ParamDefs.h:61
@ v
Definition ParamDefs.h:78
@ y
Definition ParamDefs.h:56
@ px
Definition ParamDefs.h:59
@ py
Definition ParamDefs.h:60
std::vector< double > dvect
TrackParticle_v1 TrackParticle
Reference the current persistent version:
@ FirstMeasurement
Parameter defined at the position of the 1st measurement.
double robres[vkalNTrkM]
Definition ForCFT.h:51
MsgStream & msg
Definition testRead.cxx:32