ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitCommon/VertexFit/src/MultiVertexFitterD.cc
Revision: 1.1
Committed: Thu Nov 13 16:34:29 2008 UTC (16 years, 5 months ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_008pre2, Mit_008pre1, Mit_006b, Mit_006a, Mit_006
Log Message:
Submitting the ugly three versions. But no elegant solution found.

File Contents

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