ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitCommon/VertexFit/src/MultiVertexFitter.cc
Revision: 1.4
Committed: Fri Mar 20 13:33:03 2009 UTC (16 years, 1 month ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_011, Mit_010a, Mit_010, Mit_009c, Mit_009b, Mit_009a, Mit_009, Mit_008
Changes since 1.3: +71 -59 lines
Log Message:
Cleanup

File Contents

# User Rev Content
1 loizides 1.4 // $Id: MultiVertexFitter.cc,v 1.3 2008/11/13 16:34:29 paus Exp $
2 loizides 1.1
3     #include "MitCommon/VertexFit/interface/MultiVertexFitter.h"
4     #include <algorithm>
5     #include <math.h>
6     #include <iostream>
7     #include <csignal>
8     #include <csetjmp>
9     #include <TMath.h>
10     #include <CLHEP/Matrix/Matrix.h>
11    
12     extern "C" void ctvmft_(int&, int&, int&);
13     extern "C" bool mcalc_ (int&, int*, float&, float&, double*);
14     extern "C" void dcalc_ (int&, int&, float*, float&, float&, float*);
15    
16     using namespace std;
17     using namespace mithep;
18    
19     jmp_buf env;
20 loizides 1.4
21     //--------------------------------------------------------------------------------------------------
22 loizides 1.1 extern "C" void MultiVertexFitterSetStatus(int i) {
23     cout << "Warning, you are handling a severe error in MultiVertexFitter" << endl;
24     longjmp(env,-66);
25     }
26    
27 loizides 1.4 //--------------------------------------------------------------------------------------------------
28 loizides 1.1 MultiVertexFitter::MultiVertexFitter() :
29     _bField (3.8), // default B field for running
30     _currentAllocatedVertexNumber(0), // facilitates CandNode recursion
31     _referencePoint (0,0,0), // set reference point to (0,0,0) initially
32     _primaryVertex (0,0,0), // set primary vertex to (0,0,0) initially
33     _cdfPrimaryVertex (0,0,0) // set pv in CMS coords to (0,0,0) initially
34     {
35     // Set name and email of MultiVertexFitter expert
36     _expert="Christoph Paus (paus@mit.edu)";
37    
38     // First get pointers to various FORTAN common blocks
39     _ctvmq_com = (CTVMQ*) ctvmq_address_(); //printf(" Common: _ctvmq_com %p\n", _ctvmq_com );
40     _ctvmfr_com = (CTVMFR*) ctvmfr_address_(); //printf(" Common: _ctvmfr_com %p\n", _ctvmfr_com);
41     _fiddle_com = (FIDDLE*) fiddle_address_(); //printf(" Common: _fiddle_com %p\n", _fiddle_com);
42     _trkprm_com = (TRKPRM*) trkprm_address_(); //printf(" Common: _trkprm_com %p\n", _trkprm_com);
43    
44     // Initialize various arrays
45     init();
46     // Don't bomb program on error.
47     _fiddle.excuse = 1;
48     // Don't extrapolate track errors by default
49     _extrapolateTrackErrors = false;
50     }
51    
52 loizides 1.4 //--------------------------------------------------------------------------------------------------
53 loizides 1.1 void MultiVertexFitter::init(double bField)
54     {
55     // Set internal variable which keeps track of the b field
56     _bField = bField;
57    
58     // Now initialize CTVMQ common. Run and trigger numbers are dummies - they are not used.
59     _ctvmq.runnum = 1;
60     _ctvmq.trgnum = 100;
61     // Eventually, we have to get the magnetic field from the right place.
62     // Origignal with bmag = 14.116 [kGauss]:
63     // _ctvmq.pscale = 0.000149896 * bmag;
64     // New bfield in Tesla bmag = 1.4116 [T]:
65     _ctvmq.pscale = 0.00149896 * bField;
66     // Set the default maximum chi-square per dof.
67     _fiddle.chisqmax = 225.0;
68     // We also need to get the primary vertex from the right place, but for now we put in (0,0,0).
69     setPrimaryVertex(0.0,0.0,0.0);
70     float xverr[3][3];
71     for (int j = 0; j<3; j++) {
72     for (int k = 0; k<3; k++) {
73     xverr[j][k] = 0.;
74     }
75     }
76     xverr[0][0] = 0.005;
77     xverr[1][1] = 0.005;
78     xverr[2][2] = 1.0;
79     setPrimaryVertexError(xverr);
80     // Zero number of tracks, vertices and mass constraints
81     _ctvmq.ntrack = 0;
82     _ctvmq.nvertx = 0;
83     _ctvmq.nmassc = 0;
84     // Zero track list and arrays containing the vertex and mass constraint configuration
85     for (int j=0; j<_maxtrk; ++j) {
86     _ctvmq.list[j]=0;
87     for (int jv=0; jv<_maxvtx; ++jv)
88     _ctvmq.trkvtx[jv][j] = false;
89     for (int jmc=0; jmc<_maxmcn; ++jmc)
90     _ctvmq.trkmcn[jmc][j] = false;
91     }
92     // Initialize the conversion and vertex pointing arrays
93     for (int jv=0; jv<_maxvtx; ++jv) {
94     _ctvmq.cvtx[jv] = 0;
95     _ctvmq.vtxpnt[0][jv] = -1;
96     _ctvmq.vtxpnt[1][jv] = 0;
97     }
98     _ctvmq.drmax = 2.0;
99     _ctvmq.dzmax = 20.0;
100     _ctvmq.rvmax = 70.0;
101     _ctvmq.trnmax = 0.5;
102     _ctvmq.dsmin = -100.0; // for CDF used -2, for CMS fuzzing around.
103     // Set _stat to -999 and chisqq to -1 to symbolize that no fit has yet been done.
104     _stat = -999;
105     _ctvmq.chisqr[0] = -1.0;
106    
107     _primaryVertex = Hep3Vector(0,0,0);
108     _cdfPrimaryVertex = Hep3Vector(0,0,0);
109     _referencePoint = ThreeVector(0,0,0);
110     }
111    
112 loizides 1.4 //--------------------------------------------------------------------------------------------------
113 loizides 1.1 bool MultiVertexFitter::addTrack(const HepVector &v, const HepSymMatrix &cov,
114     int trackId, float mass, vertexNumber jv)
115     {
116     // Check that this vertex number is within the allowed range.
117     if (jv<VERTEX_1 || jv>_maxvtx)
118     return false;
119     _ctvmq.nvertx = jv>_ctvmq.nvertx ? jv : _ctvmq.nvertx;
120    
121     // Add track extrapolation, if the user set it
122     HepVector w = v;
123     HepSymMatrix m = cov;
124     //if (_extrapolateTrackErrors)
125     // // This will move the reference point
126     // moveReferencePoint(w, m);
127    
128     // Check that we have not exceeded the maximum number of tracks.
129     if (_ctvmq.ntrack>=_maxtrk)
130     return false;
131    
132     // Add this track
133     _ctvmq.list [_ctvmq.ntrack] = trackId;
134     _ctvmq.tkbank[_ctvmq.ntrack][0] = 'Q';
135     _ctvmq.tkbank[_ctvmq.ntrack][1] = 'T';
136     _ctvmq.tkbank[_ctvmq.ntrack][2] = 'R';
137     _ctvmq.tkbank[_ctvmq.ntrack][3] = 'K';
138     _ctvmq.tmass [_ctvmq.ntrack] = mass;
139     _ctvmq.trkvtx[jv-1][_ctvmq.ntrack] = true;
140    
141     // Put this track's helix parameters and error matrix into a fortran common block so that they
142     // can be accessed by gettrk. This is a dummy for now.
143     _trkprm.trhelix[_ctvmq.ntrack][0] = w[0];
144     _trkprm.trhelix[_ctvmq.ntrack][1] = w[1];
145     _trkprm.trhelix[_ctvmq.ntrack][2] = w[2];
146     _trkprm.trhelix[_ctvmq.ntrack][3] = w[3];
147     _trkprm.trhelix[_ctvmq.ntrack][4] = w[4];
148    
149     for (int j=0; j<5; ++j) {
150     for (int k=0; k<5; ++k)
151     _trkprm.trem[_ctvmq.ntrack][j][k]=m[j][k];
152     }
153     _ctvmq.ntrack++;
154    
155     return true;
156     }
157    
158 loizides 1.4 //--------------------------------------------------------------------------------------------------
159 loizides 1.1 bool MultiVertexFitter::addTrack(const TVectorD &v, const TMatrixDSym &cov,
160     int trackId, double mass, vertexNumber jv)
161     {
162     // Check that this vertex number is within the allowed range.
163     if (jv<VERTEX_1 || jv>_maxvtx)
164     return false;
165     _ctvmq.nvertx = jv>_ctvmq.nvertx ? jv : _ctvmq.nvertx;
166    
167     //// Add track extrapolation, if the user set it
168     //HepVector w = v;
169     //HepSymMatrix m = cov;
170     TVectorD w(v);
171     TMatrixDSym m(cov);
172     ////if (_extrapolateTrackErrors)
173     //// // This will move the reference point
174     //// moveReferencePoint(w, m);
175    
176     // Check that we have not exceeded the maximum number of tracks.
177     if (_ctvmq.ntrack>=_maxtrk)
178     return false;
179    
180     // Add this track
181     _ctvmq.list [_ctvmq.ntrack] = trackId;
182     _ctvmq.tkbank[_ctvmq.ntrack][0] = 'Q';
183     _ctvmq.tkbank[_ctvmq.ntrack][1] = 'T';
184     _ctvmq.tkbank[_ctvmq.ntrack][2] = 'R';
185     _ctvmq.tkbank[_ctvmq.ntrack][3] = 'K';
186     _ctvmq.tmass [_ctvmq.ntrack] = mass;
187     _ctvmq.trkvtx[jv-1][_ctvmq.ntrack] = true;
188    
189     // Put this track's helix parameters and error matrix into a fortran common block so that they
190     // can be accessed by gettrk. This is a dummy for now.
191     _trkprm.trhelix[_ctvmq.ntrack][0] = w[0];
192     _trkprm.trhelix[_ctvmq.ntrack][1] = w[1];
193     _trkprm.trhelix[_ctvmq.ntrack][2] = w[2];
194     _trkprm.trhelix[_ctvmq.ntrack][3] = w[3];
195     _trkprm.trhelix[_ctvmq.ntrack][4] = w[4];
196    
197     for (int j=0; j<5; ++j) {
198     for (int k=0; k<5; ++k)
199     _trkprm.trem[_ctvmq.ntrack][j][k]=m[j][k];
200     }
201     _ctvmq.ntrack++;
202    
203     return true;
204     }
205    
206 loizides 1.4 //--------------------------------------------------------------------------------------------------
207 loizides 1.1 bool MultiVertexFitter::vertexPoint_2d(vertexNumber jv1, vertexNumber jv2)
208     {
209     // Check that these vertex numbers are within allowed range and that the vertices are unique.
210     if (jv1>_maxvtx || jv1<VERTEX_1)
211     return false;
212     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
213     return false;
214     if (jv1 <= jv2)
215     return false;
216    
217     // Setup vertex pointing.
218     _ctvmq.vtxpnt[0][jv1-1] = jv2;
219     _ctvmq.vtxpnt[1][jv1-1] = 1; // 2d pointing.
220    
221     return true;
222     }
223    
224 loizides 1.4 //--------------------------------------------------------------------------------------------------
225 loizides 1.1 bool MultiVertexFitter::vertexPoint_3d(vertexNumber jv1, vertexNumber jv2)
226     {
227     // Check that these vertex numbers are within allowed range and that the vertices are distinct
228     if (jv1>_maxvtx || jv1<VERTEX_1)
229     return false;
230     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
231     return false;
232     if (jv1 <= jv2)
233     return false;
234    
235     // Setup vertex pointing
236     _ctvmq.vtxpnt[0][jv1-1] = jv2;
237     _ctvmq.vtxpnt[1][jv1-1] = 2; // 3d pointing
238    
239     return true;
240     }
241    
242 loizides 1.4 //--------------------------------------------------------------------------------------------------
243 loizides 1.1 bool MultiVertexFitter::vertexPoint_1track(vertexNumber jv1, vertexNumber jv2)
244     {
245     // Check that these vertex numbers are within allowed range and are distinct
246     if (jv1>_maxvtx || jv1<VERTEX_1)
247     return false;
248     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
249     return false;
250     if (jv1 <= jv2)
251     return false;
252    
253     // Setup vertex pointing
254     _ctvmq.vtxpnt[0][jv1-1] = jv2;
255     _ctvmq.vtxpnt[1][jv1-1] = 3; // Point to 1 track vertex
256    
257     return true;
258     }
259    
260 loizides 1.4 //--------------------------------------------------------------------------------------------------
261 loizides 1.1 bool MultiVertexFitter::vertexPoint_0track(vertexNumber jv1, vertexNumber jv2)
262     {
263     // jv2 is the zero track vertex. jv1 is the multi track vertex which points to jv2
264    
265     // Note: You must call this routine at least twice in order for ctvmft_zerotrackvtx to work ie
266     // there must be at least 2 vertices pointed at a zero track vertex. ctvmft for this, but the
267     // error message may not make it to your log file (look at the local variables in the stack frame
268     // especially IJKERR(2). The significance of this error code is documented at the top of whatever
269     // routine chucked you out (ctvmf00 in this case)
270    
271     // see ctvmft.f source file for discussion. See especially comments at the top of subroutines:
272     // ctvmft and ctvmfa
273    
274     // Check that these vertex numbers are within allowed range and are distinct.
275     if (jv1>_maxvtx || jv1<VERTEX_1)
276     return false;
277     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
278     return false;
279     if (jv1 <= jv2)
280     return false;
281    
282     // Setup vertex pointing.
283     _ctvmq.vtxpnt[0][jv1-1] = jv2;
284     _ctvmq.vtxpnt[1][jv1-1] = 4; // Point to 0 track vertex.
285    
286     return true;
287     }
288    
289 loizides 1.4 //--------------------------------------------------------------------------------------------------
290 loizides 1.1 bool MultiVertexFitter::conversion_2d(vertexNumber jv)
291     {
292     if (jv<VERTEX_1 || jv>_ctvmq.nvertx)
293     return false;
294    
295     _ctvmq.cvtx[jv-1] = 1;
296    
297     return true;
298     }
299    
300 loizides 1.4 //--------------------------------------------------------------------------------------------------
301 loizides 1.1 bool MultiVertexFitter::conversion_3d(vertexNumber jv)
302     {
303     if (jv<VERTEX_1 || jv>_ctvmq.nvertx)
304     return false;
305    
306     _ctvmq.cvtx[jv-1] = 2;
307    
308     return true;
309     }
310    
311 loizides 1.4 //--------------------------------------------------------------------------------------------------
312 loizides 1.1 bool MultiVertexFitter::massConstrain(int ntrk, const int trkIds[], float mass)
313     {
314     // Check that we have not exceeded the allowed number of mass constraints.
315     if (_ctvmq.nmassc>=_maxmcn)
316     return false;
317    
318     // Set constraint mass
319     _ctvmq.cmass[_ctvmq.nmassc]=mass;
320    
321     // For each track in contraint, set trkmcn true. Since the number in tracks[] is the track
322     // number, we have to find each track in the list of tracks.
323     for (int jt=0; jt<ntrk; ++jt) {
324     bool found=false;
325     for (int kt=0; kt<_ctvmq.ntrack; ++kt) {
326     if (trkIds[jt] == _ctvmq.list[kt]) {
327     _ctvmq.trkmcn[_ctvmq.nmassc][kt]=true;
328     found=true;
329     }
330     }
331     if (!found)
332     return false;
333     }
334    
335     // Increment number of mass constraints.
336     _ctvmq.nmassc++;
337    
338     return true;
339     }
340    
341 loizides 1.4 //--------------------------------------------------------------------------------------------------
342 loizides 1.1 bool MultiVertexFitter::beamlineConstraint(float xb, float yb, HepSymMatrix berr,
343     float xzbslope, float yzbslope)
344     {
345     // Set beam position at z=0
346     setPrimaryVertex(xb,yb,0);
347     //if (_extrapolateTrackErrors) {
348     // float newXb = xb - _referencePoint.x() + _referencePoint.z() * xzbslope;
349     // float newYb = yb - _referencePoint.y() + _referencePoint.z() * yzbslope;
350     // setPrimaryVertex(newXb, newYb, 0);
351     //}
352    
353     bool success = setPrimaryVertexError(berr);
354    
355     // Set the beamline slope values
356     _ctvmq.xzslope = xzbslope;
357     _ctvmq.yzslope = yzbslope;
358    
359     // Turn ON beamline constraint
360     _ctvmq.vtxpnt[0][0] = -100;
361    
362     return success;
363     }
364    
365 loizides 1.4 //--------------------------------------------------------------------------------------------------
366 loizides 1.1 bool MultiVertexFitter::beamlineConstraint(Hep3Vector pv, HepSymMatrix berr, float xzbslope,
367     float yzbslope)
368     {
369     // Check if input beam position coordinates are at z=0
370     if (pv.z() != 0)
371     return false;
372    
373     return beamlineConstraint(pv.x(),pv.y(),berr,xzbslope,yzbslope);
374     }
375    
376 loizides 1.4 //--------------------------------------------------------------------------------------------------
377 loizides 1.1 void MultiVertexFitter::setPrimaryVertex(float xv, float yv, float zv)
378     {
379     // Set x,y,z position of the primary vertex.
380     _ctvmq.xyzpv0[0] = xv;
381     _ctvmq.xyzpv0[1] = yv;
382     _ctvmq.xyzpv0[2] = zv;
383    
384     _primaryVertex = Hep3Vector( xv, yv, zv );
385     }
386    
387 loizides 1.4 //--------------------------------------------------------------------------------------------------
388 loizides 1.1 void MultiVertexFitter::setPrimaryVertex(Hep3Vector pv)
389     {
390     // Set x,y,z position of the primary vertex.
391     _ctvmq.xyzpv0[0] = pv.x();
392     _ctvmq.xyzpv0[1] = pv.y();
393     _ctvmq.xyzpv0[2] = pv.z();
394    
395     _primaryVertex = pv;
396     }
397    
398 loizides 1.4 //--------------------------------------------------------------------------------------------------
399 loizides 1.1 void MultiVertexFitter::setPrimaryVertexError(const float xverr[3][3])
400     {
401     // Set the error matrix for the primary vertex.
402     for (int j=0; j<3; ++j) {
403     for (int k=0; k<3; ++k)
404     _ctvmq.exyzpv[j][k]=xverr[j][k];
405     }
406     }
407    
408 loizides 1.4 //--------------------------------------------------------------------------------------------------
409 loizides 1.1 bool MultiVertexFitter::setPrimaryVertexError(const HepSymMatrix &xverr)
410     {
411     // Set the error matrix for the primary vertex using a HepSymMatrix. First check that the matrix
412     // is the correct size.
413     if (xverr.num_row() != 3)
414     return false;
415     for (int j=0; j<3; j++) {
416     for (int k=0; k<3; k++)
417     _ctvmq.exyzpv[j][k]=xverr[j][k];
418     }
419     return true;
420     }
421    
422 loizides 1.4 //--------------------------------------------------------------------------------------------------
423 loizides 1.1 bool MultiVertexFitter::fit()
424     {
425     // Check that the diagonal elements of all the track error matrices are positive
426     bool mstat = true;
427     for (int trk=0; trk<_ctvmq.ntrack; ++trk) {
428     for (int j=0; j<5; ++j) {
429     // Check diagonal elements of error matrix.
430     if (_trkprm.trem[trk][j][j] < 0.) {
431     // The covariance matrix could not be inverted: Set the error codes and fail this fit
432     mstat = false;
433     _ctvmq.ijkerr[0] = 3;
434     _ctvmq.ijkerr[1] = 2;
435     _ctvmq.ijkerr[2] = trk + 1;
436     }
437     }
438     // Check that curvature of track is reasonable: Pt is above ~10MeV/c. If not, set the error
439     // codes and fail this fit
440     if (fabs(_trkprm.trhelix[trk][1]) > 0.1) {
441     //if (fabs(_trkprm.trhelix[trk][1]) > 0.01) {
442     mstat = false;
443     _ctvmq.ijkerr[0] = 3;
444     _ctvmq.ijkerr[1] = 5;
445     _ctvmq.ijkerr[2] = trk + 1;
446     }
447     }
448     // If there was a problem with any track, fail the fit
449     if (!mstat) {
450     _stat = 1;
451     return false;
452     }
453    
454     // First copy information into CTVMFT common blocks
455     *_ctvmq_com = _ctvmq;
456     *_ctvmfr_com = _ctvmfr;
457     *_fiddle_com = _fiddle;
458     *_trkprm_com = _trkprm;
459     // Do the vertex fit.
460     int print = 0;
461     int level = 0;
462    
463     // #if ( defined(LINUX) && defined(__USE_BSD) ) || defined(OSF1)
464     // struct sigaction myaction = {MultiVertexFitterSetStatus, 0, 0, 0}, oldaction;
465     // sigaction(SIGFPE, &myaction, &oldaction);
466     // if (setjmp(env)!=0) {
467     // sigaction(SIGFPE, &oldaction,0);
468     // return -999;
469     // }
470     // #endif
471    
472     ctvmft_(print,level,_stat);
473    
474     // #if ( defined(LINUX) && defined(__USE_BSD) ) || defined(OSF1)
475     // sigaction(SIGFPE, &oldaction,0);
476     // #endif
477    
478     // Now copy information from CTVMFT common blocks to local storage
479     _ctvmq = *_ctvmq_com;
480     _ctvmfr = *_ctvmfr_com;
481     _fiddle = *_fiddle_com;
482     _trkprm = *_trkprm_com;
483    
484     return (_stat == 0);
485     }
486    
487 loizides 1.4 //--------------------------------------------------------------------------------------------------
488 loizides 1.1 void MultiVertexFitter::print() const
489     {
490     print(cout);
491     }
492    
493 loizides 1.4 //--------------------------------------------------------------------------------------------------
494 loizides 1.1 void MultiVertexFitter::print(ostream& os) const
495     {
496     os << "****************************** MultiVertexFitter "
497     << "******************************" << endl;
498     os << "Number of tracks: " << _ctvmq.ntrack << endl;
499     os << " Tracks: ";
500     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
501     if (jt != 0) os << ", ";
502     os << _ctvmq.list[jt];
503     }
504     os << endl;
505     os << "Number of vertices: " << _ctvmq.nvertx << endl;
506     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
507     os << " Vertex " << jv+1 << " tracks: ";
508     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
509     if (_ctvmq.trkvtx[jv][jt]) {
510     os << " " << _ctvmq.list[jt];
511     }
512     }
513     os << endl;
514     }
515     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
516     if (_ctvmq.vtxpnt[0][jv]==0) {
517     os << " Vertex " << jv+1 << " points to the primary vertex ";
518     }
519     else if (_ctvmq.vtxpnt[0][jv]>0) {
520     os << " Vertex " << jv+1 << " points to vertex "
521     << _ctvmq.vtxpnt[0][jv];
522     }
523     if (_ctvmq.vtxpnt[1][jv]==1) {
524     os << " in 2 dimensions" << endl;
525     }
526     else if (_ctvmq.vtxpnt[1][jv]==2) {
527     os << " in 3 dimensions" << endl;
528     }
529     else if (_ctvmq.vtxpnt[1][jv]==3) {
530     os << ", a single track vertex" << endl;
531     }
532     if (_ctvmq.cvtx[jv]>0) {
533     os << " Vertex " << jv+1 << " is a conversion" << endl;
534     }
535     }
536     os << "Number of mass constraints: " << _ctvmq.nmassc << endl;
537     for (int jmc=0; jmc<_ctvmq.nmassc; ++jmc) {
538     os << " Tracks ";
539     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
540     if (_ctvmq.trkmcn[jmc][jt]) {
541     os << " " << _ctvmq.list[jt];
542     }
543     }
544     os << " constrained to mass " << _ctvmq.cmass[jmc]
545     << " Gev/c^2" << endl;
546     }
547     if (_stat==-999) {
548     os << "No fit has been done." << endl;
549     }
550     else {
551     os << "***** Results of Fit *****" << endl;
552     printErr(os);
553     os << " Status = " << _stat << endl;
554     os.precision(7);
555     os << " Chi-square = " << scientific << _ctvmq.chisqr[0]
556     << " for " << _ctvmq.ndof << " degrees of freedom." << endl;
557     os << " => probability = " << prob() << endl;
558     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
559     os << "Vertex " << jv+1
560     << " position: " << scientific
561     << _ctvmq.xyzvrt[jv+1][0] << " "
562     << _ctvmq.xyzvrt[jv+1][1] << " "
563     << _ctvmq.xyzvrt[jv+1][2] << endl;
564     }
565     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
566     os << "Track " << _ctvmq.list[jt]
567     << " - P4: " << scientific
568     << _ctvmq.trkp4[0][jt] << " "
569     << _ctvmq.trkp4[1][jt] << " "
570     << _ctvmq.trkp4[2][jt] << " "
571     << _ctvmq.trkp4[3][jt] << " "
572     << " - PT: " << scientific
573     << sqrt(_ctvmq.trkp4[0][jt]*_ctvmq.trkp4[0][jt]+
574     _ctvmq.trkp4[1][jt]*_ctvmq.trkp4[1][jt]) << endl;
575     }
576     }
577     os << "****************************************"
578     << "**************************" << endl;
579    
580     return;
581     }
582    
583 loizides 1.4 //--------------------------------------------------------------------------------------------------
584 loizides 1.1 void MultiVertexFitter::printErr() const
585     {
586     printErr(cout);
587     }
588    
589 loizides 1.4 //--------------------------------------------------------------------------------------------------
590 loizides 1.1 void MultiVertexFitter::printErr(ostream& os) const
591     {
592     os << "MultiVertexFitter: IJKERR = " << _ctvmq.ijkerr[0] << ", "
593     << _ctvmq.ijkerr[1] << ", "
594     << _ctvmq.ijkerr[2] << endl;
595     if (status()==0 && _ctvmq.ijkerr[0]==0) return;
596     if (_ctvmq.ijkerr[0] == -1) {
597     os << " Problem with GETTRK: track requested is not in list."
598     << endl
599     << " This should not happen - Contact MultiVertexFitter expert "
600     << _expert << "." <<endl;
601     }
602     else if (_ctvmq.ijkerr[0]==1) {
603     os << " Problem in CTVM00:" << endl;
604     if (_ctvmq.ijkerr[1]==1) {
605     os << " Number of tracks is " << _ctvmq.ntrack
606     << "." << endl;
607     if (_ctvmq.ntrack < 2) {
608     os << ", which is too few (must be at least 2)." << endl;
609     }
610     else if (_ctvmq.ntrack > _maxtrk) {
611     os << ", which is too many (maximum is " << _maxtrk
612     << ")." << endl;
613     }
614     else {
615     os << " Problem with number of tracks"
616     << " for unknown reasons." << endl;
617     }
618     }
619     else if (_ctvmq.ijkerr[1]==2) {
620     os << " Number of vertices is " << _ctvmq.nvertx
621     << "." << endl;
622     if (_ctvmq.nvertx < 1) {
623     os << ", which is too few (must be at least 1)." << endl;
624     }
625     else if (_ctvmq.nvertx > _maxvtx) {
626     os << ", which is too many (maximum is " << _maxvtx
627     << ")." << endl;
628     }
629     else {
630     os << endl << " Problem with number of vertices"
631     << " for unknown reasons." << endl;
632     }
633     }
634     else if (_ctvmq.ijkerr[1]==3) {
635     os << " Number of mass constraints is " << _ctvmq.nmassc
636     << "." << endl;
637     if (_ctvmq.nmassc < 0) {
638     os << ", which is negative." << endl;
639     }
640     else if (_ctvmq.nmassc > _maxmcn) {
641     os << ", which is too many (maximum is " << _maxmcn
642     << ")." << endl;
643     }
644     else {
645     os << endl << " Problem with number of mass"
646     << " constraints for unknown reasons." << endl;
647     }
648     }
649     else if (_ctvmq.ijkerr[1]==11) {
650     os << " Vertex " << _ctvmq.ijkerr[2]
651     << " has less than one track." << endl;
652     }
653     else if (_ctvmq.ijkerr[1]==12) {
654     os << " Vertex " << _ctvmq.ijkerr[2]
655     << " is a conversion vertex with a number of tracks"
656     << " different than two." << endl;
657     }
658     else if (_ctvmq.ijkerr[1]==13) {
659     os << " Vertex " << _ctvmq.ijkerr[2]
660     << " is a one track vertex that has no multi-track"
661     << " descendents." << endl;
662     }
663     else if (_ctvmq.ijkerr[1]==14) {
664     os << " Vertex " << _ctvmq.ijkerr[2]
665     << " does not point at a vertex with a lower number."
666     << endl;
667     }
668     else if (_ctvmq.ijkerr[1]==15) {
669     os << " Vertex " << _ctvmq.ijkerr[2]
670     << " has a parent vertex that is a conversion." << endl;
671     }
672     else if (_ctvmq.ijkerr[1]==16) {
673     os << " Vertex " << _ctvmq.ijkerr[2]
674     << " does 1 track pointing to a vertex with"
675     << " more than 1 track." << endl;
676     }
677     else if (_ctvmq.ijkerr[1]==17) {
678     os << " Vertex " << _ctvmq.ijkerr[2]
679     << " does 0 track pointing to a vertex with"
680     << " more than 0 track (?)." << endl;
681     }
682     else if (_ctvmq.ijkerr[1]==19) {
683     os << " Primary vertex error matrix is singular." << endl;
684     }
685     else if (_ctvmq.ijkerr[1]==21) {
686     os << " Track with Id " << _ctvmq.ijkerr[2]
687     << "is not in any vertex." << endl;
688     }
689     else if (_ctvmq.ijkerr[1]==22) {
690     os << " Track with Id " << _ctvmq.ijkerr[2]
691     << "is in multiple vertices." << endl;
692     }
693     else if (_ctvmq.ijkerr[1]==23) {
694     os << " Track with Id " << _ctvmq.ijkerr[2]
695     << "occurs more than once." << endl;
696     }
697     else if (_ctvmq.ijkerr[1]==31) {
698     os << " A mass constraint has less than 2 tracks." << endl;
699     }
700     else if (_ctvmq.ijkerr[1]==32) {
701     os << " The sum masses of the tracks in a mass constraint"
702     << " exceeds the constraint mass." << endl;
703     }
704     else if (_ctvmq.ijkerr[1]==33) {
705     os << " Beamline constraint. Beam covariance not set properly."
706     << " Negative diagonal elements." << endl;
707     }
708     else if (_ctvmq.ijkerr[1]==34) {
709     os << " Beamline constraint. Beam covariance not set properly."
710     << " Off-diagonal elements not zero." << endl;
711     }
712     else if (_ctvmq.ijkerr[1]==36) {
713     os << " Beamline constraint. Number of vertices = "
714     << _ctvmq.nvertx << " Should be 1." << endl;
715     }
716     }
717     else if (_ctvmq.ijkerr[0] == 2) {
718     if (_ctvmq.ijkerr[1] == 20) {
719     os << " Problem in CTVM00: " << endl;
720     os << " Track has negative Id = "
721     << _ctvmq.list[_ctvmq.ijkerr[2]-1] << "." << endl;
722     }
723     else {
724     os << " Problem in CTVMFA with vertex "
725     << _ctvmq.ijkerr[2] << ": " << endl;
726     os << " Failure in vertex first approximation." << endl;
727     if (_ctvmq.ijkerr[1] == 1) {
728     os << " Tracks are concentric circles." << endl;
729     }
730     if (_ctvmq.ijkerr[1] == 2) {
731     os << " Conversion vertex has widely separated"
732     << " exterior circles at midpoint." << endl;
733     }
734     if (_ctvmq.ijkerr[1] == 3) {
735     os << " Conversion vertex has widely separated"
736     << " interior circles at midpoint." << endl;
737     }
738     if (_ctvmq.ijkerr[1] == 4) {
739     os << " Vertex has widely separated"
740     << " exterior circles at approximate vertex." << endl;
741     }
742     if (_ctvmq.ijkerr[1] == 5) {
743     os << " Vertex has widely separated"
744     << " interior circles at approximate vertex." << endl;
745     }
746     if (_ctvmq.ijkerr[1] == 6) {
747     os << " Rv is too large at the chosen"
748     << " intersection point." << endl;
749     }
750     if (_ctvmq.ijkerr[1] == 7) {
751     os << " Delta z is too large at the chosen"
752     << " intersection point." << endl;
753     }
754     if (_ctvmq.ijkerr[1] == 8) {
755     os << " A track's turning to the chosen vertex"
756     << " is too large." << endl;
757     }
758     if (_ctvmq.ijkerr[1] == 9) {
759     os << " There is no solution with an adequately"
760     << " positive arc length." << endl;
761     }
762     if (_ctvmq.ijkerr[1] == 21) {
763     os << " zero-track vertexing: either/both vertex "
764     << " momenta are too small (<0.01 MeV)." << endl;
765     }
766     if (_ctvmq.ijkerr[1] == 22) {
767     os << " zero-track vertexing: Two lines (tracks) are "
768     << " parallel/antiparallel." << endl;
769     }
770    
771     }
772     }
773     else if (_ctvmq.ijkerr[0] == 3) {
774     os << " Problem in CTVM01 with track with Id = "
775     << _ctvmq.list[_ctvmq.ijkerr[2]-1] << ": " << endl;
776     if (_ctvmq.ijkerr[1] == 1) {
777     os << " GETTRK cannot find Id in list." << endl;
778     }
779     if (_ctvmq.ijkerr[1] == 2) {
780     os << " Covariance matrix could not be inverted." << endl
781     << " Offending track number (in order addded) is "
782     << _ctvmq.ijkerr[2] << "." << endl;
783     }
784     if (_ctvmq.ijkerr[1] == 3) {
785     os << " Track turns through too large an angle"
786     << " to the vertex." << endl;
787     }
788     if (_ctvmq.ijkerr[1] == 4) {
789     os << " Track moves too far backward to vertex." << endl;
790     }
791     if (_ctvmq.ijkerr[1] == 5) {
792     os << " Track with curvature > 0.01." << endl
793     << " Offending track number is "
794     << _ctvmq.ijkerr[2] << "." << endl;
795     }
796     }
797     else if (status() == 9) {
798     os << " General fit problem: " << endl;
799     if (_ctvmq.ijkerr[1] == 1) {
800     os << " Singular solution matrix." << endl;
801     }
802     if (_ctvmq.ijkerr[1] == 2 || _ctvmq.ijkerr[1] == 3) {
803     os << " Too many iterations ( "
804     << _ctvmq.ijkerr[2] << "(." << endl;
805     }
806     if (_ctvmq.ijkerr[1] == 4) {
807     os << " Convergence failure." << endl;
808     }
809     if (_ctvmq.ijkerr[1] == 5) {
810     os << " Bad convergence." << endl;
811     }
812     if (_ctvmq.ijkerr[1] == 9) {
813     os << " Ill-formed covariance matrix." << endl;
814     }
815     }
816     else {
817     os << " The error codes above are not recognized." << endl
818     << " Contact MultiVertexFitter expert " << _expert << "." << endl;
819     }
820     return;
821     }
822    
823 loizides 1.4 //--------------------------------------------------------------------------------------------------
824 loizides 1.1 void MultiVertexFitter::getIJKErr(int& err0, int& err1, int& err2) const
825     {
826     err0 = _ctvmq.ijkerr[0];
827     err1 = _ctvmq.ijkerr[1];
828     err2 = _ctvmq.ijkerr[2];
829     return;
830     }
831    
832 loizides 1.4 //--------------------------------------------------------------------------------------------------
833 loizides 1.1 int MultiVertexFitter::getIJKErr0() const
834     {
835     return _ctvmq.ijkerr[0];
836     }
837 loizides 1.4
838     //--------------------------------------------------------------------------------------------------
839 loizides 1.1 int MultiVertexFitter::getIJKErr1() const
840     {
841     return _ctvmq.ijkerr[1];
842     }
843 loizides 1.4
844     //--------------------------------------------------------------------------------------------------
845 loizides 1.1 int MultiVertexFitter::getIJKErr2() const
846     {
847     return _ctvmq.ijkerr[2];
848     }
849    
850 loizides 1.4 //--------------------------------------------------------------------------------------------------
851 loizides 1.1 int MultiVertexFitter::getErrTrackId() const
852     {
853     if (status() == 0) return 0;
854     int trkId = 0;
855     // Problems with track in CTVM01 or track has negative id in CTVM00 See PrintErr() for a more
856     // detailed list of error codes.
857     if ((_ctvmq.ijkerr[0] == 2 && _ctvmq.ijkerr[1] == 20) ||
858     _ctvmq.ijkerr[0] == 3) {
859     trkId = _ctvmq.list[_ctvmq.ijkerr[2]-1];
860     }
861    
862     return trkId;
863     }
864    
865 loizides 1.4 //--------------------------------------------------------------------------------------------------
866 loizides 1.1 string MultiVertexFitter::expert() const
867     {
868     return _expert;
869     }
870    
871 loizides 1.4 //--------------------------------------------------------------------------------------------------
872 loizides 1.1 int MultiVertexFitter::status() const
873     {
874     return _stat;
875     }
876    
877 loizides 1.4 //--------------------------------------------------------------------------------------------------
878 loizides 1.1 float MultiVertexFitter::chisq() const
879     {
880     // Chi-square of fit
881     return _ctvmq.chisqr[0];
882     }
883    
884 loizides 1.4 //--------------------------------------------------------------------------------------------------
885 loizides 1.1 int MultiVertexFitter::ndof() const
886     {
887     // Number of degrees of freedom of fit.
888     if (_ctvmq.chisqr[0] >= 0)
889     return _ctvmq.ndof;
890     else
891     return 0;
892     }
893    
894 loizides 1.4 //--------------------------------------------------------------------------------------------------
895 loizides 1.1 float MultiVertexFitter::prob() const
896     {
897     // Probability of chi-square of fit
898     if (_ctvmq.chisqr[0]>=0.) {
899     float chisq = _ctvmq.chisqr[0];
900     int nd = _ctvmq.ndof;
901     return TMath::Prob(chisq,nd);
902     }
903     else
904     return -999.;
905     }
906    
907 loizides 1.4 //--------------------------------------------------------------------------------------------------
908 loizides 1.1 float MultiVertexFitter::chisq(const int trkId) const
909     {
910     // This method returns the chisquare contribution for one track If fit not successful or not done
911     // yet, return -1.
912     if (_ctvmq.chisqr[0] < 0)
913     return -1.;
914     // Look for this track in the list of tracks.
915     for (int jt = 0; jt < _ctvmq.ntrack; ++jt) {
916     if (trkId == _ctvmq.list[jt]) {
917     // Found the track, so return its chisquare contribution.
918     return _ctvmq.chit[jt];
919     }
920     }
921     // If track is not in list, return -1.
922     return -1.;
923     }
924    
925 loizides 1.4 //--------------------------------------------------------------------------------------------------
926 loizides 1.1 float MultiVertexFitter::chisq_rphi() const
927     {
928     // This method returns the chisquare contribution in the r-phi plane.
929     int index[3] = {0,1,3};
930     // If fit not successful or not done yet, return -1.
931     if (_ctvmq.chisqr[0] < 0)
932     return -1.;
933     // Loop over the tracks in the event.
934     float chisq = 0.;
935     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
936     // Double loop over the relevant parameter indices.
937     for (int k1=0; k1<3; ++k1) {
938     for (int k2=0; k2<3; ++k2)
939     // Add contribution to chisquare.
940     chisq += _ctvmq.pardif[jt][index[k1]] *
941     _ctvmq.g[jt][index[k1]][index[k2]] *
942     _ctvmq.pardif[jt][index[k2]];
943     }
944     }
945     // Return the chisquare.
946     return chisq;
947     }
948    
949 loizides 1.4 //--------------------------------------------------------------------------------------------------
950 loizides 1.1 float MultiVertexFitter::chisq_z() const
951     {
952     // This method returns the chisquare contribution in the z direction.
953     int index[2] = {2,4};
954     // If fit not successful or not done yet, return -1.
955     if (_ctvmq.chisqr[0] < 0)
956     return -1.;
957     // Loop over the tracks in the event.
958     float chisq = 0.;
959     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
960     // Double loop over the relevant parameter indices.
961     for (int k1=0; k1<2; ++k1) {
962     for (int k2=0; k2<2; ++k2)
963     // Add contribution to chisquare.
964     chisq += _ctvmq.pardif[jt][index[k1]] *
965     _ctvmq.g[jt][index[k1]][index[k2]] *
966     _ctvmq.pardif[jt][index[k2]];
967     }
968     }
969     // Return the chisquare.
970     return chisq;
971     }
972    
973 loizides 1.4 //--------------------------------------------------------------------------------------------------
974 loizides 1.1 float MultiVertexFitter::chisq_rphiz() const
975     {
976     // This method returns the chisquare contribution of the cross
977     // terms in the r-phi and z directions.
978     int index1[2] = {2,4};
979     int index2[3] = {0,1,3};
980     // If fit not successful or not done yet, return -1.
981     if (_ctvmq.chisqr[0] < 0)
982     return -1.;
983     // Loop over the tracks in the event.
984     float chisq = 0.;
985     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
986     // Double loop over the relevant parameter indices.
987     for (int k1=0; k1<2; ++k1) {
988     for (int k2=0; k2<3; ++k2)
989     // Add contribution to chisquare.
990     chisq += _ctvmq.pardif[jt][index1[k1]] *
991     _ctvmq.g[jt][index1[k1]][index2[k2]] *
992     _ctvmq.pardif[jt][index2[k2]];
993     }
994     }
995    
996     // Return the chisquare.
997     return 2.0 * chisq;
998     }
999    
1000 loizides 1.4 //--------------------------------------------------------------------------------------------------
1001 loizides 1.1 float MultiVertexFitter::chisq_rphi(const int trkId) const
1002     {
1003     // This method returns the chisquare contribution in the r-phi plane.
1004     int index[3] = {0,1,3};
1005     // If fit not successful or not done yet, return -1.
1006     if (_ctvmq.chisqr[0] < 0)
1007     return -1.;
1008     // Loop over the tracks in the event, looking for the one we want
1009     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
1010     if (trkId == _ctvmq.list[jt]) {
1011     // Found the track, so calculate its chisquare contribution.
1012     float chisq = 0.;
1013     // Double loop over the relevant parameter indices.
1014     for (int k1=0; k1<3; ++k1) {
1015     for (int k2=0; k2<3; ++k2) {
1016     // Add contribution to chisquare.
1017     chisq += _ctvmq.pardif[jt][index[k1]] *
1018     _ctvmq.g[jt][index[k1]][index[k2]] *
1019     _ctvmq.pardif[jt][index[k2]];
1020     }
1021     }
1022     return chisq;
1023     }
1024     }
1025    
1026     // Track not found, return -1.
1027     return -1.;
1028     }
1029    
1030 loizides 1.4 //--------------------------------------------------------------------------------------------------
1031 loizides 1.1 float MultiVertexFitter::chisq_z(const int trkId) const
1032     {
1033     // This method returns the chisquare contribution in the z direction.
1034     int index[2] = {2,4};
1035     // If fit not successful or not done yet, return -1.
1036     if (_ctvmq.chisqr[0] < 0)
1037     return -1.;
1038     // Loop over the tracks in the event, looking for the one we want.
1039     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
1040     if (trkId == _ctvmq.list[jt]) {
1041     // Found the track, so calculate its chisquare contribution.
1042     float chisq = 0.;
1043     // Double loop over the relevant parameter indices.
1044     for (int k1=0; k1<2; ++k1) {
1045     for (int k2=0; k2<2; ++k2)
1046     // Add contribution to chisquare.
1047     chisq += _ctvmq.pardif[jt][index[k1]] *
1048     _ctvmq.g[jt][index[k1]][index[k2]] *
1049     _ctvmq.pardif[jt][index[k2]];
1050     }
1051     return chisq;
1052     }
1053     }
1054    
1055     // Track not found - return -1.
1056     return -1.;
1057     }
1058    
1059 loizides 1.4 //--------------------------------------------------------------------------------------------------
1060 loizides 1.1 float MultiVertexFitter::chisq_rphiz(const int trkId) const
1061     {
1062     // This method returns the chisquare contribution of the cross terms in the r-phi and z
1063     // directions
1064     int index1[2] = { 2,4 };
1065     int index2[3] = { 0,1,3 };
1066     // If fit not successful or not done yet, return -1
1067     if (_ctvmq.chisqr[0] < 0)
1068     return -1.;
1069     // Loop over the tracks in the event
1070     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
1071     if (trkId == _ctvmq.list[jt]) {
1072     // Found the track, so calculate its chisquare contribution
1073     float chisq = 0.;
1074     // Double loop over the relevant parameter indices
1075     for (int k1=0; k1<2; ++k1) {
1076     for (int k2=0; k2<3; ++k2)
1077     // Add contribution to chisquare
1078     chisq += _ctvmq.pardif[jt][index1[k1]] *
1079     _ctvmq.g[jt][index1[k1]][index2[k2]] *
1080     _ctvmq.pardif[jt][index2[k2]];
1081     }
1082     return 2.0 * chisq;
1083     }
1084     }
1085     // Track not found so return -1.
1086     return -1.;
1087     }
1088    
1089 loizides 1.4 //--------------------------------------------------------------------------------------------------
1090 loizides 1.1 FourVector MultiVertexFitter::getTrackP4(const int trkId) const
1091     {
1092     if (_stat != 0)
1093     return FourVector(0,0,0,0);
1094     // return four momentum of fit track
1095     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
1096     // Find which track matches this Id.
1097     if (trkId == _ctvmq.list[jt]) {
1098     FourVector p4((double)_ctvmq.trkp4[0][jt], (double)_ctvmq.trkp4[1][jt],
1099     (double)_ctvmq.trkp4[2][jt], (double)_ctvmq.trkp4[3][jt]);
1100     return p4;
1101     }
1102     }
1103     return FourVector(0,0,0,0);
1104     }
1105    
1106 loizides 1.4 //--------------------------------------------------------------------------------------------------
1107 loizides 1.1 float MultiVertexFitter::getMass(int ntrk, const int trkIds[], float &dmass) const
1108     {
1109     // #if (defined(LINUX) && defined(__USE_BSD)) || defined(OSF1)
1110     // struct sigaction myaction = {MultiVertexFitterSetStatus, 0, 0, 0}, oldaction;
1111     // sigaction(SIGFPE, &myaction, &oldaction);
1112     // if (setjmp(env)!=0) {
1113     // sigaction(SIGFPE, &oldaction,0);
1114     // return -999;
1115     // }
1116     // #endif
1117    
1118     dmass = -999.;
1119     if (_stat!=0)
1120     return -999.;
1121     // Get fit invariant mass of ntrk tracks listed in array tracks. dmass is the error on the mass.
1122     dmass=-999.;
1123     if (ntrk <= 0)
1124     return 0;
1125     int jtrks[_maxtrk];
1126     for (int jt=0; jt<ntrk; ++jt) {
1127     bool found = false;
1128     for (int kt=0; kt<_ctvmq.ntrack; ++kt) {
1129     if (trkIds[jt] == _ctvmq.list[kt]) {
1130     found = true;
1131     jtrks[jt]=kt+1;
1132     }
1133     }
1134     if (!found)
1135     return 0;
1136     }
1137     // Copy information into CTVMFT common blocks
1138     *_ctvmq_com = _ctvmq;
1139     *_ctvmfr_com = _ctvmfr;
1140     int ntr = ntrk;
1141     float mass;
1142     double p4[4];
1143     mcalc_(ntr, jtrks, mass, dmass, p4);
1144    
1145     // #if (defined(LINUX) && defined(__USE_BSD)) || defined(OSF1)
1146     // sigaction(SIGFPE, &oldaction,0);
1147     // #endif
1148    
1149     return mass;
1150     }
1151    
1152 loizides 1.4 //--------------------------------------------------------------------------------------------------
1153 loizides 1.1 float MultiVertexFitter::getDecayLength(vertexNumber nv, vertexNumber mv,
1154     const Hep3Vector& dir, float& dlerr) const
1155     {
1156     dlerr = -999.;
1157     if (_stat!=0)
1158     return -999.;
1159    
1160     // Get the signed decay length from vertex nv to vertex mv along the x-y direction defined by the
1161     // 3 vector dir. dlerr is the error on the decay length including the full error matrix. Check
1162     // that vertices are within range.
1163     if (nv<0 || nv>=_ctvmq.nvertx)
1164     return -999.;
1165     if (mv<1 || mv>_ctvmq.nvertx)
1166     return -999.;
1167     if (nv>=mv)
1168     return -999.;
1169    
1170     float dir_t = sqrt(dir.x()*dir.x()+dir.y()*dir.y());
1171     if (dir_t <= 0.)
1172     return -999.;
1173    
1174     Hep3Vector dv = getVertexHep(mv)-getVertexHep(nv);
1175     float dl = (dv.x()*dir.x()+dv.y()*dir.y())/dir_t;
1176     // Set up the column matrix of derivatives
1177     HepMatrix A(4,1);
1178     A(1,1) = dir.x()/dir_t;
1179     A(2,1) = dir.y()/dir_t;
1180     A(3,1) = -dir.x()/dir_t;
1181     A(4,1) = -dir.y()/dir_t;
1182     // Check if first vertex (nv) is primary vertex. If it is, check if it was used in the primary
1183     // vertex. If not, all of the corresponding error matrix elements are those supplied for the
1184     // primary vertex.
1185     int nvf = 0;
1186     if (nv==0) {
1187     nvf=-1;
1188     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
1189     if (_ctvmq.vtxpnt[0][jv]==0)
1190     nvf=0;
1191     }
1192     }
1193     // Get the relevant submatrix of the full error matrix.
1194     HepMatrix V(4,4,0);
1195     if (nvf < 0) {
1196     V(1,1) = getErrorMatrixHep(_ctvmq.voff[mv-1]+1,_ctvmq.voff[mv-1]+1);
1197     V(1,2) = getErrorMatrixHep(_ctvmq.voff[mv-1]+1,_ctvmq.voff[mv-1]+2);
1198     V(2,1) = getErrorMatrixHep(_ctvmq.voff[mv-1]+2,_ctvmq.voff[mv-1]+1);
1199     V(2,2) = getErrorMatrixHep(_ctvmq.voff[mv-1]+2,_ctvmq.voff[mv-1]+2);
1200     V(3,3) = _ctvmq.exyzpv[0][0];
1201     V(3,4) = _ctvmq.exyzpv[0][1];
1202     V(4,3) = _ctvmq.exyzpv[1][0];
1203     V(4,4) = _ctvmq.exyzpv[1][1];
1204     }
1205     else {
1206     // Get the indices into the error matrix vmat
1207     int index[4] = { _ctvmq.voff[mv-1]+1,_ctvmq.voff[mv-1]+2,0,0 };
1208     if (nv==0) {
1209     index[2] = 1;
1210     index[3] = 2;
1211     }
1212     else {
1213     index[2] = _ctvmq.voff[nv-1]+1;
1214     index[3] = _ctvmq.voff[nv-1]+2;
1215     }
1216     for (int j=0; j<4; ++j) {
1217     for (int k=0; k<4; ++k) {
1218     V[j][k] = getErrorMatrixHep(index[j],index[k]);
1219     }
1220     }
1221     }
1222    
1223     // Calculate square of dlerr
1224     dlerr = (A.T()*V*A)(1,1);
1225     if (dlerr >= 0.)
1226     dlerr = sqrt(dlerr);
1227     else
1228     dlerr = -sqrt(-dlerr);
1229    
1230     return dl;
1231     }
1232    
1233 loizides 1.4 //--------------------------------------------------------------------------------------------------
1234 bendavid 1.2 float MultiVertexFitter::getDecayLength(vertexNumber nv, vertexNumber mv,
1235     const ThreeVector& dir, float& dlerr) const
1236     {
1237     Hep3Vector dirHep(dir.x(),dir.y(),dir.z());
1238     return getDecayLength(nv, mv, dirHep, dlerr);
1239     }
1240    
1241 loizides 1.4 //--------------------------------------------------------------------------------------------------
1242 bendavid 1.2 float MultiVertexFitter::getZDecayLength(vertexNumber nv, vertexNumber mv,
1243     const Hep3Vector& mom, float& dlerr) const
1244     {
1245     //----------------------------------------------------------------------------
1246     // Get the signed decay length from vertex nv to vertex mv along the
1247     // z direction of the momentum vector, mom.
1248     // dlerr is the error on the decay length including the full error
1249     // matrix.
1250     //----------------------------------------------------------------------------
1251     // Start with good initialization
1252     dlerr = -999.;
1253    
1254     // Check that the fit worked
1255     if (_stat != 0)
1256     return -999.;
1257    
1258     // Check that vertices are within range.
1259     if (nv<0 || nv>=_ctvmq.nvertx)
1260     return -999.;
1261     if (mv<1 || mv> _ctvmq.nvertx)
1262     return -999.;
1263     if (nv >= mv)
1264     return -999.;
1265    
1266     // Calculate the vector length
1267     float length = fabs(mom.z());
1268     if (length <= 0.)
1269     return -999.;
1270    
1271     // Get the vector pointing from first vertex (nv) to second vertex (mv)
1272     Hep3Vector dv = getVertexHep(mv) - getVertexHep(nv);
1273    
1274     //----------------------------------------------------------------------------
1275     // Calculate the "decay distance"
1276     //----------------------------------------------------------------------------
1277     // Project the vertex vector onto the momentum vector direction
1278     float dl = (dv.z()*mom.z())/length;
1279    
1280     //----------------------------------------------------------------------------
1281     // Calculate the error on that distance
1282     //----------------------------------------------------------------------------
1283     // Set up the column matrix of derivatives
1284     HepMatrix A(2,1);
1285     A(1,1) = mom.z()/length;
1286     A(2,1) = -mom.z()/length;
1287    
1288     // Need to catch the special case if the first vertex is the primary
1289     int nvf = 0;
1290     if (nv == 0) {
1291     nvf = -1;
1292     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
1293     if (_ctvmq.vtxpnt[0][jv] == 0)
1294     nvf = 0;
1295     }
1296     }
1297     // Get the relevant submatrix of the full error matrix.
1298     HepMatrix V(2,2,0);
1299     if (nvf < 0) {
1300     // Geometric uncertainties (positions second vertex)
1301     V(1,1) = getErrorMatrixHep(_ctvmq.voff[mv-1]+3,_ctvmq.voff[mv-1]+3);
1302     // Geometric uncertainties (positions first vertex)
1303     V(2,2) = _ctvmq.exyzpv[2][2];
1304     }
1305     else {
1306     // Get the indices into the error matrix vmat
1307     int index[2] = { _ctvmq.voff[mv-1]+3,_ctvmq.voff[nv-1]+3 };
1308     // Handeling the case of the primary vertex
1309     if (nv == 0)
1310     index[1] = 3;
1311     // All right... copy
1312     for(int j=0; j<2; ++j)
1313     for(int k=0; k<2; ++k)
1314     V[j][k] = getErrorMatrixHep(index[j],index[k]);
1315     }
1316     // Calculate square of dlerr
1317     dlerr = (A.T() * V * A )(1,1);
1318     if (dlerr >= 0.)
1319     dlerr = sqrt(dlerr);
1320     else
1321     dlerr = -sqrt(-dlerr);
1322    
1323     return dl;
1324     }
1325    
1326 loizides 1.4 //--------------------------------------------------------------------------------------------------
1327 bendavid 1.2 float MultiVertexFitter::getZDecayLength(vertexNumber nv, vertexNumber mv,
1328     const ThreeVector& mom, float& dlerr) const
1329     {
1330     Hep3Vector momHep(mom.x(),mom.y(),mom.z());
1331     return getZDecayLength(nv, mv, momHep, dlerr);
1332     }
1333    
1334 loizides 1.4 //--------------------------------------------------------------------------------------------------
1335 bendavid 1.2 float MultiVertexFitter::getImpactPar(vertexNumber prdV, vertexNumber dcyV,
1336     const Hep3Vector &v, float &dxyerr) const
1337     {
1338     Hep3Vector PVtx = getVertexHep (prdV);
1339     Hep3Vector DVtx = getVertexHep (dcyV);
1340     HepSymMatrix PVtxCv = getErrorMatrixHep(prdV);
1341     HepSymMatrix DVtxCv = getErrorMatrixHep(dcyV);
1342    
1343     double norma = v.perp();
1344     if (norma <= 0) {
1345     dxyerr = -999.;
1346     return -999.0;
1347     }
1348     double dxy = ((v.cross(DVtx-PVtx)).z())/norma;
1349    
1350     // Calculate error on B impact parameter:
1351     double cosPhi = cos(v.phi());
1352     double sinPhi = sin(v.phi());
1353     dxyerr = cosPhi * cosPhi * (DVtxCv[1][1] + PVtxCv[1][1])
1354     + sinPhi * sinPhi * (DVtxCv[0][0] + PVtxCv[0][0])
1355     - 2.0 * cosPhi * sinPhi * (DVtxCv[0][1] + PVtxCv[0][1]);
1356     dxyerr = (dxyerr>0.0) ? sqrt(dxyerr) : -999.;
1357    
1358     return dxy;
1359     }
1360    
1361 loizides 1.4 //--------------------------------------------------------------------------------------------------
1362 bendavid 1.2 float MultiVertexFitter::getImpactPar(vertexNumber prdV, vertexNumber dcyV,
1363     const ThreeVector &v, float &dxyerr) const
1364     {
1365     Hep3Vector vHep(v.x(),v.y(),v.z());
1366     return getImpactPar(prdV, dcyV, vHep, dxyerr);
1367     }
1368    
1369 loizides 1.4 //--------------------------------------------------------------------------------------------------
1370 loizides 1.1 float MultiVertexFitter::get_dr(vertexNumber nv, vertexNumber mv, float& drerr) const
1371     {
1372     drerr = -999.;
1373     if (_stat!=0)
1374     return -999.;
1375     // Get the transvese distance between vertices nv and mv and return it as the function value.
1376     // drerr is the uncertainty on the transverse distance, calculated from the full error matrix
1377     // including correlations.
1378     float dxyz[3];
1379     float dr;
1380     float dz;
1381     float dl[3];
1382    
1383     int mvert = mv;
1384     int nvert = nv;
1385     // Copy information into CTVMFT common blocks
1386     *_ctvmq_com = _ctvmq;
1387     *_ctvmfr_com = _ctvmfr;
1388     // Do calculation
1389     dcalc_(mvert,nvert,dxyz,dr,dz,dl);
1390     drerr = dl[0];
1391     return dr;
1392     }
1393    
1394 loizides 1.4 //--------------------------------------------------------------------------------------------------
1395 loizides 1.1 float MultiVertexFitter::get_dz(vertexNumber nv, vertexNumber mv, float& dzerr) const
1396     {
1397     dzerr = -999.;
1398     if (_stat!=0)
1399     return -999.;
1400     // Get the longitudinal distance between vertices nv and mv and return it as the function value.
1401     // dzerr is the uncertainty on the longitudinal distance, calculated from the full error matrix
1402     // including correlations.
1403     float dxyz[3];
1404     float dr;
1405     float dz;
1406     float dl[3];
1407    
1408     int mvert = mv;
1409     int nvert = nv;
1410     // Copy information into CTVMFT common blocks
1411     *_ctvmq_com = _ctvmq;
1412     *_ctvmfr_com = _ctvmfr;
1413    
1414     // Do calculation
1415     dcalc_(mvert,nvert,dxyz,dr,dz,dl);
1416     dzerr = dl[1];
1417     return dz;
1418     }
1419    
1420 loizides 1.4 //--------------------------------------------------------------------------------------------------
1421 loizides 1.1 Hep3Vector MultiVertexFitter::getVertexHep(vertexNumber nv) const
1422     {
1423     if (_stat!=0)
1424     return Hep3Vector(-999,-999,-999);
1425     // Return x,y,z position of vertex nv.
1426     if (nv<0 || nv>_ctvmq.nvertx)
1427     return Hep3Vector(0,0,0);
1428     // Check if first vertex (nv) is primary vertex. If it is, check if it was used in the primary
1429     // vertex. If not, all of the corresponding error matrix elements are those supplied for the
1430     // primary vertex.
1431     int nvf=0;
1432     if (nv==0) {
1433     nvf=-1;
1434     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
1435     if (_ctvmq.vtxpnt[0][jv]==0)
1436     nvf=0;
1437     }
1438     }
1439     Hep3Vector vertex;
1440     // If primary vertex was not part of fit, take vertex location as supplied.
1441     if (nvf < 0) {
1442     vertex.setX(_ctvmq.xyzpv0[0]);
1443     vertex.setY(_ctvmq.xyzpv0[1]);
1444     vertex.setZ(_ctvmq.xyzpv0[2]);
1445     }
1446     else {
1447     vertex.setX(_ctvmq.xyzvrt[nv][0]);
1448     vertex.setY(_ctvmq.xyzvrt[nv][1]);
1449     vertex.setZ(_ctvmq.xyzvrt[nv][2]);
1450     }
1451     //// If we have a different reference point, need to add it back in
1452     //vertex += _referencePoint;
1453     return vertex;
1454     }
1455    
1456 loizides 1.4 //--------------------------------------------------------------------------------------------------
1457 loizides 1.1 ThreeVector MultiVertexFitter::getVertex(vertexNumber nv) const
1458     {
1459     if (_stat!=0)
1460     return ThreeVector(-999,-999,-999);
1461     // Return x,y,z position of vertex nv.
1462     if (nv<0 || nv>_ctvmq.nvertx)
1463     return ThreeVector(0,0,0);
1464     // Check if first vertex (nv) is primary vertex. If it is, check if it was used in the primary
1465     // vertex. If not, all of the corresponding error matrix elements are those supplied for the
1466     // primary vertex.
1467     int nvf=0;
1468     if (nv==0) {
1469     nvf=-1;
1470     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
1471     if (_ctvmq.vtxpnt[0][jv]==0)
1472     nvf=0;
1473     }
1474     }
1475     ThreeVector vertex;
1476     // If primary vertex was not part of fit, take vertex location as supplied.
1477     if (nvf < 0) {
1478     vertex.SetX(_ctvmq.xyzpv0[0]);
1479     vertex.SetY(_ctvmq.xyzpv0[1]);
1480     vertex.SetZ(_ctvmq.xyzpv0[2]);
1481     }
1482     else {
1483     vertex.SetX(_ctvmq.xyzvrt[nv][0]);
1484     vertex.SetY(_ctvmq.xyzvrt[nv][1]);
1485     vertex.SetZ(_ctvmq.xyzvrt[nv][2]);
1486     }
1487     // If we have a different reference point, need to add it back in
1488     vertex += _referencePoint;
1489     return vertex;
1490     }
1491    
1492 loizides 1.4 //--------------------------------------------------------------------------------------------------
1493 loizides 1.1 ThreeSymMatrix MultiVertexFitter::getErrorMatrix(MultiVertexFitter::vertexNumber nv) const
1494     {
1495     // return errors for vertex nv
1496     ThreeSymMatrix cov;
1497    
1498     // if this is the primary vertex, return the error matrix the user supplied
1499     if (nv==PRIMARY_VERTEX) {
1500     for (int j=0; j<3; j++)
1501     for (int k=0; k<3; k++)
1502     cov(j,k) = _ctvmq.exyzpv[j][k];
1503     return cov;
1504     }
1505    
1506     if (_stat!=0)
1507     return cov;
1508     // return x,y,z position of vertex nv
1509     if (nv<VERTEX_1 || nv>_ctvmq.nvertx)
1510     return cov;
1511     // get offset for vertex nv
1512     int voff = _ctvmq.voff[nv-1];
1513     // fill matrix
1514     for (int i = 0 ; i < 3 ; ++i)
1515     for (int j = i ; j < 3 ; ++j)
1516     cov(i,j) = _ctvmfr.vmat[voff+i][voff+j];
1517     return cov;
1518     }
1519    
1520 loizides 1.4 //--------------------------------------------------------------------------------------------------
1521 loizides 1.1 double MultiVertexFitter::getErrorMatrixHep(int j, int k) const
1522     {
1523     if (_stat!=0)
1524     return -999.;
1525     // Return the j,k element of the full error matrix VMAT. Since the CTVMFT documentation assumes
1526     // the indices start from 1 (ala Fortran), we will also assume this and convert the C++ indices.
1527     // Note also the that order of Fortran and C++ indices is different. We assume that j and k are
1528     // in the Fortran order.
1529     if (j<1 || k<1 || j>_maxdim+1 || k>_maxdim)
1530     return -999.;
1531    
1532     return _ctvmfr.vmat[k-1][j-1];
1533     }
1534    
1535 loizides 1.4 //--------------------------------------------------------------------------------------------------
1536 loizides 1.1 HepSymMatrix MultiVertexFitter::getErrorMatrixHep(MultiVertexFitter::vertexNumber nv) const
1537     {
1538     // return errors for vertex nv
1539     HepSymMatrix cov(3,0);
1540    
1541     // if this is the primary vertex, return the error matrix the user supplied
1542     if (nv==PRIMARY_VERTEX) {
1543     for (int j=0; j<3; j++)
1544     for (int k=0; k<3; k++)
1545     cov[j][k] = _ctvmq.exyzpv[j][k];
1546     return cov;
1547     }
1548    
1549     if (_stat!=0)
1550     return cov;
1551     // return x,y,z position of vertex nv
1552     if (nv<VERTEX_1 || nv>_ctvmq.nvertx)
1553     return cov;
1554     // get offset for vertex nv
1555     int voff = _ctvmq.voff[nv-1];
1556     // fill matrix
1557     for (int i = 0 ; i < 3 ; ++i)
1558     for (int j = i ; j < 3 ; ++j)
1559     cov[i][j] = _ctvmfr.vmat[voff+i][voff+j];
1560     return cov;
1561     }
1562    
1563 loizides 1.4 //--------------------------------------------------------------------------------------------------
1564 loizides 1.1 HepSymMatrix MultiVertexFitter::getErrorMatrixHep(const int trkId) const
1565     {
1566     HepSymMatrix cov(3,0);
1567     if (_stat != 0)
1568     return cov;
1569    
1570     for (int nt=0; nt<_ctvmq.ntrack; ++nt) {
1571    
1572     // Find which track matches this Id
1573     if (trkId == _ctvmq.list[nt]) {
1574    
1575     // Position of track nt get offset for track nt
1576     int toff = _ctvmq.toff[nt];
1577    
1578     // Fill matrix -- Crv,Phi,Ctg
1579     for (int i = 0 ; i < 3 ; ++i)
1580     for (int j = i ; j < 3 ; ++j)
1581     cov[i][j] = _ctvmfr.vmat[toff+i][toff+j];
1582     }
1583     }
1584     return cov;
1585     }
1586    
1587 loizides 1.4 //--------------------------------------------------------------------------------------------------
1588 loizides 1.1 float MultiVertexFitter::getPtError(const int trkId) const
1589     {
1590     if (_stat != 0)
1591     return 0;
1592    
1593     int toff;
1594     float pt,curv,curvErr,ptErr;
1595    
1596     for (int nt=0; nt<_ctvmq.ntrack; ++nt) {
1597    
1598     // Find which track matches this Id
1599     if (trkId == _ctvmq.list[nt]) {
1600    
1601     // Position of track nt get offset for track nt
1602     toff = _ctvmq.toff[nt];
1603    
1604     // Curvature error
1605     pt = sqrt(_ctvmq.trkp4[0][nt]*_ctvmq.trkp4[0][nt] +
1606     _ctvmq.trkp4[1][nt]*_ctvmq.trkp4[1][nt]);
1607     curv = _ctvmq.pscale/pt;
1608     curvErr = sqrt(_ctvmfr.vmat[toff+0][toff+0]);
1609     ptErr = _ctvmq.pscale/curv/curv*curvErr;
1610     return ptErr;
1611     }
1612     }
1613    
1614     return 0;
1615     }
1616    
1617 loizides 1.4 //--------------------------------------------------------------------------------------------------
1618 loizides 1.1 void MultiVertexFitter::getPosMomErr(HepMatrix& errors) const
1619     {
1620     // A c++ rewrite of the FORTRAN MASSAG function The result of this function is an error matrix in
1621     // position-momentum basis. A 7x7 matrix of errors where the rows/columns are x, y, z, px, py,
1622     // pz, e.
1623    
1624     double cosph [_maxtrk];
1625     double sinph [_maxtrk];
1626     double cosdph[_maxtrk];
1627     double sindph[_maxtrk];
1628     Hep3Vector mom3 [_maxtrk];
1629     HepLorentzVector pmom [_maxtrk];
1630     HepLorentzVector total_mom;
1631    
1632     for (int lvtx = 0; lvtx < _ctvmq.nvertx; lvtx++) {
1633     for (int ltrk = 0; ltrk < _ctvmq.ntrack; ltrk++) {
1634     if (!_ctvmq.trkvtx[lvtx][ltrk])
1635     continue;
1636     cosph[ltrk] = cos(_ctvmq.par[ltrk][1]);
1637     sinph[ltrk] = sin(_ctvmq.par[ltrk][1]);
1638     double dphi = 0;
1639     sindph[ltrk] = 2 * _ctvmq.par[ltrk][0] *
1640     (_ctvmq.xyzvrt[lvtx + 1][0] * cosph[ltrk] +
1641     _ctvmq.xyzvrt[lvtx + 1][1] * sinph[ltrk]);
1642     cosdph[ltrk] = sqrt(1.0 - sindph[ltrk] * sindph[ltrk]);
1643     if (fabs(sindph[ltrk]) <= 1.0){
1644     dphi = asin(sindph[ltrk]);
1645     }
1646     double pt = _ctvmq.pscale * fabs(1./_ctvmq.par[ltrk][0]);
1647     mom3[ltrk].setX(pt * cos(_ctvmq.par[ltrk][1] + dphi));
1648     mom3[ltrk].setY(pt * sin(_ctvmq.par[ltrk][1] + dphi));
1649     mom3[ltrk].setZ(pt * _ctvmq.par[ltrk][2]);
1650     double e = sqrt(_ctvmq.tmass[ltrk] * _ctvmq.tmass[ltrk]
1651     + mom3[ltrk].mag2());
1652     pmom[ltrk].setVect(mom3[ltrk]);
1653     pmom[ltrk].setE(e);
1654    
1655     total_mom += pmom[ltrk];
1656     }
1657     }
1658    
1659     // Easy so far, but now it gets ugly: fill dp_dpar with the derivatives of the position and
1660     // momentum with respect to the parameters
1661    
1662     int ctvmft_dim = 3 * (_ctvmq.nvertx + _ctvmq.ntrack);
1663     HepMatrix deriv(ctvmft_dim, 7, 0);
1664 paus 1.3 //HepMatrix dp_dpar[_maxvtx] = { deriv, deriv, deriv };
1665     HepMatrix dp_dpar[_maxvtx] = { deriv };
1666 loizides 1.1
1667     // Fill the x, y, z rows:
1668     for (int nvtx = 0; nvtx < _ctvmq.nvertx; nvtx++) {
1669     for (int lcomp = 0; lcomp < 3; lcomp++){
1670     dp_dpar[nvtx][(3 * nvtx) + lcomp][lcomp] = 1.0;
1671     }
1672    
1673     // Fill the px, py, pz, e rows:
1674     for (int lvtx = 0; lvtx < _ctvmq.nvertx; lvtx++){
1675     for (int ltrk = 0; ltrk < _ctvmq.ntrack; ltrk++){
1676     if (!_ctvmq.trkvtx[lvtx][ltrk])
1677     continue;
1678    
1679     // Find the derivatives of dphi with respect to x, y, curvature, and phi0:
1680     double dphi_dx = 2.0 * _ctvmq.par[ltrk][0] * cosph[ltrk]/cosdph[ltrk];
1681     double dphi_dy = 2.0 * _ctvmq.par[ltrk][0] * sinph[ltrk]/cosdph[ltrk];
1682     double dphi_dc = 2.0 *
1683     (_ctvmq.xyzvrt[lvtx + 1][0] * cosph[ltrk] +
1684     _ctvmq.xyzvrt[lvtx + 1][1] * sinph[ltrk])/cosdph[ltrk];
1685     double dphi_dp = 2.0 * _ctvmq.par[ltrk][0] *
1686     (-_ctvmq.xyzvrt[lvtx + 1][0] * sinph[ltrk] +
1687     _ctvmq.xyzvrt[lvtx + 1][1] * cosph[ltrk])/cosdph[ltrk];
1688    
1689     // Now load the derivative matrix
1690     int lvele = 3 * lvtx;
1691     // dPx/dx:
1692     dp_dpar[nvtx][lvele][3] += -pmom[ltrk].y() * dphi_dx;
1693     // dPy/dx:
1694     dp_dpar[nvtx][lvele][4] += pmom[ltrk].x() * dphi_dx;
1695     // dPz/dx:
1696     dp_dpar[nvtx][lvele][5] = 0.;
1697     // dPy/dx:
1698     dp_dpar[nvtx][lvele][6] = 0.;
1699    
1700     lvele++;
1701     // dPx/dy:
1702     dp_dpar[nvtx][lvele][3] += -pmom[ltrk].y() * dphi_dy;
1703     // dPy/dy:
1704     dp_dpar[nvtx][lvele][4] += pmom[ltrk].x() * dphi_dy;
1705     // dPz/dy:
1706     dp_dpar[nvtx][lvele][5] = 0.;
1707     // dE/dy:
1708     dp_dpar[nvtx][lvele][6] = 0.;
1709    
1710     lvele++;
1711     // dPx/dz:
1712     dp_dpar[nvtx][lvele][3] = 0.;
1713     // dPy/dz:
1714     dp_dpar[nvtx][lvele][4] = 0.;
1715     // dPz/dz:
1716     dp_dpar[nvtx][lvele][5] = 0.;
1717     // dE/dz:
1718     dp_dpar[nvtx][lvele][6] = 0.;
1719    
1720     lvele = 3 * (ltrk + _ctvmq.nvertx);
1721     // dPx/dcurv[ltrk]:
1722     dp_dpar[nvtx][lvele][3] = -(pmom[ltrk].x()/_ctvmq.par[ltrk][0])
1723     - pmom[ltrk].y() * dphi_dc;
1724     // dPy/dcurv[ltrk]:
1725     dp_dpar[nvtx][lvele][4] = -(pmom[ltrk].y()/_ctvmq.par[ltrk][0])
1726     + pmom[ltrk].x() * dphi_dc;
1727     // dPz/dcurv[ltrk]:
1728     dp_dpar[nvtx][lvele][5] = -(pmom[ltrk].z()/_ctvmq.par[ltrk][0]);
1729     // dE/dcurv[ltrk]:
1730     dp_dpar[nvtx][lvele][6] =
1731     -mom3[ltrk].mag2()/(_ctvmq.par[ltrk][0] * pmom[ltrk].e());
1732    
1733     lvele++;
1734     // dPx/dphi[ltrk]:
1735     dp_dpar[nvtx][lvele][3] = -pmom[ltrk].y() * (1.0 + dphi_dp);
1736     // dPy/dphi[ltrk]:
1737     dp_dpar[nvtx][lvele][4] = pmom[ltrk].x() * (1.0 + dphi_dp);
1738     // dPz/dphi[ltrk]:
1739     dp_dpar[nvtx][lvele][5] = 0.;
1740     // dE/dphi[ltrk]:
1741     dp_dpar[nvtx][lvele][6] = 0.;
1742    
1743     lvele++;
1744     // dPx/dcot[ltrk]:
1745     dp_dpar[nvtx][lvele][3] = 0;
1746     // dPy/dcot[ltrk]:
1747     dp_dpar[nvtx][lvele][4] = 0;
1748     // dPz/dcot[ltrk]:
1749     dp_dpar[nvtx][lvele][5] = pmom[ltrk].perp();
1750     // dE/dcot[ltrk]:
1751     dp_dpar[nvtx][lvele][6] =
1752     pmom[ltrk].perp2() * _ctvmq.par[ltrk][2] / pmom[ltrk].e();
1753     }
1754     }
1755     }
1756    
1757     // -----------------------------------------------------------------------------------------------
1758     // Now calculate the new error matrix
1759     // -----------------------------------------------------------------------------------------------
1760     // Extract the interesting bits from VMAT
1761     HepMatrix vmat(ctvmft_dim,ctvmft_dim,0);
1762     for (int lpar = 0; lpar < ctvmft_dim; lpar++){
1763     int l = lpar%3;
1764     int lvele = 0;
1765     if (lpar < 3 * _ctvmq.nvertx){
1766     int lvtx = lpar/3;
1767     lvele = _ctvmq.voff[lvtx] + l;
1768     }
1769     else {
1770     int ltrk = (lpar - 3 * _ctvmq.nvertx)/3;
1771     lvele = _ctvmq.toff[ltrk] + l;
1772     }
1773     for (int kpar = 0; kpar < ctvmft_dim; kpar++) {
1774     int k = kpar%3;
1775     int kvele = 0;
1776     if (kpar < 3 * _ctvmq.nvertx) {
1777     int kvtx = kpar/3;
1778     kvele = _ctvmq.voff[kvtx] + k;
1779     }
1780     else{
1781     int ktrk = (kpar - 3 * _ctvmq.nvertx)/3;
1782     kvele = _ctvmq.toff[ktrk] + k;
1783     }
1784     vmat[kpar][lpar] = _ctvmfr.vmat[kvele][lvele];
1785     }
1786     }
1787    
1788     // Just about done
1789     HepMatrix ans(7,7,0);
1790 paus 1.3 //HepMatrix answer[_maxvtx] = { ans, ans, ans };
1791     HepMatrix answer[_maxvtx] = { ans };
1792 loizides 1.1 for (int nvtx = 0; nvtx < _ctvmq.nvertx; nvtx++)
1793     answer[nvtx] = (dp_dpar[nvtx].T() * vmat) * dp_dpar[nvtx];
1794     errors = answer[0];
1795     }
1796    
1797 loizides 1.4 //--------------------------------------------------------------------------------------------------
1798 loizides 1.1 int MultiVertexFitter::vOff(vertexNumber jv) const
1799     {
1800     if (jv < VERTEX_1 || jv > _maxvtx)
1801     return -999;
1802     else
1803     return _ctvmq.voff[jv-1];
1804     }
1805    
1806 loizides 1.4 //--------------------------------------------------------------------------------------------------
1807 loizides 1.1 int MultiVertexFitter::tOff(const int trkId) const
1808     {
1809     for (int kt=0; kt<_ctvmq.ntrack; ++kt) {
1810     if (trkId == _ctvmq.list[kt])
1811     return _ctvmq.toff[kt];
1812     }
1813     return -999;
1814     }
1815    
1816 loizides 1.4 //--------------------------------------------------------------------------------------------------
1817 loizides 1.1 int MultiVertexFitter::pOff(int lp) const
1818     {
1819     if (lp < 1)
1820     return -999;
1821     else
1822     return _ctvmq.poff[lp-1];
1823     }
1824    
1825 loizides 1.4 //--------------------------------------------------------------------------------------------------
1826 loizides 1.1 int MultiVertexFitter::cOff(int lc) const
1827     {
1828     if (lc < 1)
1829     return -999;
1830     else
1831     return _ctvmq.coff[lc-1];
1832     }
1833    
1834 loizides 1.4 //--------------------------------------------------------------------------------------------------
1835 loizides 1.1 int MultiVertexFitter::mOff() const
1836     {
1837     return _ctvmq.moff;
1838     }
1839    
1840 loizides 1.4 //--------------------------------------------------------------------------------------------------
1841 loizides 1.1 double MultiVertexFitter::VMat(int i, int j) const
1842     {
1843     if (i <0 || j < 0)
1844     return -999.;
1845     else
1846     return _ctvmfr.vmat[i][j];
1847     }
1848    
1849     // Facilitates CandNode recursion. CandNodes have no way of deciding which vertex they are, and
1850     // these trivial functions help them do that.
1851 loizides 1.4 //--------------------------------------------------------------------------------------------------
1852 loizides 1.1 MultiVertexFitter::vertexNumber MultiVertexFitter::allocateVertexNumber()
1853     {
1854     if ((_currentAllocatedVertexNumber < PRIMARY_VERTEX) ||
1855     (_currentAllocatedVertexNumber > _maxvtx)) {
1856     cout << "MultiVertexFitter::allocateVertexNumber: out of range!" << endl;
1857     return PRIMARY_VERTEX;
1858     }
1859     return vertexNumber(++_currentAllocatedVertexNumber);
1860     }
1861    
1862 loizides 1.4 //--------------------------------------------------------------------------------------------------
1863 loizides 1.1 void MultiVertexFitter::resetAllocatedVertexNumber()
1864     {
1865     _currentAllocatedVertexNumber = 0;
1866     }
1867    
1868 loizides 1.4 //--------------------------------------------------------------------------------------------------
1869 loizides 1.1 void MultiVertexFitter::restoreFromCommons()
1870     {
1871     _stat = 0;
1872     _ctvmq_com = (CTVMQ*) ctvmq_address_();
1873     _ctvmfr_com = (CTVMFR*) ctvmfr_address_();
1874     _fiddle_com = (FIDDLE*) fiddle_address_();
1875     _trkprm_com = (TRKPRM*) trkprm_address_();
1876     _ctvmq = *_ctvmq_com;
1877     _ctvmfr = *_ctvmfr_com;
1878     _fiddle = *_fiddle_com;
1879     _trkprm = *_trkprm_com;
1880     }