ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitCommon/VertexFit/src/MultiVertexFitterD.cc
Revision: 1.3
Committed: Mon Oct 12 21:39:22 2009 UTC (15 years, 6 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_025c_branch2, Mit_025c_branch1, Mit_030, Mit_029c, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_025c_branch0, Mit_028, Mit_027a, Mit_027, Mit_026, Mit_025e, Mit_025d, Mit_025c, Mit_025b, Mit_025a, Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024, Mit_023, Mit_022a, Mit_022, Mit_020d, TMit_020d, Mit_020c, Mit_021, Mit_021pre2, Mit_021pre1, Mit_020b, Mit_020a, Mit_020, Mit_020pre1, Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, V07-05-00, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, ConvRejection-10-06-09, Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, Mit_013pre1, Mit_012i, Mit_012g, Mit_012f, Mit_012e, Mit_012d, Mit_012c, Mit_012b, Mit_012a, Mit_012, Mit_011a, HEAD
Branch point for: Mit_025c_branch
Changes since 1.2: +2 -1 lines
Log Message:
Use CLHEP namespace.

File Contents

# User Rev Content
1 loizides 1.3 // $Id: MultiVertexFitterD.cc,v 1.2 2009/03/20 13:33:03 loizides Exp $
2 paus 1.1
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 loizides 1.3 using namespace CLHEP;
19 paus 1.1
20     jmp_buf denv;
21 loizides 1.2
22     //--------------------------------------------------------------------------------------------------
23 paus 1.1 extern "C" void MultiVertexFitterDSetStatus(int i) {
24     cout << "Warning, you are handling a severe error in MultiVertexFitterD" << endl;
25     longjmp(denv,-66);
26     }
27    
28 loizides 1.2 //--------------------------------------------------------------------------------------------------
29 paus 1.1 MultiVertexFitterD::MultiVertexFitterD() :
30     _bField (3.8), // default B field for running
31     _currentAllocatedVertexNumber(0), // facilitates CandNode recursion
32     _referencePoint (0,0,0), // set reference point to (0,0,0) initially
33     _primaryVertex (0,0,0), // set primary vertex to (0,0,0) initially
34     _cdfPrimaryVertex (0,0,0) // set pv in CMS coords to (0,0,0) initially
35     {
36     // Set name and email of MultiVertexFitterD expert
37     _expert="Christoph Paus (paus@mit.edu)";
38    
39     // First get pointers to various FORTAN common blocks
40     _ctvmq_com = (CTVMQ*) dctvmq_address_(); //printf(" Common: _ctvmq_com %p\n", _ctvmq_com );
41     _ctvmfr_com = (CTVMFR*) dctvmfr_address_(); //printf(" Common: _ctvmfr_com %p\n", _ctvmfr_com);
42     _fiddle_com = (FIDDLE*) dfiddle_address_(); //printf(" Common: _fiddle_com %p\n", _fiddle_com);
43     _trkprm_com = (TRKPRM*) dtrkprm_address_(); //printf(" Common: _trkprm_com %p\n", _trkprm_com);
44    
45     // Initialize various arrays
46     init();
47     // Don't bomb program on error.
48     _fiddle.excuse = 1;
49     // Don't extrapolate track errors by default
50     _extrapolateTrackErrors = false;
51     }
52    
53 loizides 1.2 //--------------------------------------------------------------------------------------------------
54 paus 1.1 void MultiVertexFitterD::init(double bField)
55     {
56     // Set internal variable which keeps track of the b field
57     _bField = bField;
58    
59     // Now initialize CTVMQ common. Run and trigger numbers are dummies - they are not used.
60     _ctvmq.runnum = 1;
61     _ctvmq.trgnum = 100;
62     // Eventually, we have to get the magnetic field from the right place.
63     // Origignal with bmag = 14.116 [kGauss]:
64     // _ctvmq.pscale = 0.000149896 * bmag;
65     // New bfield in Tesla bmag = 1.4116 [T]:
66     _ctvmq.pscale = 0.00149896 * bField;
67     // Set the default maximum chi-square per dof.
68     _fiddle.chisqmax = 225.0;
69     // We also need to get the primary vertex from the right place, but for now we put in (0,0,0).
70     setPrimaryVertex(0.0,0.0,0.0);
71     float xverr[3][3];
72     for (int j = 0; j<3; j++) {
73     for (int k = 0; k<3; k++) {
74     xverr[j][k] = 0.;
75     }
76     }
77     xverr[0][0] = 0.005;
78     xverr[1][1] = 0.005;
79     xverr[2][2] = 1.0;
80     setPrimaryVertexError(xverr);
81     // Zero number of tracks, vertices and mass constraints
82     _ctvmq.ntrack = 0;
83     _ctvmq.nvertx = 0;
84     _ctvmq.nmassc = 0;
85     // Zero track list and arrays containing the vertex and mass constraint configuration
86     for (int j=0; j<_maxtrk; ++j) {
87     _ctvmq.list[j]=0;
88     for (int jv=0; jv<_maxvtx; ++jv)
89     _ctvmq.trkvtx[jv][j] = false;
90     for (int jmc=0; jmc<_maxmcn; ++jmc)
91     _ctvmq.trkmcn[jmc][j] = false;
92     }
93     // Initialize the conversion and vertex pointing arrays
94     for (int jv=0; jv<_maxvtx; ++jv) {
95     _ctvmq.cvtx[jv] = 0;
96     _ctvmq.vtxpnt[0][jv] = -1;
97     _ctvmq.vtxpnt[1][jv] = 0;
98     }
99     _ctvmq.drmax = 2.0;
100     _ctvmq.dzmax = 20.0;
101     _ctvmq.rvmax = 70.0;
102     _ctvmq.trnmax = 0.5;
103     _ctvmq.dsmin = -100.0; // for CDF used -2, for CMS fuzzing around.
104     // Set _stat to -999 and chisqq to -1 to symbolize that no fit has yet been done.
105     _stat = -999;
106     _ctvmq.chisqr[0] = -1.0;
107    
108     _primaryVertex = Hep3Vector(0,0,0);
109     _cdfPrimaryVertex = Hep3Vector(0,0,0);
110     _referencePoint = ThreeVector(0,0,0);
111     }
112    
113 loizides 1.2 //--------------------------------------------------------------------------------------------------
114 paus 1.1 bool MultiVertexFitterD::addTrack(const HepVector &v, const HepSymMatrix &cov,
115     int trackId, float mass, vertexNumber jv)
116     {
117     // Check that this vertex number is within the allowed range.
118     if (jv<VERTEX_1 || jv>_maxvtx)
119     return false;
120     _ctvmq.nvertx = jv>_ctvmq.nvertx ? jv : _ctvmq.nvertx;
121    
122     // Add track extrapolation, if the user set it
123     HepVector w = v;
124     HepSymMatrix m = cov;
125     //if (_extrapolateTrackErrors)
126     // // This will move the reference point
127     // moveReferencePoint(w, m);
128    
129     // Check that we have not exceeded the maximum number of tracks.
130     if (_ctvmq.ntrack>=_maxtrk)
131     return false;
132    
133     // Add this track
134     _ctvmq.list [_ctvmq.ntrack] = trackId;
135     _ctvmq.tkbank[_ctvmq.ntrack][0] = 'Q';
136     _ctvmq.tkbank[_ctvmq.ntrack][1] = 'T';
137     _ctvmq.tkbank[_ctvmq.ntrack][2] = 'R';
138     _ctvmq.tkbank[_ctvmq.ntrack][3] = 'K';
139     _ctvmq.tmass [_ctvmq.ntrack] = mass;
140     _ctvmq.trkvtx[jv-1][_ctvmq.ntrack] = true;
141    
142     // Put this track's helix parameters and error matrix into a fortran common block so that they
143     // can be accessed by gettrk. This is a dummy for now.
144     _trkprm.trhelix[_ctvmq.ntrack][0] = w[0];
145     _trkprm.trhelix[_ctvmq.ntrack][1] = w[1];
146     _trkprm.trhelix[_ctvmq.ntrack][2] = w[2];
147     _trkprm.trhelix[_ctvmq.ntrack][3] = w[3];
148     _trkprm.trhelix[_ctvmq.ntrack][4] = w[4];
149    
150     for (int j=0; j<5; ++j) {
151     for (int k=0; k<5; ++k)
152     _trkprm.trem[_ctvmq.ntrack][j][k]=m[j][k];
153     }
154     _ctvmq.ntrack++;
155    
156     return true;
157     }
158    
159 loizides 1.2 //--------------------------------------------------------------------------------------------------
160 paus 1.1 bool MultiVertexFitterD::addTrack(const TVectorD &v, const TMatrixDSym &cov,
161     int trackId, double mass, vertexNumber jv)
162     {
163     // Check that this vertex number is within the allowed range.
164     if (jv<VERTEX_1 || jv>_maxvtx)
165     return false;
166     _ctvmq.nvertx = jv>_ctvmq.nvertx ? jv : _ctvmq.nvertx;
167    
168     //// Add track extrapolation, if the user set it
169     //HepVector w = v;
170     //HepSymMatrix m = cov;
171     TVectorD w(v);
172     TMatrixDSym m(cov);
173     ////if (_extrapolateTrackErrors)
174     //// // This will move the reference point
175     //// moveReferencePoint(w, m);
176    
177     // Check that we have not exceeded the maximum number of tracks.
178     if (_ctvmq.ntrack>=_maxtrk)
179     return false;
180    
181     // Add this track
182     _ctvmq.list [_ctvmq.ntrack] = trackId;
183     _ctvmq.tkbank[_ctvmq.ntrack][0] = 'Q';
184     _ctvmq.tkbank[_ctvmq.ntrack][1] = 'T';
185     _ctvmq.tkbank[_ctvmq.ntrack][2] = 'R';
186     _ctvmq.tkbank[_ctvmq.ntrack][3] = 'K';
187     _ctvmq.tmass [_ctvmq.ntrack] = mass;
188     _ctvmq.trkvtx[jv-1][_ctvmq.ntrack] = true;
189    
190     // Put this track's helix parameters and error matrix into a fortran common block so that they
191     // can be accessed by gettrk. This is a dummy for now.
192     _trkprm.trhelix[_ctvmq.ntrack][0] = w[0];
193     _trkprm.trhelix[_ctvmq.ntrack][1] = w[1];
194     _trkprm.trhelix[_ctvmq.ntrack][2] = w[2];
195     _trkprm.trhelix[_ctvmq.ntrack][3] = w[3];
196     _trkprm.trhelix[_ctvmq.ntrack][4] = w[4];
197    
198     for (int j=0; j<5; ++j) {
199     for (int k=0; k<5; ++k)
200     _trkprm.trem[_ctvmq.ntrack][j][k]=m[j][k];
201     }
202     _ctvmq.ntrack++;
203    
204     return true;
205     }
206    
207 loizides 1.2 //--------------------------------------------------------------------------------------------------
208 paus 1.1 bool MultiVertexFitterD::vertexPoint_2d(vertexNumber jv1, vertexNumber jv2)
209     {
210     // Check that these vertex numbers are within allowed range and that the vertices are unique.
211     if (jv1>_maxvtx || jv1<VERTEX_1)
212     return false;
213     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
214     return false;
215     if (jv1 <= jv2)
216     return false;
217    
218     // Setup vertex pointing.
219     _ctvmq.vtxpnt[0][jv1-1] = jv2;
220     _ctvmq.vtxpnt[1][jv1-1] = 1; // 2d pointing.
221    
222     return true;
223     }
224    
225 loizides 1.2 //--------------------------------------------------------------------------------------------------
226 paus 1.1 bool MultiVertexFitterD::vertexPoint_3d(vertexNumber jv1, vertexNumber jv2)
227     {
228     // Check that these vertex numbers are within allowed range and that the vertices are distinct
229     if (jv1>_maxvtx || jv1<VERTEX_1)
230     return false;
231     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
232     return false;
233     if (jv1 <= jv2)
234     return false;
235    
236     // Setup vertex pointing
237     _ctvmq.vtxpnt[0][jv1-1] = jv2;
238     _ctvmq.vtxpnt[1][jv1-1] = 2; // 3d pointing
239    
240     return true;
241     }
242    
243 loizides 1.2 //--------------------------------------------------------------------------------------------------
244 paus 1.1 bool MultiVertexFitterD::vertexPoint_1track(vertexNumber jv1, vertexNumber jv2)
245     {
246     // Check that these vertex numbers are within allowed range and are distinct
247     if (jv1>_maxvtx || jv1<VERTEX_1)
248     return false;
249     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
250     return false;
251     if (jv1 <= jv2)
252     return false;
253    
254     // Setup vertex pointing
255     _ctvmq.vtxpnt[0][jv1-1] = jv2;
256     _ctvmq.vtxpnt[1][jv1-1] = 3; // Point to 1 track vertex
257    
258     return true;
259     }
260    
261 loizides 1.2 //--------------------------------------------------------------------------------------------------
262 paus 1.1 bool MultiVertexFitterD::vertexPoint_0track(vertexNumber jv1, vertexNumber jv2)
263     {
264     // jv2 is the zero track vertex. jv1 is the multi track vertex which points to jv2
265    
266     // Note: You must call this routine at least twice in order for ctvmft_zerotrackvtx to work ie
267     // there must be at least 2 vertices pointed at a zero track vertex. ctvmft for this, but the
268     // error message may not make it to your log file (look at the local variables in the stack frame
269     // especially IJKERR(2). The significance of this error code is documented at the top of whatever
270     // routine chucked you out (ctvmf00 in this case)
271    
272     // see ctvmft.f source file for discussion. See especially comments at the top of subroutines:
273     // ctvmft and ctvmfa
274    
275     // Check that these vertex numbers are within allowed range and are distinct.
276     if (jv1>_maxvtx || jv1<VERTEX_1)
277     return false;
278     if (jv2>_maxvtx || jv2<PRIMARY_VERTEX)
279     return false;
280     if (jv1 <= jv2)
281     return false;
282    
283     // Setup vertex pointing.
284     _ctvmq.vtxpnt[0][jv1-1] = jv2;
285     _ctvmq.vtxpnt[1][jv1-1] = 4; // Point to 0 track vertex.
286    
287     return true;
288     }
289    
290 loizides 1.2 //--------------------------------------------------------------------------------------------------
291 paus 1.1 bool MultiVertexFitterD::conversion_2d(vertexNumber jv)
292     {
293     if (jv<VERTEX_1 || jv>_ctvmq.nvertx)
294     return false;
295    
296     _ctvmq.cvtx[jv-1] = 1;
297    
298     return true;
299     }
300    
301 loizides 1.2 //--------------------------------------------------------------------------------------------------
302 paus 1.1 bool MultiVertexFitterD::conversion_3d(vertexNumber jv)
303     {
304     if (jv<VERTEX_1 || jv>_ctvmq.nvertx)
305     return false;
306    
307     _ctvmq.cvtx[jv-1] = 2;
308    
309     return true;
310     }
311    
312 loizides 1.2 //--------------------------------------------------------------------------------------------------
313 paus 1.1 bool MultiVertexFitterD::massConstrain(int ntrk, const int trkIds[], float mass)
314     {
315     // Check that we have not exceeded the allowed number of mass constraints.
316     if (_ctvmq.nmassc>=_maxmcn)
317     return false;
318    
319     // Set constraint mass
320     _ctvmq.cmass[_ctvmq.nmassc]=mass;
321    
322     // For each track in contraint, set trkmcn true. Since the number in tracks[] is the track
323     // number, we have to find each track in the list of tracks.
324     for (int jt=0; jt<ntrk; ++jt) {
325     bool found=false;
326     for (int kt=0; kt<_ctvmq.ntrack; ++kt) {
327     if (trkIds[jt] == _ctvmq.list[kt]) {
328     _ctvmq.trkmcn[_ctvmq.nmassc][kt]=true;
329     found=true;
330     }
331     }
332     if (!found)
333     return false;
334     }
335    
336     // Increment number of mass constraints.
337     _ctvmq.nmassc++;
338    
339     return true;
340     }
341    
342 loizides 1.2 //--------------------------------------------------------------------------------------------------
343 paus 1.1 bool MultiVertexFitterD::beamlineConstraint(float xb, float yb, HepSymMatrix berr,
344     float xzbslope, float yzbslope)
345     {
346     // Set beam position at z=0
347     setPrimaryVertex(xb,yb,0);
348     //if (_extrapolateTrackErrors) {
349     // float newXb = xb - _referencePoint.x() + _referencePoint.z() * xzbslope;
350     // float newYb = yb - _referencePoint.y() + _referencePoint.z() * yzbslope;
351     // setPrimaryVertex(newXb, newYb, 0);
352     //}
353    
354     bool success = setPrimaryVertexError(berr);
355    
356     // Set the beamline slope values
357     _ctvmq.xzslope = xzbslope;
358     _ctvmq.yzslope = yzbslope;
359    
360     // Turn ON beamline constraint
361     _ctvmq.vtxpnt[0][0] = -100;
362    
363     return success;
364     }
365    
366 loizides 1.2 //--------------------------------------------------------------------------------------------------
367 paus 1.1 bool MultiVertexFitterD::beamlineConstraint(Hep3Vector pv, HepSymMatrix berr, float xzbslope,
368     float yzbslope)
369     {
370     // Check if input beam position coordinates are at z=0
371     if (pv.z() != 0)
372     return false;
373    
374     return beamlineConstraint(pv.x(),pv.y(),berr,xzbslope,yzbslope);
375     }
376    
377 loizides 1.2 //--------------------------------------------------------------------------------------------------
378 paus 1.1 void MultiVertexFitterD::setPrimaryVertex(float xv, float yv, float zv)
379     {
380     // Set x,y,z position of the primary vertex.
381     _ctvmq.xyzpv0[0] = xv;
382     _ctvmq.xyzpv0[1] = yv;
383     _ctvmq.xyzpv0[2] = zv;
384    
385     _primaryVertex = Hep3Vector( xv, yv, zv );
386     }
387    
388 loizides 1.2 //--------------------------------------------------------------------------------------------------
389 paus 1.1 void MultiVertexFitterD::setPrimaryVertex(Hep3Vector pv)
390     {
391     // Set x,y,z position of the primary vertex.
392     _ctvmq.xyzpv0[0] = pv.x();
393     _ctvmq.xyzpv0[1] = pv.y();
394     _ctvmq.xyzpv0[2] = pv.z();
395    
396     _primaryVertex = pv;
397     }
398    
399 loizides 1.2 //--------------------------------------------------------------------------------------------------
400 paus 1.1 void MultiVertexFitterD::setPrimaryVertexError(const float xverr[3][3])
401     {
402     // Set the error matrix for the primary vertex.
403     for (int j=0; j<3; ++j) {
404     for (int k=0; k<3; ++k)
405     _ctvmq.exyzpv[j][k]=xverr[j][k];
406     }
407     }
408    
409 loizides 1.2 //--------------------------------------------------------------------------------------------------
410 paus 1.1 bool MultiVertexFitterD::setPrimaryVertexError(const HepSymMatrix &xverr)
411     {
412     // Set the error matrix for the primary vertex using a HepSymMatrix. First check that the matrix
413     // is the correct size.
414     if (xverr.num_row() != 3)
415     return false;
416     for (int j=0; j<3; j++) {
417     for (int k=0; k<3; k++)
418     _ctvmq.exyzpv[j][k]=xverr[j][k];
419     }
420     return true;
421     }
422    
423 loizides 1.2 //--------------------------------------------------------------------------------------------------
424 paus 1.1 bool MultiVertexFitterD::fit()
425     {
426     // Check that the diagonal elements of all the track error matrices are positive
427     bool mstat = true;
428     for (int trk=0; trk<_ctvmq.ntrack; ++trk) {
429     for (int j=0; j<5; ++j) {
430     // Check diagonal elements of error matrix.
431     if (_trkprm.trem[trk][j][j] < 0.) {
432     // The covariance matrix could not be inverted: Set the error codes and fail this fit
433     mstat = false;
434     _ctvmq.ijkerr[0] = 3;
435     _ctvmq.ijkerr[1] = 2;
436     _ctvmq.ijkerr[2] = trk + 1;
437     }
438     }
439     // Check that curvature of track is reasonable: Pt is above ~10MeV/c. If not, set the error
440     // codes and fail this fit
441     if (fabs(_trkprm.trhelix[trk][1]) > 0.1) {
442     //if (fabs(_trkprm.trhelix[trk][1]) > 0.01) {
443     mstat = false;
444     _ctvmq.ijkerr[0] = 3;
445     _ctvmq.ijkerr[1] = 5;
446     _ctvmq.ijkerr[2] = trk + 1;
447     }
448     }
449     // If there was a problem with any track, fail the fit
450     if (!mstat) {
451     _stat = 1;
452     return false;
453     }
454    
455     // First copy information into CTVMFT common blocks
456     *_ctvmq_com = _ctvmq;
457     *_ctvmfr_com = _ctvmfr;
458     *_fiddle_com = _fiddle;
459     *_trkprm_com = _trkprm;
460     // Do the vertex fit.
461     int print = 0;
462     int level = 0;
463    
464     // #if ( defined(LINUX) && defined(__USE_BSD) ) || defined(OSF1)
465     // struct sigaction myaction = {MultiVertexFitterDSetStatus, 0, 0, 0}, oldaction;
466     // sigaction(SIGFPE, &myaction, &oldaction);
467     // if (setjmp(denv)!=0) {
468     // sigaction(SIGFPE, &oldaction,0);
469     // return -999;
470     // }
471     // #endif
472    
473     dctvmft_(print,level,_stat);
474    
475     // #if ( defined(LINUX) && defined(__USE_BSD) ) || defined(OSF1)
476     // sigaction(SIGFPE, &oldaction,0);
477     // #endif
478    
479     // Now copy information from CTVMFT common blocks to local storage
480     _ctvmq = *_ctvmq_com;
481     _ctvmfr = *_ctvmfr_com;
482     _fiddle = *_fiddle_com;
483     _trkprm = *_trkprm_com;
484    
485     return (_stat == 0);
486     }
487    
488 loizides 1.2 //--------------------------------------------------------------------------------------------------
489 paus 1.1 void MultiVertexFitterD::print() const
490     {
491     print(cout);
492     }
493    
494 loizides 1.2 //--------------------------------------------------------------------------------------------------
495 paus 1.1 void MultiVertexFitterD::print(ostream& os) const
496     {
497     os << "****************************** MultiVertexFitterD "
498     << "******************************" << endl;
499     os << "Number of tracks: " << _ctvmq.ntrack << endl;
500     os << " Tracks: ";
501     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
502     if (jt != 0) os << ", ";
503     os << _ctvmq.list[jt];
504     }
505     os << endl;
506     os << "Number of vertices: " << _ctvmq.nvertx << endl;
507     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
508     os << " Vertex " << jv+1 << " tracks: ";
509     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
510     if (_ctvmq.trkvtx[jv][jt]) {
511     os << " " << _ctvmq.list[jt];
512     }
513     }
514     os << endl;
515     }
516     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
517     if (_ctvmq.vtxpnt[0][jv]==0) {
518     os << " Vertex " << jv+1 << " points to the primary vertex ";
519     }
520     else if (_ctvmq.vtxpnt[0][jv]>0) {
521     os << " Vertex " << jv+1 << " points to vertex "
522     << _ctvmq.vtxpnt[0][jv];
523     }
524     if (_ctvmq.vtxpnt[1][jv]==1) {
525     os << " in 2 dimensions" << endl;
526     }
527     else if (_ctvmq.vtxpnt[1][jv]==2) {
528     os << " in 3 dimensions" << endl;
529     }
530     else if (_ctvmq.vtxpnt[1][jv]==3) {
531     os << ", a single track vertex" << endl;
532     }
533     if (_ctvmq.cvtx[jv]>0) {
534     os << " Vertex " << jv+1 << " is a conversion" << endl;
535     }
536     }
537     os << "Number of mass constraints: " << _ctvmq.nmassc << endl;
538     for (int jmc=0; jmc<_ctvmq.nmassc; ++jmc) {
539     os << " Tracks ";
540     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
541     if (_ctvmq.trkmcn[jmc][jt]) {
542     os << " " << _ctvmq.list[jt];
543     }
544     }
545     os << " constrained to mass " << _ctvmq.cmass[jmc]
546     << " Gev/c^2" << endl;
547     }
548     if (_stat==-999) {
549     os << "No fit has been done." << endl;
550     }
551     else {
552     os << "***** Results of Fit *****" << endl;
553     printErr(os);
554     os << " Status = " << _stat << endl;
555     os.precision(7);
556     os << " Chi-square = " << scientific << _ctvmq.chisqr[0]
557     << " for " << _ctvmq.ndof << " degrees of freedom." << endl;
558     os << " => probability = " << prob() << endl;
559     for (int jv=0; jv<_ctvmq.nvertx; ++jv) {
560     os << "Vertex " << jv+1
561     << " position: " << scientific
562     << _ctvmq.xyzvrt[jv+1][0] << " "
563     << _ctvmq.xyzvrt[jv+1][1] << " "
564     << _ctvmq.xyzvrt[jv+1][2] << endl;
565     }
566     for (int jt=0; jt<_ctvmq.ntrack; ++jt) {
567     os << "Track " << _ctvmq.list[jt]
568     << " - P4: " << scientific
569     << _ctvmq.trkp4[0][jt] << " "
570     << _ctvmq.trkp4[1][jt] << " "
571     << _ctvmq.trkp4[2][jt] << " "
572     << _ctvmq.trkp4[3][jt] << " "
573     << " - PT: " << scientific
574     << sqrt(_ctvmq.trkp4[0][jt]*_ctvmq.trkp4[0][jt]+
575     _ctvmq.trkp4[1][jt]*_ctvmq.trkp4[1][jt]) << endl;
576     }
577     }
578     os << "****************************************"
579     << "**************************" << endl;
580    
581     return;
582     }
583    
584 loizides 1.2 //--------------------------------------------------------------------------------------------------
585 paus 1.1 void MultiVertexFitterD::printErr() const
586     {
587     printErr(cout);
588     }
589    
590 loizides 1.2 //--------------------------------------------------------------------------------------------------
591 paus 1.1 void MultiVertexFitterD::printErr(ostream& os) const
592     {
593     os << "MultiVertexFitterD: IJKERR = " << _ctvmq.ijkerr[0] << ", "
594     << _ctvmq.ijkerr[1] << ", "
595     << _ctvmq.ijkerr[2] << endl;
596     if (status()==0 && _ctvmq.ijkerr[0]==0) return;
597     if (_ctvmq.ijkerr[0] == -1) {
598     os << " Problem with GETTRK: track requested is not in list."
599     << endl
600     << " This should not happen - Contact MultiVertexFitterD expert "
601     << _expert << "." <<endl;
602     }
603     else if (_ctvmq.ijkerr[0]==1) {
604     os << " Problem in CTVM00:" << endl;
605     if (_ctvmq.ijkerr[1]==1) {
606     os << " Number of tracks is " << _ctvmq.ntrack
607     << "." << endl;
608     if (_ctvmq.ntrack < 2) {
609     os << ", which is too few (must be at least 2)." << endl;
610     }
611     else if (_ctvmq.ntrack > _maxtrk) {
612     os << ", which is too many (maximum is " << _maxtrk
613     << ")." << endl;
614     }
615     else {
616     os << " Problem with number of tracks"
617     << " for unknown reasons." << endl;
618     }
619     }
620     else if (_ctvmq.ijkerr[1]==2) {
621     os << " Number of vertices is " << _ctvmq.nvertx
622     << "." << endl;
623     if (_ctvmq.nvertx < 1) {
624     os << ", which is too few (must be at least 1)." << endl;
625     }
626     else if (_ctvmq.nvertx > _maxvtx) {
627     os << ", which is too many (maximum is " << _maxvtx
628     << ")." << endl;
629     }
630     else {
631     os << endl << " Problem with number of vertices"
632     << " for unknown reasons." << endl;
633     }
634     }
635     else if (_ctvmq.ijkerr[1]==3) {
636     os << " Number of mass constraints is " << _ctvmq.nmassc
637     << "." << endl;
638     if (_ctvmq.nmassc < 0) {
639     os << ", which is negative." << endl;
640     }
641     else if (_ctvmq.nmassc > _maxmcn) {
642     os << ", which is too many (maximum is " << _maxmcn
643     << ")." << endl;
644     }
645     else {
646     os << endl << " Problem with number of mass"
647     << " constraints for unknown reasons." << endl;
648     }
649     }
650     else if (_ctvmq.ijkerr[1]==11) {
651     os << " Vertex " << _ctvmq.ijkerr[2]
652     << " has less than one track." << endl;
653     }
654     else if (_ctvmq.ijkerr[1]==12) {
655     os << " Vertex " << _ctvmq.ijkerr[2]
656     << " is a conversion vertex with a number of tracks"
657     << " different than two." << endl;
658     }
659     else if (_ctvmq.ijkerr[1]==13) {
660     os << " Vertex " << _ctvmq.ijkerr[2]
661     << " is a one track vertex that has no multi-track"
662     << " descendents." << endl;
663     }
664     else if (_ctvmq.ijkerr[1]==14) {
665     os << " Vertex " << _ctvmq.ijkerr[2]
666     << " does not point at a vertex with a lower number."
667     << endl;
668     }
669     else if (_ctvmq.ijkerr[1]==15) {
670     os << " Vertex " << _ctvmq.ijkerr[2]
671     << " has a parent vertex that is a conversion." << endl;
672     }
673     else if (_ctvmq.ijkerr[1]==16) {
674     os << " Vertex " << _ctvmq.ijkerr[2]
675     << " does 1 track pointing to a vertex with"
676     << " more than 1 track." << endl;
677     }
678     else if (_ctvmq.ijkerr[1]==17) {
679     os << " Vertex " << _ctvmq.ijkerr[2]
680     << " does 0 track pointing to a vertex with"
681     << " more than 0 track (?)." << endl;
682     }
683     else if (_ctvmq.ijkerr[1]==19) {
684     os << " Primary vertex error matrix is singular." << endl;
685     }
686     else if (_ctvmq.ijkerr[1]==21) {
687     os << " Track with Id " << _ctvmq.ijkerr[2]
688     << "is not in any vertex." << endl;
689     }
690     else if (_ctvmq.ijkerr[1]==22) {
691     os << " Track with Id " << _ctvmq.ijkerr[2]
692     << "is in multiple vertices." << endl;
693     }
694     else if (_ctvmq.ijkerr[1]==23) {
695     os << " Track with Id " << _ctvmq.ijkerr[2]
696     << "occurs more than once." << endl;
697     }
698     else if (_ctvmq.ijkerr[1]==31) {
699     os << " A mass constraint has less than 2 tracks." << endl;
700     }
701     else if (_ctvmq.ijkerr[1]==32) {
702     os << " The sum masses of the tracks in a mass constraint"
703     << " exceeds the constraint mass." << endl;
704     }
705     else if (_ctvmq.ijkerr[1]==33) {
706     os << " Beamline constraint. Beam covariance not set properly."
707     << " Negative diagonal elements." << endl;
708     }
709     else if (_ctvmq.ijkerr[1]==34) {
710     os << " Beamline constraint. Beam covariance not set properly."
711     << " Off-diagonal elements not zero." << endl;
712     }
713     else if (_ctvmq.ijkerr[1]==36) {
714     os << " Beamline constraint. Number of vertices = "
715     << _ctvmq.nvertx << " Should be 1." << endl;
716     }
717     }
718     else if (_ctvmq.ijkerr[0] == 2) {
719     if (_ctvmq.ijkerr[1] == 20) {
720     os << " Problem in CTVM00: " << endl;
721     os << " Track has negative Id = "
722     << _ctvmq.list[_ctvmq.ijkerr[2]-1] << "." << endl;
723     }
724     else {
725     os << " Problem in CTVMFA with vertex "
726     << _ctvmq.ijkerr[2] << ": " << endl;
727     os << " Failure in vertex first approximation." << endl;
728     if (_ctvmq.ijkerr[1] == 1) {
729     os << " Tracks are concentric circles." << endl;
730     }
731     if (_ctvmq.ijkerr[1] == 2) {
732     os << " Conversion vertex has widely separated"
733     << " exterior circles at midpoint." << endl;
734     }
735     if (_ctvmq.ijkerr[1] == 3) {
736     os << " Conversion vertex has widely separated"
737     << " interior circles at midpoint." << endl;
738     }
739     if (_ctvmq.ijkerr[1] == 4) {
740     os << " Vertex has widely separated"
741     << " exterior circles at approximate vertex." << endl;
742     }
743     if (_ctvmq.ijkerr[1] == 5) {
744     os << " Vertex has widely separated"
745     << " interior circles at approximate vertex." << endl;
746     }
747     if (_ctvmq.ijkerr[1] == 6) {
748     os << " Rv is too large at the chosen"
749     << " intersection point." << endl;
750     }
751     if (_ctvmq.ijkerr[1] == 7) {
752     os << " Delta z is too large at the chosen"
753     << " intersection point." << endl;
754     }
755     if (_ctvmq.ijkerr[1] == 8) {
756     os << " A track's turning to the chosen vertex"
757     << " is too large." << endl;
758     }
759     if (_ctvmq.ijkerr[1] == 9) {
760     os << " There is no solution with an adequately"
761     << " positive arc length." << endl;
762     }
763     if (_ctvmq.ijkerr[1] == 21) {
764     os << " zero-track vertexing: either/both vertex "
765     << " momenta are too small (<0.01 MeV)." << endl;
766     }
767     if (_ctvmq.ijkerr[1] == 22) {
768     os << " zero-track vertexing: Two lines (tracks) are "
769     << " parallel/antiparallel." << endl;
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 MultiVertexFitterD expert " << _expert << "." << endl;
819     }
820     return;
821     }
822    
823 loizides 1.2 //--------------------------------------------------------------------------------------------------
824 paus 1.1 void MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
833 paus 1.1 int MultiVertexFitterD::getIJKErr0() const
834     {
835     return _ctvmq.ijkerr[0];
836     }
837 loizides 1.2
838     //--------------------------------------------------------------------------------------------------
839 paus 1.1 int MultiVertexFitterD::getIJKErr1() const
840     {
841     return _ctvmq.ijkerr[1];
842     }
843 loizides 1.2
844     //--------------------------------------------------------------------------------------------------
845 paus 1.1 int MultiVertexFitterD::getIJKErr2() const
846     {
847     return _ctvmq.ijkerr[2];
848     }
849    
850 loizides 1.2 //--------------------------------------------------------------------------------------------------
851 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
866 paus 1.1 string MultiVertexFitterD::expert() const
867     {
868     return _expert;
869     }
870    
871 loizides 1.2 //--------------------------------------------------------------------------------------------------
872 paus 1.1 int MultiVertexFitterD::status() const
873     {
874     return _stat;
875     }
876    
877 loizides 1.2 //--------------------------------------------------------------------------------------------------
878 paus 1.1 float MultiVertexFitterD::chisq() const
879     {
880     // Chi-square of fit
881     return _ctvmq.chisqr[0];
882     }
883    
884 loizides 1.2 //--------------------------------------------------------------------------------------------------
885 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
895 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
908 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
926 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
950 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
974 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1001 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1031 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1060 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1090 paus 1.1 FourVector MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1107 paus 1.1 float MultiVertexFitterD::getMass(int ntrk, const int trkIds[], float &dmass) const
1108     {
1109     // #if (defined(LINUX) && defined(__USE_BSD)) || defined(OSF1)
1110     // struct sigaction myaction = {MultiVertexFitterDSetStatus, 0, 0, 0}, oldaction;
1111     // sigaction(SIGFPE, &myaction, &oldaction);
1112     // if (setjmp(denv)!=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     dmcalc_(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.2 //--------------------------------------------------------------------------------------------------
1153 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1234 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1242 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1327 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1335 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1362 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1370 paus 1.1 float MultiVertexFitterD::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     ddcalc_(mvert,nvert,dxyz,dr,dz,dl);
1390     drerr = dl[0];
1391     return dr;
1392     }
1393    
1394 loizides 1.2 //--------------------------------------------------------------------------------------------------
1395 paus 1.1 float MultiVertexFitterD::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     ddcalc_(mvert,nvert,dxyz,dr,dz,dl);
1416     dzerr = dl[1];
1417     return dz;
1418     }
1419    
1420 loizides 1.2 //--------------------------------------------------------------------------------------------------
1421 paus 1.1 Hep3Vector MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1457 paus 1.1 ThreeVector MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1493 paus 1.1 ThreeSymMatrix MultiVertexFitterD::getErrorMatrix(MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1521 paus 1.1 double MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1536 paus 1.1 HepSymMatrix MultiVertexFitterD::getErrorMatrixHep(MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1564 paus 1.1 HepSymMatrix MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1588 paus 1.1 float MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1618 paus 1.1 void MultiVertexFitterD::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     //HepMatrix dp_dpar[_maxvtx] = { deriv, deriv, deriv };
1665     HepMatrix dp_dpar[_maxvtx] = { deriv };
1666    
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     //HepMatrix answer[_maxvtx] = { ans, ans, ans };
1791     HepMatrix answer[_maxvtx] = { ans };
1792     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.2 //--------------------------------------------------------------------------------------------------
1798 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1807 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1817 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1826 paus 1.1 int MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1835 paus 1.1 int MultiVertexFitterD::mOff() const
1836     {
1837     return _ctvmq.moff;
1838     }
1839    
1840 loizides 1.2 //--------------------------------------------------------------------------------------------------
1841 paus 1.1 double MultiVertexFitterD::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.2 //--------------------------------------------------------------------------------------------------
1852 paus 1.1 MultiVertexFitterD::vertexNumber MultiVertexFitterD::allocateVertexNumber()
1853     {
1854     if ((_currentAllocatedVertexNumber < PRIMARY_VERTEX) ||
1855     (_currentAllocatedVertexNumber > _maxvtx)) {
1856     cout << "MultiVertexFitterD::allocateVertexNumber: out of range!" << endl;
1857     return PRIMARY_VERTEX;
1858     }
1859     return vertexNumber(++_currentAllocatedVertexNumber);
1860     }
1861    
1862 loizides 1.2 //--------------------------------------------------------------------------------------------------
1863 paus 1.1 void MultiVertexFitterD::resetAllocatedVertexNumber()
1864     {
1865     _currentAllocatedVertexNumber = 0;
1866     }
1867    
1868 loizides 1.2 //--------------------------------------------------------------------------------------------------
1869 paus 1.1 void MultiVertexFitterD::restoreFromCommons()
1870     {
1871     _stat = 0;
1872     _ctvmq_com = (CTVMQ*) dctvmq_address_();
1873     _ctvmfr_com = (CTVMFR*) dctvmfr_address_();
1874     _fiddle_com = (FIDDLE*) dfiddle_address_();
1875     _trkprm_com = (TRKPRM*) dtrkprm_address_();
1876     _ctvmq = *_ctvmq_com;
1877     _ctvmfr = *_ctvmfr_com;
1878     _fiddle = *_fiddle_com;
1879     _trkprm = *_trkprm_com;
1880     }