ATLAS Offline Software
Loading...
Searching...
No Matches
CFit.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
21#include <array>
22#include <cmath>
23#include <iostream>
24
25/*----------------------------------------------------------------------------*/
26/* CONSTRAINT VERTEX FIT: */
27/* INPUT : */
28/* IFLAG : what to do */
29/* 0 - simple vertex fit */
30/* 1 - fit with mass constraint */
31/* 2 - fit with constraint on passing through vertex VRT*/
32/* 3 - fit with Z-constraint */
33/* 4 - 1+2 */
34/* 5 - 1+3 */
35/* 6 - fit with primary vertex constraint */
36/* 7 - fit with summary vector close to given vertex */
37/* 8 - 1+7 */
38/* 9 - like 7 but summary track error is not used */
39/* 10 - like 8 but summary track error is not used */
40/* 11 - Phi angle between track is 0 */
41/* 12 - Phi and Theta angles between track are 0 */
42/* 13 - 11+7 */
43/* 14 - 12+7 */
44/* IFCOVV0 : =0 only vertex error matrix(6) is returned */
45/* >0 full V0 x,p error matrix(21) is returned */
46/* (calculations are CPU consuming in case of */
47/* constraints 1:5 with big track number) */
48/* NTRK : number of tracks */
49/* ICH(JTK) : charge of tracks */
50/* XYZ0(3) : starting point */
51/* PAR0(3,JTK) : tracks parameters in XYZ */
52/* APAR(5,JTK) : measured tracks parameters */
53/* AWGT(15,JTK): tracks parameter errors */
54/* OUTPUT: XYZFIT(3) : fitted sec. vertex */
55/* PARFS(3,JTK): th,phi,1/r of tracks at vertex */
56/* PTOT(3) : total momentum in vertex */
57/* COVF(21) : covariance matrix (x,y,z,ptot(3)) */
58/* CHI2 : total chi-squared */
59/* CHI2TR(JTK) : chi-squared per track */
60/* IERR : error flag */
61/* -1 bad fit chi2>1000000 ( no fitted output) */
62/* -2 fit is ready but without covariance matrix */
63/* -3 no fit at all */
64/* -4 noninvertible vert.error matrix for IFLAG=6 */
65/* -5 Divergence of iterations */
66/* Author: V.Kostyukhin ( 1996 ) */
67/*----------------------------------------------------------------------------*/
68
69
70namespace{
71using namespace Trk;
72// Functions only used internally in this module for the implementation
73// of the external modules go to the anonymous namespace.
74// They have internal linkage.
75
76// New data model for cascade fit
77void fillVertex(VKVertex *vk, int NTRK, const long int *ich,
78 double xyz0[3], double (*par0)[3], double (*inp_Trk5)[5],
79 double (*inp_CovTrk5)[15]) {
80
81
82 ForCFT &vrtForCFT = vk->vk_fitterControl->vk_forcft;
83 double xyz[3] = {0., 0., 0.}; // Perigee point
84 //
85 // VK 13.01.2009 New strategy with ref.point left at initial position
86 vk->setRefV(xyz); // ref point
87 vk->setRefIterV(xyz0); // iteration ref. point
88 vk->tmpArr.resize(NTRK);
89 vk->TrackList.resize(NTRK);
90 for (int tk = 0; tk < NTRK; tk++) {
91 long int TrkID = tk;
92 vk->TrackList[tk] = std::make_unique<VKTrack>(
93 TrkID, &inp_Trk5[tk][0], &inp_CovTrk5[tk][0], vk, vrtForCFT.wm[tk]);
94 //printT(&inp_Trk5[tk][0], &inp_CovTrk5[tk][0] , "Input track:");
95 //std::cout<<(*vk->TrackList[tk]);
96 vk->tmpArr[tk]= std::make_unique< TWRK > ();
97 vk->TrackList[tk]->Charge = ich[tk]; // Charge coinsides with sign of curvature
98 }
99
100 vk->FVC.Charge = 0; for (int tk=0; tk<NTRK; ++tk)vk->FVC.Charge += ich[tk]; // summary charge
101 cfdcopy(vrtForCFT.vrt, vk->FVC.vrt , 3);
102 cfdcopy(vrtForCFT.covvrt, vk->FVC.covvrt , 6);
103
104 /* --------------------------------------------------------*/
105 /* Initial value */
106 /* --------------------------------------------------------*/
107 for(int tk=0; tk<NTRK; tk++){
108 VKTrack *trk = vk->TrackList[tk].get();
109 trk->iniP[0]=trk->cnstP[0]=trk->fitP[0]=par0[tk][0]; //initial guess
110 trk->iniP[1]=trk->cnstP[1]=trk->fitP[1]=par0[tk][1];
111 trk->iniP[2]=trk->cnstP[2]=trk->fitP[2]=par0[tk][2];
112 }
113 vk->setIniV(xyz); vk->setCnstV(xyz);
114}
115
116bool checkPosition(VKVertex *vk, double vertex[3]) {
117 bool insideGoodVolume = true;
118 if (vk->vk_fitterControl && vk->vk_fitterControl->vk_objProp) {
119 insideGoodVolume = vk->vk_fitterControl->vk_objProp->checkTarget(
120 vertex, *vk->vk_fitterControl->vk_istate);
121 } else {
122 insideGoodVolume = Trk::vkalPropagator::checkTarget(vertex);
123 }
124 return insideGoodVolume;
125}
126
127void protectCurvatureSign(
128 double ini, double cur,
129 double *Wgt) { // VK 15.12.2009 Protection against sign change
130 double x = cur / ini;
131 if (x >= 0.7)
132 return;
133 if (x <= 0.)
134 return;
135 Wgt[14] *= (0.7 * 0.7) / (x * x);
136}
137
138int fitVertex(VKVertex * vk)
139{
140 int i, jerr, tk, it=0;
141
142 double dparst[6];
143 double chi2min, chi21s=11., chi22s=10., vShift;
144 double aermd[6],tmpd[6]={0.}; // temporary arrays
145 double tmpPer[5],tmpCov[15], tmpWgt[15];
146 double VrtMomCov[21],PartMom[4];
147 double cnstRemnants=0., iniCnstRem=1.e-12;
148 double newVrtXYZ[3];
149
150 const double ConstraintAccuracy=1.e-6;
151
152 //
153 // New datastructure
154
155 long int IERR = 0;
156 /* Type of the constraint for the fit */
157 int NTRK = vk->TrackList.size();
158 ForCFT& vrtForCFT=vk->vk_fitterControl->vk_forcft;
159 vk->passNearVertex=false; //explicit setting - not to use
160 vk->passWithTrkCov=false; //explicit setting - not to use
161 if ( vrtForCFT.usePassNear ) vk->passNearVertex=true;
162 if ( vrtForCFT.usePassNear==2 ) vk->passWithTrkCov=true;
163 double zeroV[3]={0.,0.,0.};
164 VKTrack * trk=nullptr;
165
166 chi2min = 1e15;
167
168 if( vrtForCFT.nmcnst && vrtForCFT.useMassCnst ) { //mass constraints are present
169 std::vector<int> index;
170 for( int ic=0; ic<vkalMaxNMassCnst; ic++){
171 if(vrtForCFT.wmfit[ic]>0){ // new mass constraint
172 index.clear();
173 for(tk=0; tk<NTRK; tk++){ if( vrtForCFT.indtrkmc[ic][tk] )index.push_back(tk); }
174 vk->ConstraintList.emplace_back(std::make_unique<VKMassConstraint>( NTRK, vrtForCFT.wmfit[ic], std::move(index), vk));
175 }
176 }
177 }
178 if( vrtForCFT.usePointingCnst==1 ){ //3Dpointing
179 vk->ConstraintList.emplace_back(std::make_unique<VKPointConstraint>( NTRK, vrtForCFT.vrt, vk, false));
180 }
181 if( vrtForCFT.usePointingCnst==2 ){ //Z pointing
182 vk->ConstraintList.emplace_back(std::make_unique<VKPointConstraint>( NTRK, vrtForCFT.vrt, vk, true));
183 }
184 if ( vrtForCFT.useAprioriVrt ) {
185 cfdcopy(vrtForCFT.covvrt, tmpd, 6);
186 IERR=cfdinv(tmpd, aermd, -3); if (IERR) { IERR = -4; return IERR; }
187 vk->useApriorVertex = 1;
188 cfdcopy(vrtForCFT.vrt, vk->apriorV, 3);
189 cfdcopy( aermd, vk->apriorVWGT,6);
190 }
191 if ( vrtForCFT.usePhiCnst ) vk->ConstraintList.emplace_back(std::make_unique<VKPhiConstraint>( NTRK, vk));
192 if ( vrtForCFT.useThetaCnst )vk->ConstraintList.emplace_back(std::make_unique<VKThetaConstraint>( NTRK, vk));
193 if ( vrtForCFT.usePlaneCnst && vrtForCFT.Ap+vrtForCFT.Bp+vrtForCFT.Cp != 0.){
194 vk->ConstraintList.emplace_back(std::make_unique<VKPlaneConstraint>( NTRK, vrtForCFT.Ap, vrtForCFT.Bp, vrtForCFT.Cp, vrtForCFT.Dp, vk));
195 }
196 if ( vrtForCFT.useRadiusCnst && vrtForCFT.RC != 0.){
197 vk->ConstraintList.emplace_back(std::make_unique<VKRadiusConstraint>( NTRK, vrtForCFT.RC, vrtForCFT.radiusRefP, vk));
198 }
199
200 //-----Debug printout
201 // for(auto & cnst : vk->ConstraintList) {
202 // VKMassConstraint *ctmp=dynamic_cast<VKMassConstraint*>( cnst.get() ); if(ctmp) std::cout<<(*ctmp)<<'\n';
203 // VKPointConstraint *ptmp=dynamic_cast<VKPointConstraint*>( cnst.get() ); if(ptmp) std::cout<<(*ptmp)<<'\n';
204 // }
205 //
206 // Needed for track close to vertex constraint
207 for(i=0; i<2*6; i++)vk->FVC.cvder[i]=0.;
208 for(i=0; i<21; i++)vk->FVC.dcovf[i]=0.;
209 vk->FVC.dcovf[0]=10.;vk->FVC.dcovf[2] =10.;vk->FVC.dcovf[5] =10.;
210 vk->FVC.dcovf[9]=10.;vk->FVC.dcovf[14]=10.;vk->FVC.dcovf[20]=10.;
211
212 /* --------------------------------------------------------*/
213 /* Initial value */
214 /* --------------------------------------------------------*/
215 cfdcopy(vk->refIterV,newVrtXYZ,3); // copy initial data to start values
216 vk->setIniV(zeroV);vk->setCnstV(zeroV); // initial guess. 0 of course.
217 /* ------------------------------------------------------------*/
218 /* MAIN FITTING CYCLE */
219 /* ------------------------------------------------------------*/
220 bool extrapolationDone=false;
221 bool forcedExtrapolation=false; //To force explicit extrapolation
222 std::vector< Vect3DF > savedExtrapVertices(vrtForCFT.IterationNumber);
223 double cnstRemnantsPrev=1.e20, Chi2Prev=0.; int countBadStep=0; // For consecutive bad steps
224
225 for (it = 1; it <= vrtForCFT.IterationNumber; ++it) {
226 vShift = cfddist3D( newVrtXYZ, vk->refIterV ); // Extrapolation distance
227 savedExtrapVertices[it-1].Set(newVrtXYZ); // Save new position. Index [it] starts from 1 !!!
228 /* ---------------------------------------------------------------- */
229 /* Propagate parameters and errors to the point newVrtXYZ */
230 /* Also set up localBMAG at newVrtXYZ if nonuniform field is used */
231 /* ---------------------------------------------------------------- */
232 extrapolationDone=false;
233 bool insideGoodVolume=checkPosition(vk,newVrtXYZ);
234 if( vShift>vkalShiftToTrigExtrapolation || it==1 || forcedExtrapolation || !insideGoodVolume ){ //Propagation is needed
235 //--- Check whether propagation step must be truncated
236 if( forcedExtrapolation || vk->truncatedStep || !insideGoodVolume){
237 insideGoodVolume=false;
238 while( !insideGoodVolume ){ //check extrapolation and limit step if needed
239 newVrtXYZ[0]=(2.*savedExtrapVertices[it-1].X + savedExtrapVertices[it-2].X)/3.; //Use 2/3 of the initial distance
240 newVrtXYZ[1]=(2.*savedExtrapVertices[it-1].Y + savedExtrapVertices[it-2].Y)/3.; //for extrapolation
241 newVrtXYZ[2]=(2.*savedExtrapVertices[it-1].Z + savedExtrapVertices[it-2].Z)/3.;
242 savedExtrapVertices[it-1].Set(newVrtXYZ);
243 insideGoodVolume=checkPosition(vk,newVrtXYZ);
244 if(!insideGoodVolume && savedExtrapVertices[it-1].Dist3D(savedExtrapVertices[it-2])<5.) { return -11; }
245 } }
246 //------------------------------
247 extrapolationDone=true;
248 forcedExtrapolation=false;
249 double targV[3]={newVrtXYZ[0],newVrtXYZ[1],newVrtXYZ[2]}; //Temporary to avoid overwriting
250 for (tk = 0; tk < NTRK; ++tk) {
251 //std::cout<<__func__<<" propagate trk="<<tk<<" X,Y,Z="<<targV[0]<<","<<targV[1]<<","<<targV[2]<<'\n';
252 Trk::vkalPropagator::Propagate(vk->TrackList[tk].get(), vk->refV, targV, tmpPer, tmpCov, (vk->vk_fitterControl).get());
253 if(std::isnan(tmpCov[14])) {return -7;} // want to make sure it isn't a nan
254 if (std::abs(tmpCov[14])<1.e-20 // Zero Q/p covariance. Stop fit and return failure
255 && !vk->vk_fitterControl->m_allowUltraDisplaced) { //vertices approaching exterior of MS toroid can have wacky Q/P covariances -- that's ok
256 return -7;
257 }
258 cfTrkCovarCorr(tmpCov);
259 double eig5=cfSmallEigenvalue(tmpCov,5 );
260 if(eig5>0 && eig5<1.e-15 ){
261 tmpCov[0]+=1.e-15; tmpCov[2]+=1.e-15; tmpCov[5]+=1.e-15; tmpCov[9]+=1.e-15; tmpCov[14]+=1.e-15;
262 }else if(tmpCov[0]>1.e9 || eig5<0.) { //Bad propagation with material. Try without it.
263 Trk::vkalPropagator::Propagate(-999, vk->TrackList[tk]->Charge, vk->TrackList[tk]->refPerig,vk->TrackList[tk]->refCovar,
264 vk->refV, targV, tmpPer, tmpCov, (vk->vk_fitterControl).get());
265
266 if(cfSmallEigenvalue(tmpCov,5 )<1.e-15){ // Final protection. Zero non-diagonal terms
267 tmpCov[1]=0.;tmpCov[3]=0.;tmpCov[6]=0.;tmpCov[10]=0.; // and make positive diagonal ones
268 tmpCov[4]=0.;tmpCov[7]=0.;tmpCov[11]=0.;
269 tmpCov[8]=0.;tmpCov[12]=0.;
270 tmpCov[13]=0.;
271 tmpCov[0]=std::abs(tmpCov[0]); tmpCov[2]=std::abs(tmpCov[2]);tmpCov[5]=std::abs(tmpCov[5]);
272 tmpCov[9]=std::abs(tmpCov[9]); tmpCov[14]=std::abs(tmpCov[14]);
273 }
274 }
275 if(tmpCov[0]>1.e9){ IERR=-7; return IERR; } //extremely big A0 error !!!
276 jerr=cfInv5(tmpCov, tmpWgt);
277 if(jerr)jerr=cfdinv(tmpCov, tmpWgt,5);
278 if(jerr){ IERR=-7; return IERR; } //Nothing helps. Break fit.
279 vk->TrackList[tk]->setCurrent(tmpPer,tmpWgt);
280 }
281//
282/* Magnetic field in fitting point setting */
283 vk->setRefIterV(targV);
285 vk->setIniV(zeroV);vk->setCnstV(zeroV); // initial guess with respect to refIterV[]. 0 after extrap.
286 }else{
287 vk->setIniV( vk->fitV ); vk->setCnstV( vk->fitV ); // use previous fitV[] as start position
288 for( tk=0; tk<NTRK; tk++) vk->TrackList[tk]->restoreCurrentWgt(); // restore matrix changed by ROBTEST!!!
289 }
290/*-----------------------------------------------------------------------*/
291/* All parameters are set in new vertex position. We are ready for fit */
292/* Track WGT matrix is completely clean now */
293/* - apply protection against charge sign change */
294/*-----------------------------------------------------------------------*/
295 for (tk = 0; tk < NTRK; ++tk){
296 trk = vk->TrackList[tk].get(); protectCurvatureSign( trk->refPerig[4], trk->fitP[2] , trk->WgtM);
297 }
298/*-------------------------------- Now the fit itself -----------------*/
299 if (vrtForCFT.irob != 0) {robtest(vk, 0, it);} // ROBUSTIFICATION new data structure
300 if (vrtForCFT.irob != 0) {robtest(vk, 1, it);} // ROBUSTIFICATION new data structure
301 for( tk=0; tk<NTRK; tk++){
302 trk = vk->TrackList[tk].get();
303 trk->iniP[0]=trk->cnstP[0]=trk->fitP[0]; //use fitted track parameters as initial guess
304 trk->iniP[1]=trk->cnstP[1]=trk->fitP[1];
305 trk->iniP[2]=trk->cnstP[2]=trk->fitP[2];
306 if(extrapolationDone)trk->iniP[2]=trk->cnstP[2]=trk->Perig[4]; //Use exact value here to take into account change of mag.field
307 }
308 applyConstraints( vk );
309 if(it==1) { bool bsave = vk->passNearVertex; vk->passNearVertex=false; IERR = vtcfit( vk ); vk->passNearVertex=bsave;}
310 else IERR = vtcfit( vk ); // In first iteration data for PassNearVertex constraint don't exist!!!
311 if (IERR) return IERR;
312 vShift = cfddist3D( vk->fitV, vk->iniV );
313//----
314 if(extrapolationDone && chi21s*10.<vk->Chi2 && it>2 && vk->Chi2>NTRK*5.){ //Extrapolation produced huge degradation
315 double ddstep=savedExtrapVertices[it-1].Dist3D(savedExtrapVertices[it-2]);
316//std::cout<<" Huge degradation due to extrapolation. Limit step! it="<<it<<" step="<<ddstep<<'\n';
317 if( ddstep > 5.*vkalShiftToTrigExtrapolation) {
318 forcedExtrapolation=true;
319 continue;
320 }
321 }
322 chi21s = vk->Chi2;
323 chi22s = chi21s * 1.01 + 10.; //for safety
324 if ( vShift < 10.*vkalShiftToTrigExtrapolation) { // REASONABLE DISPLACEMENT - RECALCULATE
325/* ROBUSTIFICATION */
326 if (vrtForCFT.irob != 0) {robtest(vk, 1, it+1);} // ROBUSTIFICATION new data structure
327//Reset mag.field
328 for( i=0; i<3; i++) dparst[i]=vk->refIterV[i]+vk->fitV[i]; // fitted vertex at global frame
329 vrtForCFT.localbmag=Trk::vkalMagFld::getMagFld(dparst,(vk->vk_fitterControl).get());
330 if ( vk->passNearVertex && it>1 ) { //No necessary information at first iteration
331 jerr = afterFit(vk, vk->ader, vk->FVC.dcv, PartMom, VrtMomCov, (vk->vk_fitterControl).get());
332 if(jerr!=0) return -17; // Non-invertable matrix for combined track
333 cfdcopy( PartMom, &dparst[3], 3); //vertex part of it is filled above
334 cfdcopy(VrtMomCov,vk->FVC.dcovf,21); //Used in chi2 caclulation later...
335 cfdcopy( PartMom, vk->fitMom, 3); //save Momentum
336 cfdcopy(VrtMomCov, vk->fitCovXYZMom, 21); //save XYZMom covariance
337 vpderiv(vk->passWithTrkCov, vk->FVC.Charge, dparst, VrtMomCov, vk->FVC.vrt, vk->FVC.covvrt,
338 vk->FVC.cvder, vk->FVC.ywgt, vk->FVC.rv0, (vk->vk_fitterControl).get());
339 }
340
341 for( tk=0; tk<NTRK; tk++){
342 trk = vk->TrackList[tk].get();
343 trk->iniP[0]=trk->cnstP[0]=trk->fitP[0]; //use fitted track parameters as initial guess
344 trk->iniP[1]=trk->cnstP[1]=trk->fitP[1];
345 trk->iniP[2]=trk->cnstP[2]=trk->fitP[2];
346 }
347 vk->setIniV( vk->fitV ); vk->setCnstV( vk->fitV ); // use previous fitV[] as start position
348 applyConstraints( vk);
349 if(it==1) { bool bsave = vk->passNearVertex; vk->passNearVertex=false; IERR = vtcfit( vk ); vk->passNearVertex=bsave;}
350 else IERR = vtcfit( vk ); // In first iteration data for PassNearVertex constraint don't exist!!!
351 if ( IERR ) return IERR;
352 vShift = cfddist3D( vk->fitV, vk->iniV );
353 chi22s = vk->Chi2;
354 }
355//
356//
357 if (chi22s > 1e8) { return -1; } // TOO HIGH CHI2 - BAD FIT
358 if (chi22s != chi22s) { return -14; } // Chi2 == nan - BAD FIT
359 for( i=0; i<3; i++) newVrtXYZ[i] = vk->refIterV[i]+vk->fitV[i]; // fitted vertex in global frame
360 //std::cout.precision(11);
361 //std::cout<<"NNFIT Iter="<<it<<" Chi2ss="<< chi21s <<", "<<chi22s<<", "<<vShift<<'\n';
362 //std::cout<<"NNVertex="<<newVrtXYZ[0]<<", "<<newVrtXYZ[1]<<", "<<newVrtXYZ[2]<<'\n';
363 //std::cout<<"-----------------------------------------------"<<'\n';
364/* Test of convergence */
365 double chi2df = std::abs(chi21s - chi22s);
366 /*---------------------Normal convergence--------------------*/
367 double PrecLimit = std::min(chi22s*1.e-4, vrtForCFT.IterationPrecision);
368//std::cout<<"Convergence="<< chi2df <<"<"<<PrecLimit<<" cnst="<<cnstRemnants<<"<"<<ConstraintAccuracy<<'\n';
369 if( ( vk->vk_fitterControl->m_frozenVersionForBTagging && it>1 )
370 || (!vk->vk_fitterControl->m_frozenVersionForBTagging && (it>3||vrtForCFT.irob==0) ) ){
371 if((chi2df < PrecLimit) && (vShift < 0.001) && (cnstRemnants<ConstraintAccuracy)){
372 double dstFromExtrapPnt=sqrt(vk->fitV[0]*vk->fitV[0] + vk->fitV[1]*vk->fitV[1]+ vk->fitV[2]*vk->fitV[2]);
373 if( dstFromExtrapPnt>vkalShiftToTrigExtrapolation/2. && it < vrtForCFT.IterationNumber-15){
374 forcedExtrapolation=true;
375 continue; // Make another extrapolation exactly to found vertex position
376 }
377 break;
378 }
379 }
380 chi2min = std::min(chi2min,chi22s);
381 if ((chi2min*100. < chi22s) && (chi22s>std::max( (2*NTRK-3)*10., 100.)) && (it>5)){
382 //std::cout<<" DIVERGENCE="<<chi22s<<" Ratio="<<chi22s/chi2min<<'\n';
383 IERR = -5; return IERR; /* Divergence. Stop fit */
384 }
385
386// Track near vertex constraint recalculation for next fit
387 if ( vk->passNearVertex ) {
388 if(it==1) jerr = afterFit(vk, nullptr, vk->FVC.dcv, PartMom, VrtMomCov, (vk->vk_fitterControl).get());
389 else jerr = afterFit(vk, vk->ader, vk->FVC.dcv, PartMom, VrtMomCov, (vk->vk_fitterControl).get());
390 for( i=0; i<3; i++) dparst[i] = vk->refIterV[i]+vk->fitV[i]; // fitted vertex at global frame
391 cfdcopy( PartMom, &dparst[3], 3);
392 cfdcopy(VrtMomCov,vk->FVC.dcovf,21); //Used in chi2 caclulation later...
393 cfdcopy( PartMom, vk->fitMom, 3); //save Momentum
394 cfdcopy(VrtMomCov, vk->fitCovXYZMom, 21); //save XYZMom covariance
395 vpderiv(vk->passWithTrkCov, vk->FVC.Charge, dparst, VrtMomCov,
396 vk->FVC.vrt, vk->FVC.covvrt,
397 vk->FVC.cvder, vk->FVC.ywgt, vk->FVC.rv0, (vk->vk_fitterControl).get());
398 }
399//
400// Check constraints status
401//
402 cnstRemnants=0. ;
403 for(int ii=0; ii<(int)vk->ConstraintList.size();ii++){
404 for(int ic=0; ic<(int)vk->ConstraintList[ii]->NCDim; ic++){
405 cnstRemnants += std::abs( vk->ConstraintList[ii]->aa[ic] ); } }
406 if(it==1){ if(cnstRemnants>0.){iniCnstRem=cnstRemnants;}else{iniCnstRem=1.e-12;} }
407 if(it==10){
408 if(cnstRemnants > 0.1 && cnstRemnants/iniCnstRem > 0.1 ){
409 IERR = -6; return IERR; /* Constraints are not resolved. Stop fit */
410 }
411 }
412//
413// Check consecutive bad steps in fitting
414//
415 if(it>5 && ((cnstRemnants>cnstRemnantsPrev)||cnstRemnants==0.) && chi22s>(Chi2Prev*1.0001)){
416 if(++countBadStep>5) break; // Apparently further improvement is impossible
417 }else{countBadStep=0;}
418 Chi2Prev=chi22s; cnstRemnantsPrev=cnstRemnants;
419 }
420/*--------------------------------------*/
421/* End of cycling */
422/* Now we have a secondary vertex xyzs */
423/*--------------------------------------*/
424 if( cnstRemnants > ConstraintAccuracy ){
425 IERR = -6; return IERR; /* Constraints are not resolved. Bad fit */
426 }
427 return 0;
428
429}
430//
431} //end of anonymous
432
433
434namespace Trk {
435int CFit(VKalVrtControl *FitCONTROL, int ifCovV0, int NTRK,
436 long int *ich, double xyz0[3], double (*par0)[3],
437 double (*inp_Trk5)[5], double (*inp_CovTrk5)[15],
438 double xyzfit[3], double (*parfs)[3], double ptot[4],
439 double covf[21], double &chi2, double *chi2tr)
440{
441
442 std::fill_n(ptot,4,0.);
443 std::fill_n(covf,21,0.);
444 cfdcopy(xyz0,xyzfit,3);
445 FitCONTROL->renewFullCovariance(nullptr);
446
447//
448// Create vertex
449//
450 std::unique_ptr<VKVertex> MainVRT(new VKVertex(*FitCONTROL));
451 fillVertex(MainVRT.get(), NTRK, ich, xyz0, par0, inp_Trk5, inp_CovTrk5);
452//
453// Fit vertex
454//
455 int IERR = fitVertex( MainVRT.get());
456 if(IERR) return IERR;
457//
458// if successful - save results
459//
460 for(int i=0; i<3; i++) xyzfit[i] = MainVRT->refIterV[i]+MainVRT->fitV[i]; // fitted vertex in global frame
461 MainVRT->vk_fitterControl->vk_forcft.localbmag = Trk::vkalMagFld::getMagFld(xyzfit,(MainVRT->vk_fitterControl).get());
462 chi2 = 0.;
463 for (int tk = 0; tk < NTRK; tk++) {
464 VKTrack * trk = MainVRT->TrackList[tk].get();
465 chi2 += trk->Chi2;
466 chi2tr[tk] = trk->Chi2;
467 cfdcopy( trk->fitP, &parfs[tk][0], 3);
468 std::array<double, 4> pp = getFitParticleMom( trk, MainVRT.get() );
469 ptot[0]+=pp[0]; ptot[1]+=pp[1]; ptot[2]+=pp[2]; ptot[3]+=pp[3];
470 }
471 cfdcopy(MainVRT->fitVcov, covf , 6); //fitted vertex covariance
472 double vM2=(ptot[3]-ptot[2])*(ptot[3]+ptot[2]) - ptot[1]*ptot[1] - ptot[0]*ptot[0];
473 FitCONTROL->setVertexMass(std::sqrt(vM2>0. ? vM2 : 0. ) );
474 FitCONTROL->setVrtMassError(0.);
475//
476// If required - get full covariance matrix
477//
478 if(ifCovV0){
479 double dptot[5], VrtMomCov[21];
480 IERR = afterFit(MainVRT.get(), MainVRT->ader, MainVRT->FVC.dcv, dptot, VrtMomCov, (MainVRT->vk_fitterControl).get());
481 if (IERR) return -2; /* NONINVERTIBLE COV.MATRIX */
482 if(MainVRT->existFullCov){
483 FitCONTROL->renewFullCovariance(new double[ (3*NTRK+3)*(3*NTRK+3+1)/2 ]);
484 int out=0;
485 for(int ti=0; ti<3+3*NTRK; ti++){
486 for(int tj=0; tj<=ti; tj++){
487 FitCONTROL->getFullCovariance()[out] = ARR2D_FS(MainVRT->ader, 3*vkalNTrkM+3, ti, tj);
488 out++;
489 } }
490 int activeTrk[vkalNTrkM]={0}; std::fill_n(activeTrk,NTRK,1);
491 double vrtMass=0., vrtMassError=0.;
492 cfmasserr(MainVRT.get(), activeTrk, &vrtMass, &vrtMassError);
493 FitCONTROL->setVertexMass(vrtMass);
494 FitCONTROL->setVrtMassError(vrtMassError);
495 }
496 cfdcopy( VrtMomCov, covf, 21); //XvYvZvPxPyPz covariance
497 cfdcopy( MainVRT->vk_fitterControl->vk_forcft.robres, FitCONTROL->vk_forcft.robres, NTRK); //Track weights from robust fit
498
499 }
500
501 return 0;
502}
503
504} //End of Trk
#define vkalMaxNMassCnst
Definition CommonPars.h:27
#define vkalShiftToTrigExtrapolation
Definition CommonPars.h:26
#define ARR2D_FS(name, N, i, j)
Definition CommonPars.h:29
#define vkalNTrkM
Definition CommonPars.h:22
#define xyz
#define x
void setRefIterV(double v[]) noexcept
void setCnstV(double v[3]) noexcept
void setIniV(double v[3]) noexcept
std::vector< std::unique_ptr< VKTrack > > TrackList
double ader[(3 *vkalNTrkM+3) *(3 *vkalNTrkM+3)]
std::vector< std::unique_ptr< VKConstraintBase > > ConstraintList
double fitCovXYZMom[21]
void setRefV(double v[3]) noexcept
std::vector< std::unique_ptr< TWRK > > tmpArr
std::unique_ptr< VKalVrtControl > vk_fitterControl
void renewFullCovariance(double *)
void setVertexMass(double mass)
void setVrtMassError(double error)
const double * getFullCovariance() const
static void getMagFld(const double, const double, const double, double &, double &, double &, const VKalVrtControlBase *)
static bool checkTarget(double *RefEnd)
static void Propagate(long int TrkID, long int Charge, double *ParOld, double *CovOld, double *RefStart, double *RefEnd, double *ParNew, double *CovNew, VKalVrtControlBase *FitControl=0)
double chi2(TH1 *h0, TH1 *h1)
Ensure that the ATLAS eigen extensions are properly loaded.
void vpderiv(bool UseTrackErr, long int Charge, const double *pari0, double *covi, double *vrtref, double *covvrtref, double *drdpar, double *dwgt, double *rv0, VKalVrtControl *FitCONTROL)
Definition VtDeriv.cxx:16
void cfdcopy(double *source, double *target, int n)
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
int cfInv5(double *cov, double *wgt)
double cfddist3D(double *V1, double *V2)
int cfdinv(double *cov, double *wgt, long int NI)
void cfmasserr(VKVertex *vk, const int *list, double *MASS, double *sigM)
Definition cfMassErr.cxx:16
int vtcfit(VKVertex *vk)
Definition VtCFit.cxx:293
std::array< double, 4 > getFitParticleMom(const VKTrack *trk, const VKVertex *vk)
void applyConstraints(VKVertex *vk)
Definition stCnst.cxx:22
void robtest(VKVertex *vk, int ifl, int nIteration)
Definition RobTest.cxx:14
int afterFit(VKVertex *vk, double *ader, double *dcv, double *ptot, double *VrtMomCov, const VKalVrtControlBase *CONTROL)
Definition cfTotCov.cxx:26
double cfSmallEigenvalue(double *cov, long int n)
int ic
Definition grepfile.py:33
Definition index.py:1
int useMassCnst
Definition ForCFT.h:15
int usePhiCnst
Definition ForCFT.h:16
int useRadiusCnst
Definition ForCFT.h:22
double robres[vkalNTrkM]
Definition ForCFT.h:51
int indtrkmc[vkalMaxNMassCnst][vkalNTrkM]
Definition ForCFT.h:52
double Dp
Definition ForCFT.h:34
int useThetaCnst
Definition ForCFT.h:17
double RC
Definition ForCFT.h:35
int usePlaneCnst
Definition ForCFT.h:21
double IterationPrecision
Definition ForCFT.h:54
double covvrt[6]
Definition ForCFT.h:45
double Cp
Definition ForCFT.h:33
int irob
Definition ForCFT.h:49
double radiusRefP[2]
Definition ForCFT.h:36
double Ap
Definition ForCFT.h:31
int useAprioriVrt
Definition ForCFT.h:19
int IterationNumber
Definition ForCFT.h:53
double localbmag
Definition ForCFT.h:29
double wm[vkalNTrkM]
Definition ForCFT.h:27
int nmcnst
Definition ForCFT.h:26
int usePointingCnst
Definition ForCFT.h:18
int usePassNear
Definition ForCFT.h:20
double wmfit[vkalMaxNMassCnst]
Definition ForCFT.h:28
double Bp
Definition ForCFT.h:32
double vrt[3]
Definition ForCFT.h:45
double covvrt[6]
Definition ForVrtClose.h:23
double cvder[12]
Definition ForVrtClose.h:30
double dcovf[21]
Definition ForVrtClose.h:24
double dcv[6 *(vkalNTrkM *3+3)]
Definition ForVrtClose.h:31