Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DReferenceTrajectory.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DReferenceTrajectory.cc
4 // Created: Wed Jul 19 13:42:58 EDT 2006
5 // Creator: davidl (on Darwin swire-b241.jlab.org 8.7.0 powerpc)
6 //
7 
8 #include <signal.h>
9 #include <memory>
10 #include <cmath>
11 
12 #include <DVector3.h>
13 using namespace std;
14 #include <math.h>
15 #include <algorithm>
16 
17 #include "DReferenceTrajectory.h"
18 #include "DTrackCandidate.h"
19 #include "DMagneticFieldStepper.h"
20 #include <TMatrix.h>
21 #include "HDGEOMETRY/DRootGeom.h"
22 #define ONE_THIRD 0.33333333333333333
23 #define TWO_THIRD 0.66666666666666667
24 #define EPS 1e-8
25 #define QuietNaN std::numeric_limits<double>::quiet_NaN()
26 
28 
29 //---------------------------------
30 // DReferenceTrajectory (Constructor)
31 //---------------------------------
33  , double q
34  , swim_step_t *swim_steps
35  , int max_swim_steps
36  , double step_size)
37 {
38  // Copy some values into data members
39  this->q = q;
40  this->step_size = step_size;
41  this->bfield = bfield;
42  this->Nswim_steps = 0;
43  this->dist_to_rt_depth = 0;
44  this->mass = 0.13957; // assume pion mass until otherwise specified
45  this->mass_sq=this->mass*this->mass;
46  this->hit_cdc_endplate = false;
47  this->RootGeom=NULL;
48  this->geom = NULL;
49  this->ploss_direction = kForward;
50  this->check_material_boundaries = true;
51  this->zmin_track_boundary = -100.0; // boundary at which to stop swimming
52  this->zmax_track_boundary = 670.0; // boundary at which to stop swimming
53  this->Rsqmax_interior = 65.0*65.0; // Maximum radius (in cm) corresponding to inside of BCAL
54  this->Rsqmax_exterior = 88.0*88.0; // Maximum radius (in cm) corresponding to outside of BCAL
55 
56  this->last_phi = 0.0;
57  this->last_swim_step = NULL;
58  this->last_dist_along_wire = 0.0;
59  this->last_dz_dphi = 0.0;
60 
61  this->debug_level = 0;
62 
63  // Initialize some values from configuration parameters
64  BOUNDARY_STEP_FRACTION = 0.80;
65  MIN_STEP_SIZE = 0.1; // cm
66  MAX_STEP_SIZE = 3.0; // cm
67  int MAX_SWIM_STEPS = 2500;
68 
69  gPARMS->SetDefaultParameter("TRK:BOUNDARY_STEP_FRACTION" , BOUNDARY_STEP_FRACTION, "Fraction of estimated distance to boundary to use as step size");
70  gPARMS->SetDefaultParameter("TRK:MIN_STEP_SIZE" , MIN_STEP_SIZE, "Minimum step size in cm to take when swimming a track with adaptive step sizes");
71  gPARMS->SetDefaultParameter("TRK:MAX_STEP_SIZE" , MAX_STEP_SIZE, "Maximum step size in cm to take when swimming a track with adaptive step sizes");
72  gPARMS->SetDefaultParameter("TRK:MAX_SWIM_STEPS" , MAX_SWIM_STEPS, "Number of swim steps for DReferenceTrajectory to allocate memory for (when not using external buffer)");
73  gPARMS->SetDefaultParameter("TRK:DEBUG_LEVEL" , this->debug_level);
74 
75  // It turns out that the greatest bottleneck in speed here comes from
76  // allocating/deallocating the large block of memory required to hold
77  // all of the trajectory info. The preferred way of calling this is
78  // with a pointer allocated once at program startup. This code block
79  // though allows it to be allocated here if necessary.
80  if(!swim_steps){
81  own_swim_steps = true;
82  this->max_swim_steps = MAX_SWIM_STEPS;
83  this->swim_steps = new swim_step_t[this->max_swim_steps];
84  }else{
85  own_swim_steps = false;
86  this->max_swim_steps = max_swim_steps;
87  this->swim_steps = swim_steps;
88  }
89 }
90 
91 //---------------------------------
92 // DReferenceTrajectory (Copy Constructor)
93 //---------------------------------
95 {
96  /// The copy constructor will always allocate its own memory for the
97  /// swim steps and set its internal flag to indicate that is owns them
98  /// regardless of the owner of the source trajectory's.
99 
100  this->Nswim_steps = rt.Nswim_steps;
101  this->q = rt.q;
102  this->max_swim_steps = rt.max_swim_steps;
103  this->own_swim_steps = true;
104  this->step_size = rt.step_size;
105  this->bfield = rt.bfield;
106  this->last_phi = rt.last_phi;
107  this->last_dist_along_wire = rt.last_dist_along_wire;
108  this->last_dz_dphi = rt.last_dz_dphi;
109  this->RootGeom = rt.RootGeom;
110  this->geom = rt.geom;
111  this->dist_to_rt_depth = 0;
112  this->mass = rt.GetMass();
113  this->mass_sq=this->mass*this->mass;
114  this->ploss_direction = rt.ploss_direction;
115  this->check_material_boundaries = rt.GetCheckMaterialBoundaries();
116  this->BOUNDARY_STEP_FRACTION = rt.GetBoundaryStepFraction();
117  this->MIN_STEP_SIZE = rt.GetMinStepSize();
118  this->MAX_STEP_SIZE = rt.GetMaxStepSize();
119  this->debug_level=rt.debug_level;
120  this->zmin_track_boundary = -100.0; // boundary at which to stop swimming
121  this->zmax_track_boundary = 670.0; // boundary at which to stop swimming
122  this->Rsqmax_interior = 65.0*65.0; // Maximum radius (in cm) corresponding to inside of BCAL
123  this->Rsqmax_exterior = 88.0*88.0; // Maximum radius (in cm) corresponding to outside of BCAL
124 
125 
126  this->swim_steps = new swim_step_t[this->max_swim_steps];
127  this->last_swim_step = NULL;
128  for(int i=0; i<Nswim_steps; i++)
129  {
130  swim_steps[i] = rt.swim_steps[i];
131  if(&(rt.swim_steps[i]) == rt.last_swim_step)
132  this->last_swim_step = &(swim_steps[i]);
133  }
134 
135 }
136 
137 //---------------------------------
138 // operator= (Assignment operator)
139 //---------------------------------
141 {
142  /// The assignment operator will always make sure the memory allocated
143  /// for the swim_steps is owned by the object being copied into.
144  /// If it already owns memory of sufficient size, then it will be
145  /// reused. If it owns memory that is too small, it will be freed and
146  /// a new block allocated. If it does not own its swim_steps coming
147  /// in, then it will allocate memory so that it does own it on the
148  /// way out.
149 
150  if(&rt == this)return *this; // protect against self copies
151 
152  // Free memory if block is too small
153  if(own_swim_steps==true && max_swim_steps<rt.Nswim_steps){
154  delete[] swim_steps;
155  swim_steps=NULL;
156  }
157 
158  // Forget memory block if we don't currently own it
159  if(!own_swim_steps){
160  swim_steps=NULL;
161  }
162 
163  this->Nswim_steps = rt.Nswim_steps;
164  this->q = rt.q;
165  this->max_swim_steps = rt.max_swim_steps;
166  this->own_swim_steps = true;
167  this->step_size = rt.step_size;
168  this->bfield = rt.bfield;
169  this->last_phi = rt.last_phi;
170  this->last_dist_along_wire = rt.last_dist_along_wire;
171  this->last_dz_dphi = rt.last_dz_dphi;
172  this->RootGeom = rt.RootGeom;
173  this->geom = rt.geom;
174  this->dist_to_rt_depth = rt.dist_to_rt_depth;
175  this->mass = rt.GetMass();
176  this->mass_sq=this->mass*this->mass;
177  this->ploss_direction = rt.ploss_direction;
178  this->check_material_boundaries = rt.GetCheckMaterialBoundaries();
179  this->BOUNDARY_STEP_FRACTION = rt.GetBoundaryStepFraction();
180  this->MIN_STEP_SIZE = rt.GetMinStepSize();
181  this->MAX_STEP_SIZE = rt.GetMaxStepSize();
182 
183  // Allocate memory if needed
184  if(swim_steps==NULL)this->swim_steps = new swim_step_t[this->max_swim_steps];
185 
186  // Copy swim steps
187  this->last_swim_step = NULL;
188  for(int i=0; i<Nswim_steps; i++)
189  {
190  swim_steps[i] = rt.swim_steps[i];
191  if(&(rt.swim_steps[i]) == rt.last_swim_step)
192  this->last_swim_step = &(swim_steps[i]);
193  }
194 
195 
196  return *this;
197 }
198 
199 //---------------------------------
200 // ~DReferenceTrajectory (Destructor)
201 //---------------------------------
203 {
204  if(own_swim_steps){
205  delete[] swim_steps;
206  }
207 }
208 
209 //---------------------------------
210 // CopyWithShift
211 //---------------------------------
213 {
214  // First, do a straight copy
215  *this = *rt;
216 
217  // Second, shift all positions
218  for(int i=0; i<Nswim_steps; i++)swim_steps[i].origin += shift;
219 }
220 
221 
222 //---------------------------------
223 // Reset
224 //---------------------------------
226  //reset DReferenceTrajectory for re-use
227  this->Nswim_steps = 0;
228  this->ploss_direction = kForward;
229  this->mass = 0.13957; // assume pion mass until otherwise specified
230  this->mass_sq=this->mass*this->mass;
231  this->hit_cdc_endplate = false;
232  this->last_phi = 0.0;
233  this->last_swim_step = NULL;
234  this->last_dist_along_wire = 0.0;
235  this->last_dz_dphi = 0.0;
236  this->dist_to_rt_depth = 0;
237  this->check_material_boundaries = true;
238 }
239 
240 //---------------------------------
241 // FastSwim -- light-weight swim to a wire that does not treat multiple
242 // scattering but does handle energy loss.
243 // No checks for distance to boundaries are done.
244 //---------------------------------
246  DVector3 &last_pos,DVector3 &last_mom,
247  double q,double smax,
248  const DCoordinateSystem *wire){
249  const DVector3 wire_origin=wire->origin;
250  const DVector3 wire_dir=wire->udir;
251  FastSwim(pos,mom,last_pos,last_mom,q,wire_origin,wire_dir,smax);
252 }
254  DVector3 &last_pos,DVector3 &last_mom,
255  double q,
256  const DVector3 &origin,
257  const DVector3 &dir,double smax){
258 
259  DVector3 mypos(pos);
260  DVector3 mymom(mom);
261 
262  // Initialize the stepper
263  DMagneticFieldStepper stepper(bfield, q, &pos, &mom);
264  double s=0,doca=1000.,old_doca=1000.,dP_dx=0.;
265  // double mass=GetMass();
266  while (s<smax){
267  // Save old value of doca
268  old_doca=doca;
269 
270  // Adjust step size to take smaller steps in regions of high momentum loss
271  if(geom){
272  double KrhoZ_overA=0.0;
273  double rhoZ_overA=0.0;
274  double LogI=0.0;
275  double X0=0.0;
276  if (geom->FindMatALT1(mypos,mymom,KrhoZ_overA,rhoZ_overA,LogI,X0)
277  ==NOERROR){
278  // Calculate momentum loss due to ionization
279  dP_dx = dPdx(mymom.Mag(), KrhoZ_overA, rhoZ_overA,LogI);
280  double my_step_size = 0.0001/fabs(dP_dx);
281 
282  if(my_step_size>MAX_STEP_SIZE)my_step_size=MAX_STEP_SIZE; // maximum step size in cm
283  if(my_step_size<MIN_STEP_SIZE)my_step_size=MIN_STEP_SIZE; // minimum step size in cm
284 
285  if (ploss_direction==kBackward) my_step_size*=-1.;
286  stepper.SetStepSize(my_step_size);
287  }
288  }
289  // Swim to next
290  double ds=stepper.Step(NULL);
291  s+=fabs(ds);
292  stepper.GetPosMom(mypos,mymom);
293 
294  // Take into account energy loss
295  double ptot=mymom.Mag();
296  // if (ploss_direction==kForward) ptot+=dP_dx*ds;
297  //else ptot-=dP_dx*ds;
298  ptot+=dP_dx*ds;
299 
300  mymom.SetMag(ptot);
301  stepper.SetStartingParams(q, &mypos, &mymom);
302 
303  // Break if we have passed the line to which we are trying to find the doca
304  DVector3 linepos=origin;
305  if (fabs(dir.z())>0.){
306  linepos+=((mypos.z()-origin.z())/dir.z())*dir;
307  }
308  doca=(linepos-mypos).Mag();
309  if (doca>old_doca) break;
310 
311  // Store the position and momentum for this step
312  last_pos=mypos;
313  last_mom=mymom;
314  }
315 }
316 
317 // Faster version of swimmer that only deals with the region inside the bore of
318 // the magnet. This is intended to be used with the hit selector, so is a
319 // paired-down version of the usual trajectory.
320 void DReferenceTrajectory::FastSwimForHitSelection(const DVector3 &pos, const DVector3 &mom, double q){
321  DMagneticFieldStepper stepper(bfield, q, &pos, &mom);
322  if(step_size>0.0)stepper.SetStepSize(step_size);
323 
324  // Step until we leave the active tracking region
325  swim_step_t *swim_step = this->swim_steps;
326  Nswim_steps = 0;
327  double itheta02 = 0.0;
328  double itheta02s = 0.0;
329  double itheta02s2 = 0.0;
330  double X0sum=0.0;
331  swim_step_t *last_step=NULL;
332 
333  for(double s=0; fabs(s)<1000.; Nswim_steps++, swim_step++){
334 
335  if(Nswim_steps>=this->max_swim_steps){
336  if (debug_level>0){
337  jerr<<__FILE__<<":"<<__LINE__<<" Too many steps in trajectory. Truncating..."<<endl;
338  }
339  break;
340  }
341 
342  stepper.GetDirs(swim_step->sdir, swim_step->tdir, swim_step->udir);
343  stepper.GetPosMom(swim_step->origin, swim_step->mom);
344  stepper.GetBField(swim_step->B);
345  swim_step->Ro = stepper.GetRo();
346  swim_step->s = s;
347  swim_step->t = 0.; // we do not need the time for our purpose...
348 
349  double Rsq=swim_step->origin.Perp2();
350  double z=swim_step->origin.Z();
351  if (z>345. || z<0. || Rsq>60.*60.){
352  Nswim_steps++; break;
353  }
354 
355  //magnitude of momentum and beta
356  double p_sq=swim_step->mom.Mag2();
357  double one_over_beta_sq=1.+mass_sq/p_sq;
358 
359  // Add material if geom is not NULL
360  double dP = 0.0;
361  double dP_dx=0.0;
362  if(geom){
363  double KrhoZ_overA=0.0;
364  double rhoZ_overA=0.0;
365  double LogI=0.0;
366  double X0=0.0;
367  jerror_t err = geom->FindMatALT1(swim_step->origin, swim_step->mom, KrhoZ_overA, rhoZ_overA,LogI, X0);
368  if(err == NOERROR){
369  if(X0>0.0){
370  double p=sqrt(p_sq);
371  double delta_s = s;
372  if(last_step)delta_s -= last_step->s;
373  double radlen = delta_s/X0;
374 
375  if(radlen>1.0E-5){ // PDG 2008 pg 271, second to last paragraph
376 
377  // double theta0 = 0.0136*sqrt(one_over_beta_sq)/p*sqrt(radlen)*(1.0+0.038*log(radlen)); // From PDG 2008 eq 27.12
378  //double theta02 = theta0*theta0;
379  double factor=1.0+0.038*log(radlen);
380  double theta02=1.8496e-4*factor*factor*radlen*one_over_beta_sq/p_sq;
381 
382  itheta02 += theta02;
383  itheta02s += s*theta02;
384  itheta02s2 += s*s*theta02;
385  X0sum+=X0;
386  }
387 
388  // Calculate momentum loss due to ionization
389  dP_dx = dPdx(p, KrhoZ_overA, rhoZ_overA,LogI);
390  }
391  }
392  last_step = swim_step;
393  }
394  swim_step->itheta02 = itheta02;
395  swim_step->itheta02s = itheta02s;
396  swim_step->itheta02s2 = itheta02s2;
397  swim_step->invX0=Nswim_steps/X0sum;
398 
399  if(step_size<0.0){ // step_size<0 indicates auto-calculated step size
400  // Adjust step size to take smaller steps in regions of high momentum loss
401  double my_step_size = 0.0001/fabs(dP_dx);
402  if(my_step_size>MAX_STEP_SIZE)my_step_size=MAX_STEP_SIZE; // maximum step size in cm
403  if(my_step_size<MIN_STEP_SIZE)my_step_size=MIN_STEP_SIZE; // minimum step size in cm
404 
405  stepper.SetStepSize(my_step_size);
406  }
407 
408  // Swim to next
409  double ds=stepper.FastStep();
410 
411  // Calculate momentum loss due to the step we're about to take
412  dP = ds*dP_dx;
413 
414  // Adjust momentum due to ionization losses
415  if(dP!=0.0){
416  DVector3 pos, mom;
417  stepper.GetPosMom(pos, mom);
418  double ptot = mom.Mag() - dP; // correct for energy loss
419  if (ptot<0.005) {Nswim_steps++; break;}
420  mom.SetMag(ptot);
421  stepper.SetStartingParams(q, &pos, &mom);
422  }
423  s += ds;
424  }
425 }
426 
427 // Faster version of the swimmer that uses an alternate stepper and does not
428 // check for material boundaries.
429 void DReferenceTrajectory::FastSwim(const DVector3 &pos, const DVector3 &mom, double q,double smax, double zmin,double zmax){
430 
431  /// (Re)Swim the trajectory starting from pos with momentum mom.
432  /// This will use the charge and step size (if given) passed to
433  /// the constructor when the object was created. It will also
434  /// (re)use the swim_step buffer, replacing it's contents.
435 
436  // If the charged passed to us is greater that 10, it means use the charge
437  // already stored in the class. Otherwise, use what was passed to us.
438  if(fabs(q)>10)
439  q = this->q;
440  else
441  this->q = q;
442 
443  DMagneticFieldStepper stepper(bfield, q, &pos, &mom);
444  if(step_size>0.0)stepper.SetStepSize(step_size);
445 
446  // Step until we hit a boundary (don't track more than 20 meters)
447  swim_step_t *swim_step = this->swim_steps;
448  double t=0.;
449  Nswim_steps = 0;
450  double itheta02 = 0.0;
451  double itheta02s = 0.0;
452  double itheta02s2 = 0.0;
453  double X0sum=0.0;
454  swim_step_t *last_step=NULL;
455  double old_radius_sq=1e6;
456 
457  // Variables used to tag the step at which the track passes into one one of
458  // the outer detectors
459  index_at_bcal=-1;
460  index_at_tof=-1;
461  index_at_fcal=-1;
462  bool hit_bcal=false,hit_fcal=false,hit_tof=false;
463 
464  for(double s=0; fabs(s)<smax; Nswim_steps++, swim_step++){
465 
466  if(Nswim_steps>=this->max_swim_steps){
467  if (debug_level>0){
468  jerr<<__FILE__<<":"<<__LINE__<<" Too many steps in trajectory. Truncating..."<<endl;
469  }
470  break;
471  }
472 
473  stepper.GetDirs(swim_step->sdir, swim_step->tdir, swim_step->udir);
474  stepper.GetPosMom(swim_step->origin, swim_step->mom);
475  stepper.GetBField(swim_step->B);
476  swim_step->Ro = stepper.GetRo();
477  swim_step->s = s;
478  swim_step->t = t;
479 
480  //magnitude of momentum and beta
481  double p_sq=swim_step->mom.Mag2();
482  double one_over_beta_sq=1.+mass_sq/p_sq;
483 
484  // Add material if geom or RootGeom is not NULL
485  // If both are non-NULL, then use RootGeom
486  double dP = 0.0;
487  double dP_dx=0.0;
488  if(RootGeom || geom){
489  double KrhoZ_overA=0.0;
490  double rhoZ_overA=0.0;
491  double LogI=0.0;
492  double X0=0.0;
493  jerror_t err;
494  if(RootGeom){
495  double rhoZ_overA,rhoZ_overA_logI;
496  err = RootGeom->FindMatLL(swim_step->origin,
497  rhoZ_overA,
498  rhoZ_overA_logI,
499  X0);
500  KrhoZ_overA=0.1535e-3*rhoZ_overA;
501  LogI=rhoZ_overA_logI/rhoZ_overA;
502  }else{
503  err = geom->FindMatALT1(swim_step->origin, swim_step->mom, KrhoZ_overA, rhoZ_overA,LogI, X0);
504  }
505  if(err == NOERROR){
506  if(X0>0.0){
507  double p=sqrt(p_sq);
508  double delta_s = s;
509  if(last_step)delta_s -= last_step->s;
510  double radlen = delta_s/X0;
511 
512  if(radlen>1.0E-5){ // PDG 2008 pg 271, second to last paragraph
513 
514  // double theta0 = 0.0136*sqrt(one_over_beta_sq)/p*sqrt(radlen)*(1.0+0.038*log(radlen)); // From PDG 2008 eq 27.12
515  //double theta02 = theta0*theta0;
516  double factor=1.0+0.038*log(radlen);
517  double theta02=1.8496e-4*factor*factor*radlen*one_over_beta_sq/p_sq;
518 
519  itheta02 += theta02;
520  itheta02s += s*theta02;
521  itheta02s2 += s*s*theta02;
522  X0sum+=X0;
523  }
524 
525  // Calculate momentum loss due to ionization
526  dP_dx = dPdx(p, KrhoZ_overA, rhoZ_overA,LogI);
527  }
528  }
529  last_step = swim_step;
530  }
531  swim_step->itheta02 = itheta02;
532  swim_step->itheta02s = itheta02s;
533  swim_step->itheta02s2 = itheta02s2;
534  swim_step->invX0=Nswim_steps/X0sum;
535 
536 
537  if(step_size<0.0){ // step_size<0 indicates auto-calculated step size
538  // Adjust step size to take smaller steps in regions of high momentum loss
539  double my_step_size = 0.0001/fabs(dP_dx);
540  if(my_step_size>MAX_STEP_SIZE)my_step_size=MAX_STEP_SIZE; // maximum step size in cm
541  if(my_step_size<MIN_STEP_SIZE)my_step_size=MIN_STEP_SIZE; // minimum step size in cm
542 
543  stepper.SetStepSize(my_step_size);
544  }
545 
546  // Swim to next
547  double ds=stepper.FastStep();
548 
549  // Calculate momentum loss due to the step we're about to take
550  dP = ds*dP_dx;
551 
552  // Adjust momentum due to ionization losses
553  if(dP!=0.0){
554  DVector3 pos, mom;
555  stepper.GetPosMom(pos, mom);
556  double ptot = mom.Mag() - dP; // correct for energy loss
557  if (ptot<0.005) {Nswim_steps++; break;}
558  mom.SetMag(ptot);
559  stepper.SetStartingParams(q, &pos, &mom);
560  }
561 
562  // update flight time
563  t+=ds*sqrt(one_over_beta_sq)/SPEED_OF_LIGHT;
564  s += ds;
565 
566  // Mark places along trajectory where we pass into one of the
567  // main detectors
568  double Rsq=swim_step->origin.Perp2();
569  double z=swim_step->origin.Z();
570  if (hit_bcal==false && Rsq>Rsqmax_interior && z<407 &&z>0){
571  index_at_bcal=Nswim_steps-1;
572  hit_bcal=true;
573  }
574  if (hit_tof==false && z>606.){
575  index_at_tof=Nswim_steps-1;
576  hit_tof=true;
577  }
578  if (hit_fcal==false && z>625.){
579  index_at_fcal=Nswim_steps-1;
580  hit_fcal=true;
581  }
582 
583 
584  // Exit the loop if we are already inside the volume of the BCAL
585  // and the radius is decreasing
586  if (Rsq<old_radius_sq && Rsq>Rsqmax_interior && z<407.0 && z>-100.0){
587  Nswim_steps++; break;
588  }
589 
590  // Exit loop if we leave the tracking volume
591  if (z>zmax){Nswim_steps++; break;}
592  if(Rsq>Rsqmax_exterior && z<407.0){Nswim_steps++; break;} // ran into BCAL
593  if (fabs(swim_step->origin.X())>160.0
594  || fabs(swim_step->origin.Y())>160.0)
595  {Nswim_steps++; break;} // left extent of TOF + 31cm Buffer
596  if(z>670.0){Nswim_steps++; break;} // ran into FCAL
597  if(z<zmin){Nswim_steps++; break;} // exit upstream
598 
599  old_radius_sq=Rsq;
600  }
601 
602  // OK. At this point the positions of the trajectory in the lab
603  // frame have been recorded along with the momentum of the
604  // particle and the directions of reference trajectory
605  // coordinate system at each point.
606 }
607 
608 //---------------------------------
609 // Swim
610 //---------------------------------
611 void DReferenceTrajectory::Swim(const DVector3 &pos, const DVector3 &mom, double q, const TMatrixFSym *cov,double smax, const DCoordinateSystem *wire)
612 {
613  /// (Re)Swim the trajectory starting from pos with momentum mom.
614  /// This will use the charge and step size (if given) passed to
615  /// the constructor when the object was created. It will also
616  /// (re)use the sim_step buffer, replacing it's contents.
617 
618  // If the charged passed to us is greater that 10, it means use the charge
619  // already stored in the class. Otherwise, use what was passed to us.
620  if(fabs(q)>10)
621  q = this->q;
622  else
623  this->q = q;
624 
625  DMagneticFieldStepper stepper(bfield, q, &pos, &mom);
626  if(step_size>0.0)stepper.SetStepSize(step_size);
627 
628  // Step until we hit a boundary (don't track more than 20 meters)
629  swim_step_t *swim_step = this->swim_steps;
630  double t=0.;
631  Nswim_steps = 0;
632  double itheta02 = 0.0;
633  double itheta02s = 0.0;
634  double itheta02s2 = 0.0;
635  double X0sum=0.0;
636  swim_step_t *last_step=NULL;
637  double old_radius_sq=1e6;
638 
639  TMatrixFSym mycov(7);
640  if (cov!=NULL){
641  mycov=*cov;
642  }
643 
644  // Reset flag indicating whether we hit the CDC endplate
645  // and get the parameters of the endplate so we can check
646  // if we hit it while swimming.
647  //hit_cdc_endplate = false;
648  /*
649 #if 0 // The GetCDCEndplate call goes all the way back to the XML and slows down
650  // overall tracking by a factor of 20. Therefore, we skip finding it
651  // and just hard-code the values instead. 1/28/2011 DL
652  double cdc_endplate_z=150+17; // roughly, from memory
653  double cdc_endplate_dz=5.0; // roughly, from memory
654  double cdc_endplate_rmin=10.0; // roughly, from memory
655  double cdc_endplate_rmax=55.0; // roughly, from memory
656  if(geom)geom->GetCDCEndplate(cdc_endplate_z, cdc_endplate_dz, cdc_endplate_rmin, cdc_endplate_rmax);
657  double cdc_endplate_zmin = cdc_endplate_z - cdc_endplate_dz/2.0;
658  double cdc_endplate_zmax = cdc_endplate_zmin + cdc_endplate_dz;
659 #else
660  double cdc_endplate_rmin=10.0; // roughly, from memory
661  double cdc_endplate_rmax=55.0; // roughly, from memory
662  double cdc_endplate_zmin = 167.6;
663  double cdc_endplate_zmax = 168.2;
664 #endif
665  */
666 
667 #if 0
668  // Get Bfield from stepper to initialize Bz_old
669  DVector3 B;
670  stepper.GetBField(B);
671  double Bz_old = B.z();
672 #endif
673 
674  // Variables used to tag the step at which the track passes into one
675  // one of the outer detectors
676  index_at_bcal=-1;
677  index_at_tof=-1;
678  index_at_fcal=-1;
679  bool hit_bcal=false,hit_fcal=false,hit_tof=false;
680 
681  for(double s=0; fabs(s)<smax; Nswim_steps++, swim_step++){
682 
683  if(Nswim_steps>=this->max_swim_steps){
684  if (debug_level>0){
685  jerr<<__FILE__<<":"<<__LINE__<<" Too many steps in trajectory. Truncating..."<<endl;
686  }
687  break;
688  }
689 
690  stepper.GetDirs(swim_step->sdir, swim_step->tdir, swim_step->udir);
691  stepper.GetPosMom(swim_step->origin, swim_step->mom);
692  swim_step->Ro = stepper.GetRo();
693  swim_step->s = s;
694  swim_step->t = t;
695 
696  //magnitude of momentum and beta
697  double p_sq=swim_step->mom.Mag2();
698  double one_over_beta_sq=1.+mass_sq/p_sq;
699 
700  // Add material if geom or RootGeom is not NULL
701  // If both are non-NULL, then use RootGeom
702  double dP = 0.0;
703  double dP_dx=0.0;
704  double s_to_boundary=1.0E6; // initialize to "infinity" in case we don't set this below
705  if(RootGeom || geom){
706  double KrhoZ_overA=0.0;
707  double rhoZ_overA=0.0;
708  double LogI=0.0;
709  double X0=0.0;
710  jerror_t err;
711  if(RootGeom){
712  double rhoZ_overA,rhoZ_overA_logI;
713  err = RootGeom->FindMatLL(swim_step->origin,
714  rhoZ_overA,
715  rhoZ_overA_logI,
716  X0);
717  KrhoZ_overA=0.1535e-3*rhoZ_overA;
718  LogI=rhoZ_overA_logI/rhoZ_overA;
719  }else{
720  if(check_material_boundaries){
721  err = geom->FindMatALT1(swim_step->origin, swim_step->mom, KrhoZ_overA, rhoZ_overA,LogI, X0, &s_to_boundary);
722  }else{
723  err = geom->FindMatALT1(swim_step->origin, swim_step->mom, KrhoZ_overA, rhoZ_overA,LogI, X0);
724  }
725 
726  // Check if we hit the CDC endplate
727  //double z = swim_step->origin.Z();
728  //if(z>=cdc_endplate_zmin && z<=cdc_endplate_zmax){
729  // double r = swim_step->origin.Perp();
730  // if(r>=cdc_endplate_rmin && r<=cdc_endplate_rmax){
731  // hit_cdc_endplate = true;
732  //}
733  //}
734  }
735 
736  if(err == NOERROR){
737  if(X0>0.0){
738  double p=sqrt(p_sq);
739  double delta_s = s;
740  if(last_step)delta_s -= last_step->s;
741  double radlen = delta_s/X0;
742 
743  if(radlen>1.0E-5){ // PDG 2008 pg 271, second to last paragraph
744 
745  // double theta0 = 0.0136*sqrt(one_over_beta_sq)/p*sqrt(radlen)*(1.0+0.038*log(radlen)); // From PDG 2008 eq 27.12
746  //double theta02 = theta0*theta0;
747  double factor=1.0+0.038*log(radlen);
748  double theta02=1.8496e-4*factor*factor*radlen*one_over_beta_sq/p_sq;
749 
750  itheta02 += theta02;
751  itheta02s += s*theta02;
752  itheta02s2 += s*s*theta02;
753  X0sum+=X0;
754 
755  if (cov){
756 
757  }
758  }
759 
760  // Calculate momentum loss due to ionization
761  dP_dx = dPdx(p, KrhoZ_overA, rhoZ_overA,LogI);
762  }
763  }
764  last_step = swim_step;
765  }
766  swim_step->itheta02 = itheta02;
767  swim_step->itheta02s = itheta02s;
768  swim_step->itheta02s2 = itheta02s2;
769  swim_step->invX0=Nswim_steps/X0sum;
770 
771  // Adjust step size to take smaller steps in regions of high momentum loss or field gradient
772  if(step_size<0.0){ // step_size<0 indicates auto-calculated step size
773  // Take step so as to change momentum by 100keV
774  //double my_step_size=p/fabs(dP_dx)*0.01;
775  double my_step_size = 0.0001/fabs(dP_dx);
776 
777  // Now check the field gradient
778 #if 0
779  stepper.GetBField(B);
780  double Bz = B.z();
781  if (fabs(Bz-Bz_old)>EPS){
782  double my_step_size_B=0.01*my_step_size
783  *fabs(Bz/(Bz_old-Bz));
784  if (my_step_size_B<my_step_size)
785  my_step_size=my_step_size_B;
786  }
787  Bz_old=Bz; // Save old z-component of B-field
788 #endif
789  // Use the estimated distance to the boundary to make sure we don't overstep
790  // into a high density region and miss some material. Use half the estimated
791  // distance since it's only an estimate. Note that even though this would lead
792  // to infinitely small steps, there is a minimum step size imposed below to
793  // ensure the step size is reasonable.
794  /*
795  double step_size_to_boundary = BOUNDARY_STEP_FRACTION*s_to_boundary;
796  if(step_size_to_boundary < my_step_size)my_step_size = step_size_to_boundary;
797  */
798 
799  if(my_step_size>MAX_STEP_SIZE)my_step_size=MAX_STEP_SIZE; // maximum step size in cm
800  if(my_step_size<MIN_STEP_SIZE)my_step_size=MIN_STEP_SIZE; // minimum step size in cm
801 
802  stepper.SetStepSize(my_step_size);
803  }
804 
805  // Swim to next
806  double ds=stepper.Step(NULL,&swim_step->B);
807  if (cov){
808  PropagateCovariance(ds,q,mass_sq,mom,pos,swim_step->B,mycov);
809  swim_step->cov_t_t=mycov(6,6);
810  swim_step->cov_px_t=mycov(6,0);
811  swim_step->cov_py_t=mycov(6,1);
812  swim_step->cov_pz_t=mycov(6,2);
813  }
814 
815  // Calculate momentum loss due to the step we're about to take
816  dP = ds*dP_dx;
817 
818  // Adjust momentum due to ionization losses
819  if(dP!=0.0){
820  DVector3 pos, mom;
821  stepper.GetPosMom(pos, mom);
822  double ptot = mom.Mag() - dP; // correct for energy loss
823  bool ranged_out = false;
824  /*
825  if (ptot<0.05){
826  swim_step->origin.Print();
827  cout<<"N: " << Nswim_steps <<" x " << pos.x() <<" y " <<pos.y() <<" z " << pos.z() <<" r " << pos.Perp()<< " s " << s << " p " << ptot << endl;
828  }
829  */
830  if(ptot<0.005)ranged_out=true;
831  if(dP<0.0 && ploss_direction==kForward)ranged_out=true;
832  if(dP>0.0 && ploss_direction==kBackward)ranged_out=true;
833  //if(mom.Mag()==0.0)ranged_out=true;
834  if(ranged_out){
835  Nswim_steps++; // This will at least allow for very low momentum particles to have 1 swim step
836  break;
837  }
838  mom.SetMag(ptot);
839  stepper.SetStartingParams(q, &pos, &mom);
840  }
841 
842  // update flight time
843  t+=ds*sqrt(one_over_beta_sq)/SPEED_OF_LIGHT;
844  s += ds;
845 
846  // Mark places along trajectory where we pass into one of the
847  // main detectors
848  double Rsq=swim_step->origin.Perp2();
849  double z=swim_step->origin.Z();
850  if (hit_bcal==false && Rsq>Rsqmax_interior && z<407 &&z>0){
851  index_at_bcal=Nswim_steps-1;
852  hit_bcal=true;
853  }
854  if (hit_tof==false && z>618.){
855  index_at_tof=Nswim_steps-1;
856  hit_tof=true;
857  }
858  if (hit_fcal==false && z>625.){
859  index_at_fcal=Nswim_steps-1;
860  hit_fcal=true;
861  }
862 
863  // Exit the loop if we are already inside the volume of the BCAL
864  // and the radius is decreasing
865  if (Rsq<old_radius_sq && Rsq>Rsqmax_interior && z<407.0 && z>-100.0){
866  Nswim_steps++; break;
867  }
868 
869 
870  // Exit loop if we leave the tracking volume
871  if(Rsq>Rsqmax_exterior && z<407.0){Nswim_steps++; break;} // ran into BCAL
872  if (fabs(swim_step->origin.X())>160.0
873  || fabs(swim_step->origin.Y())>160.0)
874  {Nswim_steps++; break;} // left extent of TOF + 31cm Buffer
875  if(z>zmax_track_boundary){Nswim_steps++; break;} // ran into FCAL
876  if(z<zmin_track_boundary){Nswim_steps++; break;} // exit upstream
877  if(wire && Nswim_steps>0){ // optionally check if we passed a wire we're supposed to be swimming to
878  swim_step_t *closest_step = FindClosestSwimStep(wire);
879  if(++closest_step!=swim_step){Nswim_steps++; break;}
880  }
881 
882  old_radius_sq=Rsq;
883  }
884 
885  // OK. At this point the positions of the trajectory in the lab
886  // frame have been recorded along with the momentum of the
887  // particle and the directions of reference trajectory
888  // coordinate system at each point.
889 }
890 
891 // Routine to find position on the trajectory where the track crosses a radial
892 // position R. Also returns the path length to this position.
894  DVector3 &mypos,
895  double *s,
896  double *t,
897  DVector3 *p_at_intersection) const{
898  mypos.SetXYZ(QuietNaN,QuietNaN,QuietNaN);
899  if(p_at_intersection)
900  p_at_intersection->SetXYZ(QuietNaN,QuietNaN,QuietNaN);
901 
902  if(Nswim_steps<1){
903  _DBG_<<"No swim steps! You must \"Swim\" the track before calling GetIntersectionWithRadius(...)"<<endl;
904  }
905  // Return early if the radius at the end of the reference trajectory is still less than R
906  double outer_radius=swim_steps[Nswim_steps-1].origin.Perp();
907  if (outer_radius<R){
908  if (s) *s=0.;
909  if (t) *t=0.;
910  return VALUE_OUT_OF_RANGE;
911  }
912  // Return early if the radius at the beginning of the trajectory is outside
913  // the radius to which we are trying to match
914  double inner_radius=swim_steps[0].origin.Perp();
915  if (inner_radius>R){
916  if (s) *s=0.;
917  if (t) *t=0.;
918  return VALUE_OUT_OF_RANGE;
919  }
920 
921  // Loop over swim steps and find the one that crosses the radius
922  swim_step_t *swim_step = swim_steps;
923  swim_step_t *step=NULL;
924  swim_step_t *last_step=NULL;
925 
926  // double inner_radius=swim_step->origin.Perp();
927  for(int i=0; i<Nswim_steps; i++, swim_step++){
928  if (swim_step->origin.Perp()>R){
929  step=swim_step;
930  break;
931  }
932  if (swim_step->origin.Z()>407.0) return VALUE_OUT_OF_RANGE;
933  last_step=swim_step;
934  }
935  if (step==NULL||last_step==NULL) return VALUE_OUT_OF_RANGE;
936  if (p_at_intersection!=NULL){
937  *p_at_intersection=last_step->mom;
938  }
939 
940  // At this point, the location where the track intersects the cyclinder
941  // is somewhere between last_step and step. For simplicity, we're going
942  // to just find the intersection of the cylinder with the line that joins
943  // the 2 positions. We do this by working in the X/Y plane only and
944  // finding the value of "alpha" which is the fractional distance the
945  // intersection point is between last_pos and mypos. We'll then apply
946  // the alpha found in the 2D X/Y space to the 3D x/y/Z space to find
947  // the actual intersection point.
948  DVector2 x1(last_step->origin.X(), last_step->origin.Y());
949  DVector2 x2(step->origin.X(), step->origin.Y());
950  DVector2 dx = x2-x1;
951  double A = dx.Mod2();
952  double B = 2.0*(x1.X()*dx.X() + x1.Y()*dx.Y());
953  double C = x1.Mod2() - R*R;
954 
955  double sqrt_D=sqrt(B*B-4.0*A*C);
956  double one_over_denom=0.5/A;
957  double alpha1 = (-B + sqrt_D)*one_over_denom;
958  double alpha2 = (-B - sqrt_D)*one_over_denom;
959  double alpha = alpha1;
960  if(alpha1<0.0 || alpha1>1.0)alpha=alpha2;
961  if(!isfinite(alpha))return VALUE_OUT_OF_RANGE;
962 
963  DVector3 delta = step->origin - last_step->origin;
964  mypos = last_step->origin + alpha*delta;
965 
966  // The value of s actually represents the pathlength
967  // to the outside point. Adjust it back to the
968  // intersection point (approximately).
969  if (s) *s = step->s-(1.0-alpha)*delta.Mag();
970 
971  // flight time
972  if (t){
973  double p_sq=step->mom.Mag2();
974  double one_over_beta=sqrt(1.+mass_sq/p_sq);
975  *t = step->t-(1.0-alpha)*delta.Mag()*one_over_beta/SPEED_OF_LIGHT;
976  }
977 
978  return NOERROR;
979 }
980 
981 //---------------------------------
982 // GetIntersectionWithPlane
983 //---------------------------------
984 jerror_t DReferenceTrajectory::GetIntersectionWithPlane(const DVector3 &origin, const DVector3 &norm, DVector3 &pos, double *s,double *t,double *var_t,DetectorSystem_t detector) const{
985  DVector3 dummy;
986  return GetIntersectionWithPlane(origin,norm,pos,dummy,s,t,var_t,detector);
987 }
988 jerror_t DReferenceTrajectory::GetIntersectionWithPlane(const DVector3 &origin, const DVector3 &norm, DVector3 &pos, DVector3 &p_at_intersection, double *s,
989  double *t,double *var_t,DetectorSystem_t detector) const
990 {
991  /// Get the intersection point of this trajectory with a plane.
992  /// The plane is specified by <i>origin</i> and <i>norm</i>. The
993  /// <i>origin</i> vector should give the coordinates of any point
994  /// on the plane and <i>norm</i> should give a vector normal to
995  /// the plane. The <i>norm</i> vector will be copied and normalized
996  /// so it can be of any magnitude upon entry.
997  ///
998  /// The coordinates of the intersection point will copied into
999  /// the supplied <i>pos</i> vector. If a non-NULL pointer for <i>s</i>
1000  /// is passed in, the pathlength of the trajectory from its begining
1001  /// to the intersection point is copied into location pointed to.
1002 
1003  // Set reasonable defaults
1004  pos.SetXYZ(0,0,0);
1005  if(s)*s=0.0;
1006  if(t)*t=0.0;
1007  p_at_intersection.SetXYZ(0,0,0);
1008 
1009  // Return early if the z-position of the plane with which we are
1010  // intersecting is beyong the reference trajectory.
1011  if (origin.z()>swim_steps[Nswim_steps-1].origin.z()){
1012  return VALUE_OUT_OF_RANGE;
1013  }
1014  // Find the closest swim step to the position where the track crosses
1015  // the plane
1016  int first_i=0;
1017  switch(detector){
1018  case SYS_FCAL:
1019  if (index_at_fcal<0) return VALUE_OUT_OF_RANGE;
1020  first_i=index_at_fcal;
1021  break;
1022  case SYS_TOF:
1023  if (index_at_tof<0) return VALUE_OUT_OF_RANGE;
1024  first_i=index_at_tof;
1025  break;
1026  default:
1027  break;
1028  }
1029  swim_step_t *step=FindPlaneCrossing(origin,norm,first_i,detector);
1030  if(!step){
1031  _DBG_<<"Could not find closest swim step!"<<endl;
1032  return RESOURCE_UNAVAILABLE;
1033  }
1034 
1035  // Here we follow a scheme described in more detail in the
1036  // DistToRT(DVector3 hit) method below. The basic idea is to
1037  // express a point on the helix in terms of a single variable
1038  // and then solve for that variable by setting the distance
1039  // to zero.
1040  //
1041  // x = Ro*(cos(phi) - 1)
1042  // y = Ro*sin(phi)
1043  // z = phi*(dz/dphi)
1044  //
1045  // As is done below, we work in the RT coordinate system. Well,
1046  // sort of. The distance to the plane is given by:
1047  //
1048  // d = ( x - c ).n
1049  //
1050  // where x is a point on the helix, c is the "origin" point
1051  // which lies somewhere in the plane and n is the "norm"
1052  // vector. Since we want a point in the plane, we set d=0
1053  // and solve for phi (with the components of x expressed in
1054  // terms of phi as given in the DistToRT method below).
1055  //
1056  // Thus, the equation we need to solve is:
1057  //
1058  // x.n - c.n = 0
1059  //
1060  // notice that "c" only gets dotted into "n" so that
1061  // dot product can occur in any coordinate system. Therefore,
1062  // we do that in the lab coordinate system to avoid the
1063  // overhead of transforming "c" to the RT system. The "n"
1064  // vector, however, still must be transformed.
1065  //
1066  // Expanding the trig functions to 2nd order in phi, performing
1067  // the x.n dot product, and gathering equal powers of phi
1068  // leads us to he following:
1069  //
1070  // (-Ro*nx/2)*phi^2 + (Ro*ny+dz_dphi*nz)*phi - c.n = 0
1071  //
1072  // which is quadratic in phi. We solve for both roots, but use
1073  // the one with the smller absolute value (if both are finite).
1074 
1075  double &Ro = step->Ro;
1076 
1077  // OK, having said all of that, it turns out that the above
1078  // mechanism will tend to fail in regions of low or no
1079  // field because the value of Ro is very large. Thus, we need to
1080  // use a straight line projection in such cases. We also
1081  // want to use a straight line projection if the helical intersection
1082  // fails for some other reason.
1083  //
1084  // The algorthim is then to only try the helical calculation
1085  // for small (<10m) values of Ro and then do the straight line
1086  // if R is larger than that OR the helical calculation fails.
1087 
1088  // Try helical calculation
1089  if(Ro<1000.0){
1090  double nx = norm.Dot(step->sdir);
1091  double ny = norm.Dot(step->tdir);
1092  double nz = norm.Dot(step->udir);
1093 
1094  double delta_z = step->mom.Dot(step->udir);
1095  double delta_phi = step->mom.Dot(step->tdir)/Ro;
1096  double dz_dphi = delta_z/delta_phi;
1097 
1098  double A = -Ro*nx/2.0;
1099  double B = Ro*ny + dz_dphi*nz;
1100  double C = norm.Dot(step->origin-origin);
1101  double sqroot=sqrt(B*B-4.0*A*C);
1102  double twoA=2.0*A;
1103 
1104  double phi_1 = (-B + sqroot)/(twoA);
1105  double phi_2 = (-B - sqroot)/(twoA);
1106 
1107  double phi = fabs(phi_1)<fabs(phi_2) ? phi_1:phi_2;
1108  if(!isfinite(phi_1))phi = phi_2;
1109  if(!isfinite(phi_2))phi = phi_1;
1110  if(isfinite(phi)){
1111 
1112  double my_s = -Ro/2.0 * phi*phi;
1113  double my_t = Ro * phi;
1114  double my_u = dz_dphi * phi;
1115 
1116  pos = step->origin + my_s*step->sdir + my_t*step->tdir + my_u*step->udir;
1117  p_at_intersection = step->mom;
1118  if(s){
1119  double delta_s = sqrt(my_t*my_t + my_u*my_u);
1120  *s = step->s + (phi>0 ? +delta_s:-delta_s);
1121  }
1122  // flight time
1123  if (t){
1124  double delta_s = sqrt(my_t*my_t + my_u*my_u);
1125  double ds=(phi>0 ? +delta_s:-delta_s);
1126  double p_sq=step->mom.Mag2();
1127  double one_over_beta=sqrt(1.+mass_sq/p_sq);
1128  *t = step->t+ds*one_over_beta/SPEED_OF_LIGHT;
1129  }
1130  if (var_t){
1131  *var_t=step->cov_t_t;
1132  }
1133 
1134  // Success. Go ahead and return
1135  return NOERROR;
1136  }
1137  }
1138 
1139  // If we got here then we need to try a straight line calculation
1140  double p_sq=step->mom.Mag2();
1141  double dz_over_pz=(origin.z()-step->origin.z())/step->mom.z();
1142  double ds=sqrt(p_sq)*dz_over_pz;
1143  pos.SetXYZ(step->origin.x()+dz_over_pz*step->mom.x(),
1144  step->origin.y()+dz_over_pz*step->mom.y(),
1145  origin.z());
1146  p_at_intersection = step->mom;
1147 
1148  if(s){
1149  *s = step->s + ds;
1150  }
1151  // flight time
1152  if (t){
1153  double one_over_beta=sqrt(1.+mass_sq/p_sq);
1154  *t = step->t+ds*one_over_beta/SPEED_OF_LIGHT;
1155  }
1156  // Flight time variance
1157  if (var_t){
1158  *var_t=step->cov_t_t;
1159  }
1160 
1161  return NOERROR;
1162 }
1163 
1164 //---------------------------------
1165 // InsertSteps
1166 //---------------------------------
1167 int DReferenceTrajectory::InsertSteps(const swim_step_t *start_step, double delta_s, double step_size)
1168 {
1169  /// Insert additional steps into the reference trajectory starting
1170  /// at start_step and swimming for at least delta_s by step_size
1171  /// sized steps. Both delta_s and step_size are in centimeters.
1172  /// If the value of delta_s is negative then the particle's momentum
1173  /// and charge are reversed before swimming. This could be a problem
1174  /// if energy loss is implemented.
1175 
1176  if(!start_step)return -1;
1177 
1178  // We do this by creating another, temporary DReferenceTrajectory object
1179  // on the stack and swimming it.
1180  DVector3 pos = start_step->origin;
1181  DVector3 mom = start_step->mom;
1182  double my_q = q;
1183  int direction = +1;
1184  if(delta_s<0.0){
1185  mom *= -1.0;
1186  my_q = -q;
1187  direction = -1;
1188  }
1189 
1190  // Here I allocate the steps using an auto_ptr so I don't have to mess with
1191  // deleting them at all of the possible exits. The problem with auto_ptr
1192  // is it can't handle arrays so it has to be wrapped in a struct.
1193  auto_ptr<StepStruct> steps_aptr(new StepStruct);
1194  DReferenceTrajectory::swim_step_t *steps = steps_aptr->steps;
1195  DReferenceTrajectory rt(bfield , my_q , steps , 256);
1196  rt.SetStepSize(step_size);
1197  rt.Swim(pos, mom, my_q,NULL,fabs(delta_s));
1198  if(rt.Nswim_steps==0)return 1;
1199 
1200  // Check that there is enough space to add these points
1201  if((Nswim_steps+rt.Nswim_steps)>max_swim_steps){
1202  //_DBG_<<"Not enough swim steps available to add new ones! Max="<<max_swim_steps<<" had="<<Nswim_steps<<" new="<<rt.Nswim_steps<<endl;
1203  return 2;
1204  }
1205 
1206  // At this point, we may have swum forward or backwards so the points
1207  // will need to be added either before start_step or after it. We also
1208  // may want to replace an old step that overlaps our high density steps
1209  // since they are presumably more accurate. Find the indexes of the
1210  // existing steps that the new steps will be inserted between.
1211  double sdiff = rt.swim_steps[rt.Nswim_steps-1].s;
1212  double s1 = start_step->s;
1213  double s2 = start_step->s + (double)direction*sdiff;
1214  double smin = s1<s2 ? s1:s2;
1215  double smax = s1<s2 ? s2:s1;
1216  int istep_start = 0;
1217  int istep_end = 0;
1218  for(int i=0; i<Nswim_steps; i++){
1219  if(swim_steps[i].s < smin)istep_start = i;
1220  if(swim_steps[i].s <= smax)istep_end = i+1;
1221  }
1222 
1223  // istep_start and istep_end now point to the steps we want to keep.
1224  // All steps between them (exclusive) will be overwritten. Note that
1225  // the original start_step should be in the "overwrite" range since
1226  // it is included already in the new trajectory.
1227  int steps_to_overwrite = istep_end - istep_start - 1;
1228  int steps_to_shift = rt.Nswim_steps - steps_to_overwrite;
1229 
1230  // Shift the steps down (or is it up?) starting with istep_end.
1231  for(int i=Nswim_steps-1; i>=istep_end; i--)swim_steps[i+steps_to_shift] = swim_steps[i];
1232 
1233  // Copy the new steps into this object
1234  double s_0 = start_step->s;
1235  double itheta02_0 = start_step->itheta02;
1236  double itheta02s_0 = start_step->itheta02s;
1237  double itheta02s2_0 = start_step->itheta02s2;
1238  for(int i=0; i<rt.Nswim_steps; i++){
1239  int index = direction>0 ? (istep_start+1+i):(istep_start+1+rt.Nswim_steps-1-i);
1240  swim_steps[index] = rt.swim_steps[i];
1241  swim_steps[index].s = s_0 + (double)direction*swim_steps[index].s;
1242  swim_steps[index].itheta02 = itheta02_0 + (double)direction*swim_steps[index].itheta02;
1243  swim_steps[index].itheta02s = itheta02s_0 + (double)direction*swim_steps[index].itheta02s;
1244  swim_steps[index].itheta02s2 = itheta02s2_0 + (double)direction*swim_steps[index].itheta02s2;
1245  if(direction<0.0){
1246  swim_steps[index].sdir *= -1.0;
1247  swim_steps[index].tdir *= -1.0;
1248  }
1249  }
1250  Nswim_steps += rt.Nswim_steps-steps_to_overwrite;
1251 
1252  // Note that the above procedure may leave us with "kinks" in the itheta0
1253  // variables. It may be that we need to recalculate those for all of the
1254  // new points and the ones after them by making one more pass. I'm hoping
1255  // it is a realitively small correction though so we can skip it here.
1256  return 0;
1257 }
1258 
1259 //---------------------------------
1260 // DistToRTwithTime
1261 //---------------------------------
1262 double DReferenceTrajectory::DistToRTwithTime(DVector3 hit, double *s,double *t,
1263  double *var_t,
1264  DetectorSystem_t detector) const{
1265  double dist=DistToRT(hit,s,detector);
1266  if (s!=NULL && t!=NULL)
1267  {
1268  if(last_swim_step==NULL)
1269  {
1270  *s = 1.0E6;
1271  *t = 1.0E6;
1272  if (var_t!=NULL){
1273  *var_t=1.0E6;
1274  }
1275  }
1276  else
1277  {
1278  double p_sq=last_swim_step->mom.Mag2();
1279  double one_over_beta=sqrt(1.+mass_sq/p_sq);
1280  *t=last_swim_step->t+(*s-last_swim_step->s)*one_over_beta/SPEED_OF_LIGHT;
1281  if (var_t!=NULL){
1282  *var_t=last_swim_step->cov_t_t;
1283  }
1284  }
1285  }
1286  return dist;
1287 }
1288 
1289 //---------------------------------
1290 // DistToRT
1291 //---------------------------------
1293  DetectorSystem_t detector) const
1294 {
1295  last_swim_step=NULL;
1296  if(Nswim_steps<1)_DBG__;
1297 
1298  int start_index=0;
1299  switch(detector){
1300  case SYS_BCAL:
1301  if (index_at_bcal<0) return numeric_limits<double>::quiet_NaN();
1302  start_index=index_at_bcal;
1303  break;
1304  case SYS_FCAL:
1305  if (index_at_fcal<0) return numeric_limits<double>::quiet_NaN();
1306  start_index=index_at_fcal;
1307  break;
1308  case SYS_TOF:
1309  if (index_at_tof<0) return numeric_limits<double>::quiet_NaN();
1310  start_index=index_at_tof;
1311  break;
1312  default:
1313  break;
1314  }
1315 
1316  // First, find closest step to point
1317  swim_step_t *swim_step = &swim_steps[start_index];
1318  swim_step_t *step=NULL;
1319 
1320  //double min_delta2 = 1.0E6;
1321  double old_delta2=10.e6,delta2=1.0e6;
1322 
1323  // Check if we should start at the end of the reference trajectory
1324  // or the beginning...
1325  int last_index=Nswim_steps-1;
1326  double forward_delta2=(swim_step->origin - hit).Mag2();
1327  double backward_delta2=(swim_steps[last_index].origin-hit).Mag2();
1328 
1329  if (forward_delta2<backward_delta2){ // start at the beginning
1330  for(int i=start_index; i<Nswim_steps; i++, swim_step++){
1331 
1332  DVector3 pos_diff = swim_step->origin - hit;
1333  delta2 = pos_diff.Mag2();
1334 
1335  if (delta2>old_delta2){
1336  break;
1337  }
1338 
1339  //if(delta2 < min_delta2){
1340  //min_delta2 = delta2;
1341 
1342  step = swim_step;
1343  old_delta2=delta2;
1344  //}
1345  }
1346  }
1347  else{// start at the end
1348  for(int i=last_index; i>=start_index; i--){
1349  swim_step=&swim_steps[i];
1350  DVector3 pos_diff = swim_step->origin - hit;
1351  delta2 = pos_diff.Mag2();
1352  if (delta2>old_delta2) break;
1353 
1354  //if(delta2 < min_delta2){
1355  //min_delta2 = delta2;
1356 
1357  step = swim_step;
1358  old_delta2=delta2;
1359  //}
1360  }
1361 
1362  }
1363 
1364  if(step==NULL){
1365  // It seems to occasionally occur that we have 1 swim step
1366  // and it's values are invalid. Supress warning messages
1367  // for these as they are "known" (even if not fully understood!)
1368  if(s != NULL)
1369  *s = 1.0E6;
1370  if(Nswim_steps>1){
1371  _DBG_<<"\"hit\" passed to DistToRT(DVector3) out of range!"<<endl;
1372  _DBG_<<"hit x,y,z = "<<hit.x()<<", "<<hit.y()<<", "<<hit.z()<<" Nswim_steps="<<Nswim_steps<<" min_delta2="<<delta2<<endl;
1373  }
1374  return 1.0E6;
1375  }
1376 
1377  // store last step
1378  last_swim_step=step;
1379 
1380 
1381  // Next, define a point on the helical segment defined by the
1382  // swim step it the RT coordinate system. The directions of
1383  // the RT coordinate system are defined by step->xdir, step->ydir,
1384  // and step->zdir. The coordinates of a point on the helix
1385  // in this coordinate system are:
1386  //
1387  // x = Ro*(cos(phi) - 1)
1388  // y = Ro*sin(phi)
1389  // z = phi*(dz/dphi)
1390  //
1391  // where phi is the phi angle of the point in this coordinate system.
1392  // phi=0 corresponds to the swim step point itself
1393  //
1394  // Transform the given coordinates to the RT coordinate system
1395  // and call these x0,y0,z0. Then, the distance of point to a
1396  // point on the helical segment is given by:
1397  //
1398  // d^2 = (x0-x)^2 + (y0-y)^2 + (z0-z)^2
1399  //
1400  // where x,y,z are all functions of phi as given above.
1401  //
1402  // writing out d^2 in terms of phi, but using the small angle
1403  // approximation for the trig functions, an equation for the
1404  // distance in only phi is obtained. Taking the derivative
1405  // and setting it equal to zero leaves a 3rd order polynomial
1406  // in phi whose root corresponds to the minimum distance.
1407  // Skipping some math, this equation has the form:
1408  //
1409  // d(d^2)/dphi = 0 = Ro^2*phi^3 + 2*alpha*phi + beta
1410  //
1411  // where:
1412  // alpha = x0*Ro + Ro^2 + (dz/dphi)^2
1413  //
1414  // beta = -2*y0*Ro - 2*z0*(dz/dphi)
1415  //
1416  // The above 3rd order poly is convenient in that it does not
1417  // contain a phi^2 term. This means we can skip the step
1418  // done in the general case where a change of variables is
1419  // made such that the 2nd order term disappears.
1420  //
1421  // In general, an equation of the form
1422  //
1423  // w^3 + 3.0*b*w + 2*c = 0
1424  //
1425  // has one real root:
1426  //
1427  // w0 = q - p
1428  //
1429  // where:
1430  // q^3 = d - c
1431  // p^3 = d + c
1432  // d^2 = b^3 + c^2 (don't confuse with d^2 above!)
1433  //
1434  // So for us ...
1435  //
1436  // 3b = 2*alpha/(Ro^2)
1437  // 2c = beta/(Ro^2)
1438 
1439  hit -= step->origin;
1440  double x0 = hit.Dot(step->sdir);
1441  double y0 = hit.Dot(step->tdir);
1442  double z0 = hit.Dot(step->udir);
1443  double &Ro = step->Ro;
1444  double Ro2 = Ro*Ro;
1445  double delta_z = step->mom.Dot(step->udir);
1446  double delta_phi = step->mom.Dot(step->tdir)/Ro;
1447  double dz_dphi = delta_z/delta_phi;
1448 
1449  // double alpha = x0*Ro + Ro2 + pow(dz_dphi,2.0);
1450  double alpha=x0*Ro + Ro2 +dz_dphi*dz_dphi;
1451  // double beta = -2.0*y0*Ro - 2.0*z0*dz_dphi;
1452  double beta = -2.0*(y0*Ro + z0*dz_dphi);
1453  // double b = (2.0*alpha/Ro2)/3.0;
1454  double b= TWO_THIRD*alpha/Ro2;
1455  // double c = (beta/Ro2)/2.0;
1456  double c = 0.5*(beta/Ro2);
1457  // double d = sqrt(pow(b,3.0) + pow(c,2.0));
1458  double d2=b*b*b+c*c;
1459  double phi=0.,dist2=1e8;
1460  if (d2>=0){
1461  double d=sqrt(d2);
1462  //double q = pow(d-c, ONE_THIRD);
1463  //double p = pow(d+c, ONE_THIRD);
1464  double p=cbrt(d+c);
1465  double q=cbrt(d-c);
1466  phi = q - p;
1467  if (fabs(phi)<0.2){ // check small angle approximation
1468  double phisq=phi*phi;
1469 
1470  dist2 = 0.25*Ro2*phisq*phisq + alpha*phisq + beta*phi
1471  + x0*x0 + y0*y0 + z0*z0;
1472  }
1473  else{
1474  return numeric_limits<double>::quiet_NaN();
1475  }
1476  }
1477  else{
1478  // Use DeMoivre's theorem to find the cube root of a complex
1479  // number. In this case there are three real solutions.
1480  double d=sqrt(-d2);
1481  c*=-1.;
1482  double temp=sqrt(cbrt(c*c+d*d));
1483  double theta1=ONE_THIRD*atan2(d,c);
1484  double sum_over_2=temp*cos(theta1);
1485  double diff_over_2=-temp*sin(theta1);
1486 
1487  double phi0=2.*sum_over_2;
1488  double phi0sq=phi0*phi0;
1489  double phi1=-sum_over_2+sqrt(3)*diff_over_2;
1490  double phi1sq=phi1*phi1;
1491  double phi2=-sum_over_2-sqrt(3)*diff_over_2;
1492  double phi2sq=phi2*phi2;
1493  double d2_2 = 0.25*Ro2*phi2sq*phi2sq + alpha*phi2sq + beta*phi2 + x0*x0 + y0*y0 + z0*z0;
1494  double d2_1 = 0.25*Ro2*phi1sq*phi1sq + alpha*phi1sq + beta*phi1 + x0*x0 + y0*y0 + z0*z0;
1495  double d2_0 = 0.25*Ro2*phi0sq*phi0sq + alpha*phi0sq + beta*phi0 + x0*x0 + y0*y0 + z0*z0;
1496 
1497  if (d2_0<d2_1 && d2_0<d2_2){
1498  phi=phi0;
1499  dist2=d2_0;
1500  }
1501  else if (d2_1<d2_0 && d2_1<d2_2){
1502  phi=phi1;
1503  dist2=d2_1;
1504  }
1505  else{
1506  phi=phi2;
1507  dist2=d2_2;
1508  }
1509  if (fabs(phi)<0.2){ // Check small angle approximation
1510  return numeric_limits<double>::quiet_NaN();
1511  }
1512 
1513  }
1514 
1515  // Calculate distance along track ("s")
1516  if(s!=NULL){
1517  double dz = dz_dphi*phi;
1518  double Rodphi = Ro*phi;
1519  double ds = sqrt(dz*dz + Rodphi*Rodphi);
1520  *s = step->s + (phi>0.0 ? ds:-ds);
1521  }
1522 
1523  this->last_phi = phi;
1524  this->last_swim_step = step;
1525  this->last_dz_dphi = dz_dphi;
1526 
1527  return sqrt(dist2);
1528 }
1529 
1530 //---------------------------------
1531 // FindClosestSwimStep
1532 //---------------------------------
1534 {
1535  /// Find the closest swim step to the given wire. The value of
1536  /// "L" should be the active wire length. The coordinate system
1537  /// defined by "wire" should have its origin at the center of
1538  /// the wire with the wire running in the direction of udir.
1539 
1540  if(istep_ptr)*istep_ptr=-1;
1541 
1542  if(Nswim_steps<1){
1543  _DBG_<<"No swim steps! You must \"Swim\" the track before calling FindClosestSwimStep(...)"<<endl;
1544  }
1545 
1546  // Make sure we have a wire first!
1547  if(!wire)return NULL;
1548 
1549  // Loop over swim steps and find the one closest to the wire
1550  swim_step_t *swim_step = swim_steps;
1551  swim_step_t *step=NULL;
1552  //double min_delta2 = 1.0E6;
1553  double old_delta2=1.0e6;
1554  double L_over_2 = wire->L/2.0; // half-length of wire in cm
1555  int istep=-1;
1556 
1557  double dx, dy, dz;
1558 
1559  // w is a vector to the origin of the wire
1560  // u is a unit vector along the wire
1561 
1562  double wx, wy, wz;
1563  double ux, uy, uz;
1564 
1565  wx = wire->origin.X();
1566  wy = wire->origin.Y();
1567  wz = wire->origin.Z();
1568 
1569  ux = wire->udir.X();
1570  uy = wire->udir.Y();
1571  uz = wire->udir.Z();
1572 
1573  int i;
1574  for(i=0; i<Nswim_steps; i++, swim_step++){
1575  // Find the point's position along the wire. If the point
1576  // is past the end of the wire, calculate the distance
1577  // from the end of the wire.
1578  // DVector3 pos_diff = swim_step->origin - wire->origin;
1579 
1580  dx = swim_step->origin.X() - wx;
1581  dy = swim_step->origin.Y() - wy;
1582  dz = swim_step->origin.Z() - wz;
1583 
1584  // double u = wire->udir.Dot(pos_diff);
1585  double u = ux * dx + uy * dy + uz * dz;
1586 
1587  // Find distance perpendicular to wire
1588  // double delta2 = pos_diff.Mag2() - u*u;
1589  double delta2 = dx*dx + dy*dy + dz*dz - u*u;
1590 
1591  // If point is past end of wire, calculate distance
1592  // from wire's end by adding on distance along wire direction.
1593  if( fabs(u)>L_over_2){
1594  // delta2 += pow(fabs(u)-L_over_2, 2.0);
1595  double u_minus_L_over_2=fabs(u)-L_over_2;
1596  delta2 += ( u_minus_L_over_2*u_minus_L_over_2 );
1597  // printf("step %d\n",i);
1598  }
1599 
1600  if(debug_level>3)_DBG_<<"delta2="<<delta2<<" old_delta2="<<old_delta2<<endl;
1601  if (delta2>old_delta2) break;
1602 
1603  //if(delta2 < min_delta2){
1604  // min_delta2 = delta2;
1605  step = swim_step;
1606  istep=i;
1607 
1608  //}
1609  //printf("%d delta %f min %f\n",i,delta2,min_delta2);
1610  old_delta2=delta2;
1611  }
1612 
1613  if(istep_ptr)*istep_ptr=istep;
1614 
1615  if(debug_level>3)_DBG_<<"found closest step at i="<<i<<" istep_ptr="<<istep_ptr<<endl;
1616 
1617  return step;
1618 }
1619 
1620 //---------------------------------
1621 // FindClosestSwimStep
1622 //---------------------------------
1624 {
1625  /// Find the closest swim step to the plane specified by origin
1626  /// and norm. origin should indicate any point in the plane and
1627  /// norm a vector normal to the plane.
1628  if(istep_ptr)*istep_ptr=-1;
1629 
1630  if(Nswim_steps<1){
1631  _DBG_<<"No swim steps! You must \"Swim\" the track before calling FindClosestSwimStep(...)"<<endl;
1632  }
1633 
1634  // Make sure normal vector is unit lenght
1635  norm.SetMag(1.0);
1636 
1637  // Loop over swim steps and find the one closest to the plane
1638  swim_step_t *swim_step = swim_steps;
1639  swim_step_t *step=NULL;
1640  //double min_dist = 1.0E6;
1641  double old_dist=1.0e6;
1642  int istep=-1;
1643 
1644  for(int i=0; i<Nswim_steps; i++, swim_step++){
1645 
1646  // Distance to plane is dot product of normal vector with any
1647  // vector pointing from the current step to a point in the plane
1648  double dist = fabs(norm.Dot(swim_step->origin-origin));
1649 
1650  if (dist>old_dist) break;
1651 
1652  // Check if we're the closest step
1653  //if(dist < min_dist){
1654  //min_dist = dist;
1655 
1656  step = swim_step;
1657  istep=i;
1658  //}
1659  old_dist=dist;
1660 
1661  // We should probably have a break condition here so we don't
1662  // waste time looking all the way to the end of the track after
1663  // we've passed the plane.
1664  }
1665 
1666  if(istep_ptr)*istep_ptr=istep;
1667 
1668  return step;
1669 }
1670 
1671 
1672 //---------------------------------
1673 // FindPlaneCrossing
1674 //---------------------------------
1676 {
1677  /// Find the closest swim step to the position where the track crosses
1678  /// the plane specified by origin
1679  /// and norm. origin should indicate any point in the plane and
1680  /// norm a vector normal to the plane.
1681 
1682  if(Nswim_steps<1){
1683  _DBG_<<"No swim steps! You must \"Swim\" the track before calling FindPlaneCrossing(...)"<<endl;
1684  raise(SIGSEGV);// force seg. fault
1685  }
1686 
1687  // Make sure normal vector is unit length
1688  norm.SetMag(1.0);
1689 
1690  // Loop over swim steps and find the one closest to the plane
1691  swim_step_t *swim_step = &swim_steps[first_i];
1692  swim_step_t *step=NULL;
1693  double old_dist=1.0e6;
1694 
1695  // Check if we should start from the beginning of the reference
1696  // trajectory or the end
1697  int last_index=Nswim_steps-1;
1698  double forward_dist= norm.Dot(swim_step->origin-origin);
1699  if( forward_dist == 0.0 ) return swim_step;
1700  double backward_dist= norm.Dot(swim_steps[last_index].origin-origin);
1701  if( backward_dist ==0.0 ) return &swim_steps[last_index];
1702  if (detector==SYS_START || fabs(forward_dist)<fabs(backward_dist)){ // start at beginning
1703  for(int i=first_i; i<Nswim_steps; i++, swim_step++){
1704 
1705  // Distance to plane is dot product of normal vector with any
1706  // vector pointing from the current step to a point in the plane
1707  //double dist = fabs(norm.Dot(swim_step->origin-origin));
1708  double dist = norm.Dot(swim_step->origin-origin);
1709 
1710  // We've crossed the plane when the sign of dist changes
1711  if (dist*old_dist<0 && i>0) {
1712  if (fabs(dist)<fabs(old_dist)){
1713  step=swim_step;
1714  }
1715  break;
1716  }
1717  step = swim_step;
1718  old_dist=dist;
1719  }
1720  }
1721  else{ // start at end
1722  for(int i=last_index; i>=0; i--){
1723  swim_step=&swim_steps[i];
1724  double dist = norm.Dot(swim_step->origin-origin);
1725  // We've crossed the plane when the sign of dist changes
1726  if (dist*old_dist<0 && i<last_index) {
1727  if (fabs(dist)<fabs(old_dist)){
1728  step=swim_step;
1729  }
1730  break;
1731  }
1732  step = swim_step;
1733  old_dist=dist;
1734  }
1735 
1736  }
1737 
1738  return step;
1739 }
1740 
1741 
1742 
1743 
1744 //---------------------------------
1745 // DistToRT
1746 //---------------------------------
1747 double DReferenceTrajectory::DistToRT(const DCoordinateSystem *wire, double *s) const
1748 {
1749  /// Find the closest distance to the given wire in cm. The value of
1750  /// "L" should be the active wire length (in cm). The coordinate system
1751  /// defined by "wire" should have its origin at the center of
1752  /// the wire with the wire running in the direction of udir.
1753  swim_step_t *step=FindClosestSwimStep(wire);
1754 
1755  return (step && step->s>0) ? DistToRT(wire, step, s):std::numeric_limits<double>::quiet_NaN();
1756 }
1757 
1758 //---------------------------------
1759 // DistToRTBruteForce
1760 //---------------------------------
1762 {
1763  /// Find the closest distance to the given wire in cm. The value of
1764  /// "L" should be the active wire length (in cm). The coordinate system
1765  /// defined by "wire" should have its origin at the center of
1766  /// the wire with the wire running in the direction of udir.
1767  swim_step_t *step=FindClosestSwimStep(wire);
1768 
1769  return step ? DistToRTBruteForce(wire, step, s):std::numeric_limits<double>::quiet_NaN();
1770 }
1771 
1772 //------------------
1773 // DistToRT
1774 //------------------
1775 double DReferenceTrajectory::DistToRT(const DCoordinateSystem *wire, const swim_step_t *step, double *s) const
1776 {
1777  /// Calculate the distance of the given wire(in the lab
1778  /// reference frame) to the Reference Trajectory which the
1779  /// given swim step belongs to. This uses the momentum directions
1780  /// and positions of the swim step
1781  /// to define a curve and calculate the distance of the hit
1782  /// from it. The swim step should be the closest one to the wire.
1783  /// IMPORTANT: This approximates the helix locally by a parabola.
1784  /// This means the swim step should be fairly close
1785  /// to the wire so that this approximation is valid. If the
1786  /// reference trajectory from which the swim step came is too
1787  /// sparse, the results will not be nearly as good.
1788 
1789  // Interestingly enough, this is one of the harder things to figure
1790  // out in the tracking code which is why the explanations may be
1791  // a bit long.
1792 
1793  // The general idea is to define the helix in a coordinate system
1794  // in which the wire runs along the z-axis. The distance to the
1795  // wire is then defined just in the X/Y plane of this coord. system.
1796  // The distance is expressed as a function of the phi angle in the
1797  // natural coordinate system of the helix. This way, phi=0 corresponds
1798  // to the swim step point itself and the DOCA point should be
1799  // at a small phi angle.
1800  //
1801  // The minimum distance between the helical segment and the wire
1802  // will be a function of sin(phi), cos(phi) and phi. Approximating
1803  // sin(phi) by phi and cos(phi) by (1-phi^2) leaves a 4th order
1804  // polynomial in phi. Taking the derivative leaves a 3rd order
1805  // polynomial whose root is the phi corresponding to the
1806  // Distance Of Closest Approach(DOCA) point on the helix. Plugging
1807  // that value of phi back into the distance formula gives
1808  // us the minimum distance between the track and the wire.
1809 
1810  // First, we need to define the coordinate system in which the
1811  // wire runs along the z-axis. This is actually done already
1812  // in the CDC package for each wire once, at program start.
1813  // The directions of the axes are defined in wire->sdir,
1814  // wire->tdir, and wire->udir.
1815 
1816  // Next, define a point on the helical segment defined by the
1817  // swim step it the RT coordinate system. The directions of
1818  // the RT coordinate system are defined by step->xdir, step->ydir,
1819  // and step->zdir. The coordinates of a point on the helix
1820  // in this coordinate system are:
1821  //
1822  // x = Ro*(cos(phi) - 1)
1823  // y = Ro*sin(phi)
1824  // z = phi*(dz/dphi)
1825  //
1826  // where phi is the phi angle of the point in this coordinate system.
1827 
1828  // Now, a vector describing the helical point in the LAB coordinate
1829  // system is:
1830  //
1831  // h = x*xdir + y*ydir + z*zdir + pos
1832  //
1833  // where h,xdir,ydir,zdir and pos are all 3-vectors.
1834  // xdir,ydir,zdir are unit vectors defining the directions
1835  // of the RT coord. system axes in the lab coord. system.
1836  // pos is a vector defining the position of the swim step
1837  // in the lab coord.system
1838 
1839  // Now we just need to find the extent of "h" in the wire's
1840  // coordinate system (period . means dot product):
1841  //
1842  // s = (h-wpos).sdir
1843  // t = (h-wpos).tdir
1844  // u = (h-wpos).udir
1845  //
1846  // where wpos is the position of the center of the wire in
1847  // the lab coord. system and is given by wire->wpos.
1848 
1849  // At this point, the values of s,t, and u repesent a point
1850  // on the helix in the coord. system of the wire with the
1851  // wire in the "u" direction and positioned at the origin.
1852  // The distance(squared) from the wire to the point on the helix
1853  // is given by:
1854  //
1855  // d^2 = s^2 + t^2
1856  //
1857  // where s and t are both functions of phi.
1858 
1859  // So, we'll define the values of "s" and "t" above as:
1860  //
1861  // s = A*x + B*y + C*z + D
1862  // t = E*x + F*y + G*z + H
1863  //
1864  // where A,B,C,D,E,F,G, and H are constants defined below
1865  // and x,y,z are all functions of phi defined above.
1866  // (period . means dot product)
1867  //
1868  // A = sdir.xdir
1869  // B = sdir.ydir
1870  // C = sdir.zdir
1871  // D = sdir.(pos-wpos)
1872  //
1873  // E = tdir.xdir
1874  // F = tdir.ydir
1875  // G = tdir.zdir
1876  // H = tdir.(pos-wpos)
1877  const DVector3 &xdir = step->sdir;
1878  const DVector3 &ydir = step->tdir;
1879  const DVector3 &zdir = step->udir;
1880  const DVector3 &sdir = wire->sdir;
1881  const DVector3 &tdir = wire->tdir;
1882  const DVector3 &udir = wire->udir;
1883  DVector3 pos_diff = step->origin - wire->origin;
1884 
1885  double A = sdir.Dot(xdir);
1886  double B = sdir.Dot(ydir);
1887  double C = sdir.Dot(zdir);
1888  double D = sdir.Dot(pos_diff);
1889 
1890  double E = tdir.Dot(xdir);
1891  double F = tdir.Dot(ydir);
1892  double G = tdir.Dot(zdir);
1893  double H = tdir.Dot(pos_diff);
1894 
1895  // OK, here is the dirty part. Using the approximations given above
1896  // to write the x and y functions in terms of phi^2 and phi (instead
1897  // of cos and sin) we put them into the equations for s and t above.
1898  // Then, inserting those into the equation for d^2 above that, we
1899  // get a very long equation in terms of the constants A,...H and
1900  // phi up to 4th order. Combining coefficients for similar powers
1901  // of phi yields an equation of the form:
1902  //
1903  // d^2 = Q*phi^4 + R*phi^3 + S*phi^2 + T*phi + U
1904  //
1905  // The dirty part is that it takes the better part of a sheet of
1906  // paper to work out the relations for Q,...U in terms of
1907  // A,...H, and Ro, dz/dphi. You can work it out yourself on
1908  // paper to verify that the equations below are correct.
1909  double Ro = step->Ro;
1910  double Ro2 = Ro*Ro;
1911  double delta_z = step->mom.Dot(step->udir);
1912  double delta_phi = step->mom.Dot(step->tdir)/Ro;
1913  double dz_dphi = delta_z/delta_phi;
1914  double dz_dphi2=dz_dphi*dz_dphi;
1915  double Ro_dz_dphi=Ro*dz_dphi;
1916 
1917  // double Q = pow(A*Ro/2.0, 2.0) + pow(E*Ro/2.0, 2.0);
1918  double Q=0.25*Ro2*(A*A+E*E);
1919  // double R = -(2.0*A*B*Ro2 + 2.0*A*C*Ro_dz_dphi + 2.0*E*F*Ro2 + 2.0*E*G*Ro_dz_dphi)/2.0;
1920  double R = -((A*B+E*F)*Ro2 + (A*C+E*G)*Ro_dz_dphi);
1921  // double S = pow(B*Ro, 2.0) + pow(C*dz_dphi,2.0) + 2.0*B*C*Ro_dz_dphi - 2.0*A*D*Ro/2.0
1922  //+ pow(F*Ro, 2.0) + pow(G*dz_dphi,2.0) + 2.0*F*G*Ro_dz_dphi - 2.0*E*H*Ro/2.0;
1923  double S= (B*B+F*F)*Ro2+(C*C+G*G)*dz_dphi2+2.0*(B*C+F*G)*Ro_dz_dphi
1924  -(A*D+E*H)*Ro;
1925  // double T = 2.0*B*D*Ro + 2.0*C*D*dz_dphi + 2.0*F*H*Ro + 2.0*G*H*dz_dphi;
1926  double T = 2.0*((B*D+F*H)*Ro + (C*D+G*H)*dz_dphi);
1927  double U = D*D + H*H;
1928 
1929  // Aaarghh! my fingers hurt just from typing all of that!
1930  //
1931  // OK, now we differentiate the above equation for d^2 to get:
1932  //
1933  // d(d^2)/dphi = 4*Q*phi^3 + 3*R*phi^2 + 2*S*phi + T
1934  //
1935  // NOTE: don't confuse "R" with "Ro" in the above equations!
1936  //
1937  // Now we have to solve the 3rd order polynomial for the phi value of
1938  // the point of closest approach on the RT. This is a well documented
1939  // procedure. Essentially, when you have an equation of the form:
1940  //
1941  // x^3 + a2*x^2 + a1*x + a0 = 0;
1942  //
1943  // a change of variables is made such that w = x + a2/3 which leads
1944  // to a third order poly with no w^2 term:
1945  //
1946  // w^3 + 3.0*b*w + 2*c = 0
1947  //
1948  // where:
1949  // b = a1/3 - (a2^2)/9
1950  // c = a0/2 - a1*a2/6 + (a2^3)/27
1951  //
1952  // The one real root of this is:
1953  //
1954  // w0 = q - p
1955  //
1956  // where:
1957  // q^3 = d - c
1958  // p^3 = d + c
1959  // d^2 = b^3 + c^2 (don't confuse with d^2 above!)
1960  //
1961  // For us this means that:
1962  // a2 = 3*R/(4*Q)
1963  // a1 = 2*S/(4*Q)
1964  // a0 = T/(4*Q)
1965  //
1966  // A potential problem could occur if Q is at or very close to zero.
1967  // This situation occurs when both A and E are zero. This would mean
1968  // that both sdir and tdir are perpendicular to xdir which means
1969  // xdir is in the same direction as udir (got that?). Physically,
1970  // this corresponds to the situation when both the momentum and
1971  // the magnetic field are perpendicular to the wire (though not
1972  // necessarily perpendicular to each other). This situation can't
1973  // really occur in the CDC detector where the chambers are well
1974  // contained in a region where the field is essentially along z as
1975  // are the wires.
1976  //
1977  // Just to be safe, we check that Q is greater than
1978  // some minimum before solving for phi. If it is too small, we fall
1979  // back to solving the quadratic equation for phi.
1980  double phi =0.0;
1981  if(fabs(Q)>1.0E-6){
1982  /*
1983  double fourQ = 4.0*Q;
1984  double a2 = 3.0*R/fourQ;
1985  double a1 = 2.0*S/fourQ;
1986  double a0 = T/fourQ;
1987  */
1988  double one_over_fourQ=0.25/Q;
1989  double a2=3.0*R*one_over_fourQ;
1990  double a1=2.0*S*one_over_fourQ;
1991  double a0=T*one_over_fourQ;
1992  double a2sq=a2*a2;
1993  /*
1994  double b = a1/3.0 - a2*a2/9.0;
1995  double c = a0/2.0 - a1*a2/6.0 + a2*a2*a2/27.0;
1996  */
1997  double b=ONE_THIRD*(a1-ONE_THIRD*a2sq);
1998  double c=0.5*(a0-ONE_THIRD*a1*a2)+a2*a2sq/27.0;
1999  double my_d2=b*b*b+c*c;
2000  if (my_d2>0){
2001  //double d = sqrt(pow(b, 3.0) + pow(c, 2.0)); // occasionally, this is zero. See below
2002  double d=sqrt(my_d2);
2003  //double q = pow(d - c, ONE_THIRD);
2004  //double p = pow(d + c, ONE_THIRD);
2005  double q=cbrt(d-c);
2006  double p=cbrt(d+c);
2007 
2008  double w0 = q - p;
2009  //phi = w0 - a2/3.0;
2010  phi = w0 - ONE_THIRD*a2;
2011  }
2012  else{
2013  // Use DeMoivre's theorem to find the cube root of a complex
2014  // number. In this case there are three real solutions.
2015  double d=sqrt(-my_d2);
2016  c*=-1.;
2017  double temp=sqrt(cbrt(c*c+d*d));
2018  double theta1=ONE_THIRD*atan2(d,c);
2019  double sum_over_2=temp*cos(theta1);
2020  double diff_over_2=-temp*sin(theta1);
2021 
2022  double phi0=-a2/3+2.*sum_over_2;
2023  double phi1=-a2/3-sum_over_2+sqrt(3)*diff_over_2;
2024  double phi2=-a2/3-sum_over_2-sqrt(3)*diff_over_2;
2025 
2026  double d2_0 = U + phi0*(T + phi0*(S + phi0*(R + phi0*Q)));
2027  double d2_1 = U + phi1*(T + phi1*(S + phi1*(R + phi1*Q)));
2028  double d2_2 = U + phi2*(T + phi2*(S + phi2*(R + phi2*Q)));
2029 
2030  if (d2_0<d2_1 && d2_0<d2_2){
2031  phi=phi0;
2032  }
2033  else if (d2_1<d2_0 && d2_1<d2_2){
2034  phi=phi1;
2035  }
2036  else{
2037  phi=phi2;
2038  }
2039  }
2040  }
2041 
2042  if(fabs(Q)<=1.0E-6 || !isfinite(phi)){
2043  double a = 3.0*R;
2044  double b = 2.0*S;
2045  double c = 1.0*T;
2046  phi = (-b + sqrt(b*b - 4.0*a*c))/(2.0*a);
2047  }
2048 
2049  // The accuracy of this method is limited by how close the step is to the
2050  // actual minimum. If the value of phi is large then the step size is
2051  // not too close and we should add another couple of steps in the right
2052  // place in order to get a more accurate value. Note that while this will
2053  // increase the time it takes this round, presumably the fitter will be
2054  // calling this often for each wire and having a high density of points
2055  // near the wires will just make subsequent calls go quicker. This also
2056  // allows larger initial step sizes with the high density regions getting
2057  // filled in as needed leading to overall faster tracking.
2058 #if 0
2059  if(isfinite(phi) && fabs(phi)>2.0E-4){
2060  if(dist_to_rt_depth>=3){
2061  _DBG_<<"3 or more recursive calls to DistToRT(). Something is wrong! bailing ..."<<endl;
2062  //for(int k=0; k<Nswim_steps; k++){
2063  // DVector3 &v = swim_steps[k].origin;
2064  // _DBG_<<" "<<k<<": "<<v.X()<<", "<<v.Y()<<", "<<v.Z()<<endl;
2065  //}
2066  //exit(-1);
2067  return std::numeric_limits<double>::quiet_NaN();
2068  }
2069  double scale_step = 1.0;
2070  double s_range = 1.0*scale_step;
2071  double step_size = 0.02*scale_step;
2072  int err = InsertSteps(step, phi>0.0 ? +s_range:-s_range, step_size); // Add new steps near this step by swimming in the direction of phi
2073  if(!err){
2074  step=FindClosestSwimStep(wire); // Find the new closest step
2075  if(!step)return std::numeric_limits<double>::quiet_NaN();
2076  dist_to_rt_depth++;
2077  double doca = DistToRT(wire, step, s); // re-call ourself with the new step
2078  dist_to_rt_depth--;
2079  return doca;
2080  }else{
2081  if(err<0)return std::numeric_limits<double>::quiet_NaN();
2082 
2083  // If InsertSteps() returns an error > 0 then it indicates that it
2084  // was unable to add additional steps (perhaps because there
2085  // aren't enough spaces available). In that case, we just go ahead
2086  // and use the phi we have and make the best estimate possible.
2087  }
2088  }
2089 #endif
2090 
2091  // It is possible at this point that the value of phi corresponds to
2092  // a point past the end of the wire. We should check for this here and
2093  // recalculate, if necessary, the DOCA at the end of the wire. First,
2094  // calculate h (the vector defined way up above) and dot it into the
2095  // wire's u-direction to get the position of the DOCA point along the
2096  // wire.
2097  double x = -0.5*Ro*phi*phi;
2098  double y = Ro*phi;
2099  double z = dz_dphi*phi;
2100  DVector3 h = pos_diff + x*xdir + y*ydir + z*zdir;
2101  double u = h.Dot(udir);
2102  if(fabs(u) > wire->L/2.0){
2103  // Looks like our DOCA point is past the end of the wire.
2104  // Find phi corresponding to the end of the wire.
2105  double L_over_2 = u>0.0 ? wire->L/2.0:-wire->L/2.0;
2106  double a = -0.5*Ro*udir.Dot(xdir);
2107  double b = Ro*udir.Dot(ydir) + dz_dphi*udir.Dot(zdir);
2108  double c = udir.Dot(pos_diff) - L_over_2;
2109  double twoa=2.0*a;
2110  double sqroot=sqrt(b*b-4.0*a*c);
2111  double phi1 = (-b + sqroot)/(twoa);
2112  double phi2 = (-b - sqroot)/(twoa);
2113  phi = fabs(phi1)<fabs(phi2) ? phi1:phi2;
2114  u=L_over_2;
2115  }
2116  this->last_dist_along_wire = u;
2117 
2118  // Use phi to calculate DOCA
2119  double d2 = U + phi*(T + phi*(S + phi*(R + phi*Q)));
2120  double d = sqrt(d2);
2121 
2122  // Calculate distance along track ("s")
2123  double dz = dz_dphi*phi;
2124  double Rodphi = Ro*phi;
2125  double ds = sqrt(dz*dz + Rodphi*Rodphi);
2126  if(s)*s=step->s + (phi>0.0 ? ds:-ds);
2127  if(debug_level>3){
2128  _DBG_<<"distance to rt: "<< step->s + (phi>0.0 ? ds:-ds) <<" from step at "<<step->s<<" with ds="<<ds<<" d="<<d<<" dz="<<dz<<" Rodphi="<<Rodphi<<endl;
2129  _DBG_<<"phi="<<phi<<" U="<<U<<" u="<<u<<endl;
2130  }
2131 
2132  // Remember phi and step so additional info on the point can be obtained
2133  this->last_phi = phi;
2134  this->last_swim_step = step;
2135  this->last_dz_dphi = dz_dphi;
2136 
2137  return d; // WARNING: This could return nan!
2138 }
2139 
2140 //------------------
2141 // DistToRTBruteForce
2142 //------------------
2143 double DReferenceTrajectory::DistToRTBruteForce(const DCoordinateSystem *wire, const swim_step_t *step, double *s) const
2144 {
2145  /// Calculate the distance of the given wire(in the lab
2146  /// reference frame) to the Reference Trajectory which the
2147  /// given swim step belongs to. This uses the momentum directions
2148  /// and positions of the swim step
2149  /// to define a curve and calculate the distance of the hit
2150  /// from it. The swim step should be the closest one to the wire.
2151  /// IMPORTANT: This calculates the distance using a "brute force"
2152  /// method of taking tiny swim steps to find the minimum distance.
2153  /// It is vey SLOW and you should be using DistToRT(...) instead.
2154  /// This is only here to provide an independent check of DistToRT(...).
2155 
2156  const DVector3 &xdir = step->sdir;
2157  const DVector3 &ydir = step->tdir;
2158  const DVector3 &zdir = step->udir;
2159  const DVector3 &sdir = wire->sdir;
2160  const DVector3 &tdir = wire->tdir;
2161  DVector3 pos_diff = step->origin - wire->origin;
2162 
2163  double Ro = step->Ro;
2164  double delta_z = step->mom.Dot(step->udir);
2165  double delta_phi = step->mom.Dot(step->tdir)/Ro;
2166  double dz_dphi = delta_z/delta_phi;
2167 
2168  // Brute force
2169  double min_d2 = 1.0E6;
2170  double phi=M_PI;
2171  for(int i=-2000; i<2000; i++){
2172  double myphi=(double)i*0.000005;
2173  DVector3 d = Ro*(cos(myphi)-1.0)*xdir
2174  + Ro*sin(myphi)*ydir
2175  + dz_dphi*myphi*zdir
2176  + pos_diff;
2177 
2178  double d2 = pow(d.Dot(sdir),2.0) + pow(d.Dot(tdir),2.0);
2179  if(d2<min_d2){
2180  min_d2 = d2;
2181  phi = myphi;
2182  this->last_phi = myphi;
2183  }
2184  }
2185  double d2 = min_d2;
2186  double d = sqrt(d2);
2187  this->last_phi = phi;
2188  this->last_swim_step = step;
2189  this->last_dz_dphi = dz_dphi;
2190 
2191  // Calculate distance along track ("s")
2192  double dz = dz_dphi*phi;
2193  double Rodphi = Ro*phi;
2194  double ds = sqrt(dz*dz + Rodphi*Rodphi);
2195  if(s)*s=step->s + (phi>0.0 ? ds:-ds);
2196 
2197  return d;
2198 }
2199 
2200 //------------------
2201 // Straw_dx
2202 //------------------
2203 double DReferenceTrajectory::Straw_dx(const DCoordinateSystem *wire, double radius) const
2204 {
2205  /// Find the distance traveled within the specified radius of the
2206  /// specified wire. This will give the "dx" component of a dE/dx
2207  /// measurement for cylindrical geometry as we have with straw tubes.
2208  ///
2209  /// At this point, the estimate is done using a simple linear
2210  /// extrapolation from the DOCA point in the direction of the momentum
2211  /// to the 2 points at which it itersects the given radius. Segments
2212  /// which extend past the end of the wire will be clipped to the end
2213  /// of the wire before calculating the total dx.
2214 
2215  // First, find the DOCA point for this wire
2216  double s;
2217  double doca = DistToRT(wire, &s);
2218  if(!isfinite(doca))
2219  return 0.0;
2220 
2221  // If doca is outside of the given radius, then we're done
2222  if(doca>=radius)return 0.0;
2223 
2224  // Get the location and momentum direction of the DOCA point
2225  DVector3 pos, momdir;
2226  GetLastDOCAPoint(pos, momdir);
2227  if(momdir.Mag()!=0.0)momdir.SetMag(1.0);
2228 
2229  // Get wire direction
2230  const DVector3 &udir = wire->udir;
2231 
2232  // Calculate vectors used to form quadratic equation for "alpha"
2233  // the distance along the mometum direction from the DOCA point
2234  // to the intersection with a cylinder of the given radius.
2235  DVector3 A = udir.Cross(pos-wire->origin);
2236  DVector3 B = udir.Cross(momdir);
2237 
2238  // If the magnitude of B is zero at this point, it means the momentum
2239  // direction is parallel to the wire. In this case, this method will
2240  // not work. Return NaN.
2241  if(B.Mag()<1.0E-10)return std::numeric_limits<double>::quiet_NaN();
2242 
2243  double a = B.Mag();
2244  double b = A.Dot(B);
2245  double c = A.Mag() - radius;
2246  double d = sqrt(b*b - 4.0*a*c);
2247 
2248  // The 2 roots should correspond to the 2 intersection points.
2249  double alpha1 = (-b + d)/(2.0*a);
2250  double alpha2 = (-b - d)/(2.0*a);
2251 
2252  DVector3 int1 = pos + alpha1*momdir;
2253  DVector3 int2 = pos + alpha2*momdir;
2254 
2255  // Check if point1 is past the end of the wire
2256  double q = udir.Dot(int1 - wire->origin);
2257  if(fabs(q) > wire->L/2.0){
2258  double gamma = udir.Dot(wire->origin - pos) + (q>0.0 ? +1.0:-1.0)*wire->L/2.0;
2259  gamma /= momdir.Dot(udir);
2260  int1 = pos + gamma*momdir;
2261  }
2262 
2263  // Check if point2 is past the end of the wire
2264  q = udir.Dot(int2 - wire->origin);
2265  if(fabs(q) > wire->L/2.0){
2266  double gamma = udir.Dot(wire->origin - pos) + (q>0.0 ? +1.0:-1.0)*wire->L/2.0;
2267  gamma /= momdir.Dot(udir);
2268  int2 = pos + gamma*momdir;
2269  }
2270 
2271  // Calculate distance
2272  DVector3 delta = int1 - int2;
2273 
2274  return delta.Mag();
2275 }
2276 
2277 //------------------
2278 // GetLastDOCAPoint
2279 //------------------
2281 {
2282  /// Use values saved by the last call to one of the DistToRT functions
2283  /// to calculate the 3-D DOCA position in lab coordinates and momentum
2284  /// in GeV/c.
2285 
2286  if(last_swim_step==NULL){
2287  if(Nswim_steps>0){
2288  last_swim_step = &swim_steps[0];
2289  last_phi = 0.0;
2290  }else{
2291  pos.SetXYZ(QuietNaN,QuietNaN,QuietNaN);
2292  mom.SetXYZ(QuietNaN,QuietNaN,QuietNaN);
2293  return;
2294  }
2295  }
2296 
2297  // If last_phi is not finite, set it to 0 as a last resort
2298  if(!isfinite(last_phi))last_phi = 0.0;
2299 
2300  const DVector3 &xdir = last_swim_step->sdir;
2301  const DVector3 &ydir = last_swim_step->tdir;
2302  const DVector3 &zdir = last_swim_step->udir;
2303 
2304  double x = -(last_swim_step->Ro/2.0)*last_phi*last_phi;
2305  double y = last_swim_step->Ro*last_phi;
2306  double z = last_dz_dphi*last_phi;
2307 
2308  pos = last_swim_step->origin + x*xdir + y*ydir + z*zdir;
2309  mom = last_swim_step->mom;
2310 
2311  mom.Rotate(-last_phi, zdir);
2312 }
2313 
2314 //------------------
2315 // GetLastDOCAPoint
2316 //------------------
2318 {
2319  /// Use values saved by the last call to one of the DistToRT functions
2320  /// to calculate the 3-D DOCA position in lab coordinates. This is
2321  /// mainly intended for debugging.
2322  if(last_swim_step==NULL){
2323  if(Nswim_steps>0){
2324  last_swim_step = &swim_steps[0];
2325  last_phi = 0.0;
2326  }else{
2328  }
2329  }
2330  const DVector3 &xdir = last_swim_step->sdir;
2331  const DVector3 &ydir = last_swim_step->tdir;
2332  const DVector3 &zdir = last_swim_step->udir;
2333  double Ro = last_swim_step->Ro;
2334  double delta_z = last_swim_step->mom.Dot(zdir);
2335  double delta_phi = last_swim_step->mom.Dot(ydir)/Ro;
2336  double dz_dphi = delta_z/delta_phi;
2337 
2338  double x = -(Ro/2.0)*last_phi*last_phi;
2339  double y = Ro*last_phi;
2340  double z = dz_dphi*last_phi;
2341 
2342  return last_swim_step->origin + x*xdir + y*ydir + z*zdir;
2343 }
2344 
2345 //------------------
2346 // dPdx
2347 //------------------
2348 double DReferenceTrajectory::dPdx_from_A_Z_rho(double ptot, double A, double Z, double density) const
2349 {
2350  double I = (Z*12.0 + 7.0)*1.0E-9; // From Leo 2nd ed. pg 25.
2351  if (Z>=13) I=(9.76*Z+58.8*pow(Z,-0.19))*1.0e-9;
2352  double rhoZ_overA=density*Z/A;
2353  double KrhoZ_overA = 0.1535e-3*rhoZ_overA;
2354 
2355  return dPdx(ptot, KrhoZ_overA,rhoZ_overA,log(I));
2356 }
2357 
2358 //------------------
2359 // dPdx
2360 //------------------
2361 double DReferenceTrajectory::dPdx(double ptot, double KrhoZ_overA,
2362  double rhoZ_overA,double LogI) const
2363 {
2364  /// Calculate the momentum loss per unit distance traversed of the material with
2365  /// the given A, Z, and density. Value returned is in GeV/c per cm
2366  /// This follows the July 2008 PDG section 27.2 ppg 268-270.
2367  if(mass==0.0)return 0.0; // no ionization losses for neutrals
2368 
2369  double gammabeta = ptot/mass;
2370  double gammabeta2=gammabeta*gammabeta;
2371  double gamma = sqrt(gammabeta2+1.);
2372  double beta = gammabeta/gamma;
2373  double beta2=beta*beta;
2374  double me = 0.511E-3;
2375  double m_ratio=me/mass;
2376  double two_me_gammabeta2=2.*me*gammabeta2;
2377 
2378  double Tmax = two_me_gammabeta2/(1.0+2.0*gamma*m_ratio+m_ratio*m_ratio);
2379  //double K = 0.307075E-3; // GeV gm^-1 cm^2
2380  // Density effect
2381  double delta=0.;
2382  double X=log10(gammabeta);
2383  double X0,X1;
2384  double Cbar=2.*(LogI-log(28.816e-9*sqrt(rhoZ_overA)))+1.;
2385  if (rhoZ_overA>0.01){ // not a gas
2386  if (LogI<-1.6118){ // I<100
2387  if (Cbar<=3.681) X0=0.2;
2388  else X0=0.326*Cbar-1.;
2389  X1=2.;
2390  }
2391  else{
2392  if (Cbar<=5.215) X0=0.2;
2393  else X0=0.326*Cbar-1.5;
2394  X1=3.;
2395  }
2396  }
2397  else { // gases
2398  X1=4.;
2399  if (Cbar<=9.5) X0=1.6;
2400  else if (Cbar>9.5 && Cbar<=10.) X0=1.7;
2401  else if (Cbar>10 && Cbar<=10.5) X0=1.8;
2402  else if (Cbar>10.5 && Cbar<=11.) X0=1.9;
2403  else if (Cbar>11.0 && Cbar<=12.25) X0=2.;
2404  else if (Cbar>12.25 && Cbar<=13.804){
2405  X0=2.;
2406  X1=5.;
2407  }
2408  else {
2409  X0=0.326*Cbar-2.5;
2410  X1=5.;
2411  }
2412  }
2413  if (X>=X0 && X<X1)
2414  delta=4.606*X-Cbar+(Cbar-4.606*X0)*pow((X1-X)/(X1-X0),3.);
2415  else if (X>=X1)
2416  delta= 4.606*X-Cbar;
2417 
2418  double dEdx = KrhoZ_overA/beta2*(log(two_me_gammabeta2*Tmax)
2419  -2.*LogI - 2.0*beta2 -delta);
2420 
2421  double dP_dx = dEdx/beta;
2422 
2423  //double g = 0.350/sqrt(-log(0.06));
2424  //dP_dx *= 1.0 + exp(-pow(ptot/g,2.0)); // empirical for really low momentum particles
2425 
2426 
2427  if(ploss_direction==kBackward)dP_dx = -dP_dx;
2428 
2429  return dP_dx;
2430 }
2431 
2432 //------------------
2433 // Dump
2434 //------------------
2436 {
2437  swim_step_t *step = swim_steps;
2438  for(int i=0; i<Nswim_steps; i++, step++){
2439  vector<pair<string,string> > item;
2440  double x = step->origin.X();
2441  double y = step->origin.Y();
2442  double z = step->origin.Z();
2443  if(z<zmin || z>zmax)continue;
2444 
2445  double px = step->mom.X();
2446  double py = step->mom.Y();
2447  double pz = step->mom.Z();
2448 
2449  cout<<i<<": ";
2450  cout<<"(x,y,z)=("<<x<<","<<y<<","<<z<<") ";
2451  cout<<"(px,py,pz)=("<<px<<","<<py<<","<<pz<<") ";
2452  cout<<"(Ro,s,t)=("<<step->Ro<<","<<step->s<<","<<step->t<<") ";
2453  cout<<endl;
2454  }
2455 
2456 }
2457 
2458 // Propagate the covariance matrix for {px,py,pz,x,y,z,t} along the step ds
2459 jerror_t DReferenceTrajectory::PropagateCovariance(double ds,double q,
2460  double mass_sq,
2461  const DVector3 &mom,
2462  const DVector3 &pos,
2463  const DVector3 &B,
2464  TMatrixFSym &C) const{
2465  DMatrix J(7,7);
2466 
2467  double one_over_p_sq=1./mom.Mag2();
2468  double one_over_p=sqrt(one_over_p_sq);
2469  double px=mom.X();
2470  double py=mom.Y();
2471  double pz=mom.Z();
2472  double Bx=B.x(),By=B.y(),Bz=B.z();
2473 
2474  double ds_over_p=ds*one_over_p;
2475  double factor=0.003*q*ds_over_p;
2476  double temp=(Bz*py-Bx*pz)*one_over_p_sq;
2477  J(0,0)=1-factor*px*temp;
2478  J(0,1)=factor*(Bz-py*temp);
2479  J(0,2)=-factor*(By+pz*temp);
2480 
2481  temp=(Bx*pz-Bz*px)*one_over_p_sq;
2482  J(1,0)=-factor*(Bz+px*temp);
2483  J(1,1)=1-factor*py*temp;
2484  J(1,2)=factor*(Bx-pz*temp);
2485 
2486  temp=(By*px-Bx*py)*one_over_p_sq;
2487  J(2,0)=factor*(By-px*temp);
2488  J(2,1)=-factor*(Bx+py*temp);
2489  J(2,2)=1-factor*pz*temp;
2490 
2491  J(3,3)=1.;
2492  double ds_over_p3=one_over_p_sq*ds_over_p;
2493  J(3,0)=ds_over_p*(1-px*px*one_over_p_sq);
2494  J(3,1)=-px*py*ds_over_p3;
2495  J(3,2)=-px*pz*ds_over_p3;
2496 
2497  J(4,4)=1.;
2498  J(4,0)=J(3,1);
2499  J(4,1)=ds_over_p*(1-py*py*one_over_p_sq);
2500  J(4,2)=-py*pz*ds_over_p3;
2501 
2502  J(5,5)=1.;
2503  J(5,0)=J(3,2);
2504  J(5,1)=J(4,2);
2505  J(5,2)=ds_over_p*(1-pz*pz*one_over_p_sq);
2506 
2507  J(6,6)=1.;
2508 
2509  double fac2=(-ds/SPEED_OF_LIGHT)*mass_sq*one_over_p_sq*one_over_p_sq
2510  /sqrt(1.+mass_sq*one_over_p_sq);
2511  J(6,0)=fac2*px;
2512  J(6,1)=fac2*py;
2513  J(6,2)=fac2*pz;
2514 
2515  C=C.Similarity(J);
2516 
2517  return NOERROR;
2518 }
2519 
2520 // Find the position along a reference trajectory closest to a line.
2521 // The error matrix for the line can also be input via a pointer. The error
2522 // matrix is expected to be 7x7 with the order {Px,Py,Pz,X,Y,Z,T}.
2523 // Outputs the kinematic data object (including the covariance) at this
2524 // position, and the doca and the variance on the doca.
2526  const DVector3 &dir,
2527  const DMatrixDSym *covline,
2528  DKinematicData *track_kd,
2529  DVector3 &commonpos, double &doca, double &var_doca) const{
2530  const swim_step_t *swim_step=this->swim_steps;
2531 
2532  shared_ptr<TMatrixFSym> cov = (track_kd!=NULL) ? dResourcePool_TMatrixFSym->Get_SharedResource() : nullptr;
2533  if(track_kd!=NULL)
2534  {
2535  cov->ResizeTo(7, 7);
2536  *cov = *(track_kd->errorMatrix());
2537  }
2538  doca=1000.;
2539  double tflight=0.;
2540  double mass_sq=this->mass_sq;
2541  double q=this->q;
2542  double step_size=1.0,s=-step_size;
2543  DVector3 oldpos,oldmom;
2544  DVector3 point=origin;
2545 
2546  // Find the magnitude of the direction vector
2547  double pscale=dir.Mag();
2548  // If the magnitude of the direction vector is zero, don't bother to propagate
2549  // along a line from the input origin...
2550  bool move_along_line=(pscale>0)?true:false;
2551 
2552  // Propagate along the reference trajectory, comparing to the line at each
2553  // step
2554  for (int i=0;i<this->Nswim_steps-1; i++, swim_step++){
2555  DVector3 pos=swim_step->origin;
2556  DVector3 diff=pos-point;
2557  double new_doca=diff.Mag();
2558  if (new_doca>doca){
2559  if (i==1){ // backtrack to find the true doca
2560  tflight=0.;
2561 
2562  swim_step=this->swim_steps;
2563  if(track_kd!=NULL)
2564  *cov=*track_kd->errorMatrix();
2565 
2566  pos=swim_step->origin;
2567  DVector3 mom=swim_step->mom;
2568  DMagneticFieldStepper stepper(this->bfield, this->q, &pos, &mom);
2569 
2570  int inew=0;
2571  while (inew<100){
2572  DVector3 B;
2573  double ds=stepper.Step(&pos,&B,-0.5);
2574  // Compute the revised estimate for the doca
2575  diff=pos-point;
2576  new_doca=diff.Mag();
2577 
2578  if(new_doca > doca) break;
2579 
2580  // Propagate the covariance matrix of the track along the trajectory
2581  if(track_kd!=NULL){
2582  this->PropagateCovariance(ds,q,mass_sq,mom,oldpos,B,*cov);
2583  }
2584 
2585  // Store the current positions, doca and adjust flight times
2586  oldpos=pos;
2587  doca=new_doca;
2588 
2589  double one_over_p_sq=1./mom.Mag2();
2590  tflight+=ds*sqrt(1.+mass_sq*one_over_p_sq)/SPEED_OF_LIGHT;
2591 
2592  // New momentum
2593  stepper.GetMomentum(mom);
2594 
2595  oldmom=/*(-1.)*/mom;
2596  inew++;
2597 
2598  // New point on line
2599  if (move_along_line){
2600  point-=(step_size/pscale)*dir;
2601  s-=step_size;
2602  }
2603  }
2604  }
2605  if(track_kd!=NULL)
2606  {
2607  track_kd->setErrorMatrix(cov);
2608  track_kd->setMomentum(oldmom);
2609  track_kd->setPosition(oldpos);
2610  track_kd->setTime(track_kd->time() + tflight);
2611  }
2612 
2613  // Compute the variance on the doca
2614  diff=oldpos-point;
2615  double dx=diff.x();
2616  double dy=diff.y();
2617  double dz=diff.z();
2618 
2619  if(track_kd==NULL)
2620  break;
2621  //calculate var_doca
2622  if (covline==NULL){
2623  var_doca=(dx*dx*((*cov)(kX,kX))+dy*dy*((*cov)(kY,kY))
2624  +dz*dz*((*cov)(kZ,kZ))+2.*dx*dy*((*cov)(kX,kY))
2625  +2.*dx*dz*((*cov)(kX,kZ))+2.*dy*dz*((*cov)(kY,kZ)))
2626  /(doca*doca);
2627  }
2628  else{
2629  DMatrixDSym cov2(*covline);
2630  if (move_along_line){
2631  double two_s=2.*s;
2632  double s_sq=s*s;
2633  cov2(kX,kX)+=two_s*cov2(kPx,kX)+s_sq*cov2(kPx,kPx);
2634  cov2(kY,kY)+=two_s*cov2(kPy,kY)+s_sq*cov2(kPy,kPy);
2635  cov2(kZ,kZ)+=two_s*cov2(kPz,kZ)+s_sq*cov2(kPz,kPz);
2636  }
2637  var_doca=(dx*dx*((*cov)(kX,kX)+cov2(kX,kX))
2638  +dy*dy*((*cov)(kY,kY)+cov2(kY,kY))
2639  +dz*dz*((*cov)(kZ,kZ)+cov2(kZ,kZ))
2640  +2.*dx*dy*((*cov)(kX,kY)+cov2(kX,kY))
2641  +2.*dx*dz*((*cov)(kX,kZ)+cov2(kX,kZ))
2642  +2.*dy*dz*((*cov)(kY,kZ)+cov2(kY,kZ)))
2643  /(doca*doca);
2644  }
2645  break;
2646  }
2647  // New point on line
2648  if (move_along_line){
2649  point+=(step_size/pscale)*dir;
2650  s+=step_size;
2651  }
2652 
2653  // Propagate the covariance matrix of the track along the trajectory
2654  if(track_kd!=NULL)
2655  this->PropagateCovariance(this->swim_steps[i+1].s-swim_step->s,q,mass_sq,swim_step->mom,swim_step->origin,swim_step->B,*cov);
2656 
2657  // Store the current position and doca
2658  oldpos=pos;
2659  oldmom=swim_step->mom;
2660  tflight=swim_step->t;
2661  doca=new_doca;
2662  }
2663 
2664  // "Vertex" is mid-point of line connecting the positions of closest
2665  // approach of the two tracks
2666  commonpos = 0.5*(oldpos + point);
2667 
2668  return NOERROR;
2669 }
2670 
2671 // Find the position along a reference trajectory closest to a given point.
2672 // The error matrix for the point can also be input via a pointer. The error
2673 // matrix is expected to be 7x7, with the order {Px,Py,Pz,X,Y,Z,T}.
2674 // Outputs the kinematic data object (including the covariance) at this
2675 // position,and the doca and the variance on the doca.
2677  const DMatrixDSym *covpoint,
2678  DKinematicData *track_kd,
2679  double &doca, double &var_doca) const{
2680  if (track_kd==NULL) return RESOURCE_UNAVAILABLE;
2681 
2682  DVector3 dir, commonpos;
2683  return FindPOCAtoLine(point,dir,covpoint,track_kd,commonpos,doca,var_doca);
2684 }
2685 
2686 // Find the mid-point of the line connecting the points of closest approach of the
2687 // trajectories of two tracks. Return the positions, momenta, and error matrices
2688 // at these points for the two tracks.
2689 jerror_t DReferenceTrajectory::IntersectTracks(const DReferenceTrajectory *rt2, DKinematicData *track1_kd, DKinematicData *track2_kd, DVector3 &pos, double &doca, double &var_doca, double &vertex_chi2, bool DoFitVertex) const {
2690  const swim_step_t *swim_step1=this->swim_steps;
2691  const swim_step_t *swim_step2=rt2->swim_steps;
2692 
2693  TMatrixFSym cov1(7), cov2(7);
2694  shared_ptr<TMatrixFSym> locCovarianceMatrix1 = (track1_kd != NULL) ? dResourcePool_TMatrixFSym->Get_SharedResource() : nullptr;
2695  shared_ptr<TMatrixFSym> locCovarianceMatrix2 = (track2_kd != NULL) ? dResourcePool_TMatrixFSym->Get_SharedResource() : nullptr;
2696 
2697  if((track1_kd != NULL) && (track2_kd != NULL)){
2698  locCovarianceMatrix1->ResizeTo(7, 7);
2699  locCovarianceMatrix2->ResizeTo(7, 7);
2700  cov1=*track1_kd->errorMatrix();
2701  cov2=*track2_kd->errorMatrix();
2702  }
2703 
2704  double q1=this->q;
2705  double q2=rt2->q;
2706  double mass_sq1=this->mass_sq;
2707  double mass_sq2=rt2->mass_sq;
2708  vertex_chi2=0.;
2709 
2710  // Initialize the doca and traverse both particles' trajectories
2711  doca=1000.;
2712  double tflight1=0.,tflight2=0.;
2713  for (int i=0;i<this->Nswim_steps-1&&i<rt2->Nswim_steps-1; i++, swim_step1++, swim_step2++){
2714  DVector3 pos1=swim_step1->origin;
2715  DVector3 pos2=swim_step2->origin;
2716  DVector3 diff=pos1-pos2;
2717  double new_doca=diff.Mag();
2718 
2719  if (new_doca>doca){
2720  int prev_i=i-1;
2721  // positions and momenta of tracks at the center of the
2722  // bracketed region
2723  pos1=this->swim_steps[prev_i].origin;
2724  DVector3 mom1=this->swim_steps[prev_i].mom;
2725  pos2=rt2->swim_steps[prev_i].origin;
2726  DVector3 mom2=rt2->swim_steps[prev_i].mom;
2727 
2728  // If we break out of the loop immediately, we have not bracketed the
2729  // doca yet...
2730  if (i==1) { // backtrack to find the true doca
2731  tflight1=tflight2=0.;
2732  if((track1_kd != NULL) && (track2_kd != NULL)){
2733  cov1=*track1_kd->errorMatrix();
2734  cov2=*track2_kd->errorMatrix();
2735  }
2736  // Initialize the steppers
2737  DMagneticFieldStepper stepper1(this->bfield, q1, &pos1, &mom1);
2738  DMagneticFieldStepper stepper2(this->bfield, q2, &pos2, &mom2);
2739 
2740  // Do the backtracking...
2741  int inew=0;
2742  DVector3 oldpos1=pos1;
2743  DVector3 oldpos2=pos2;
2744  while (inew<20){
2745  if (pos1.z()<0. || pos2.z()<0. || pos1.z()>400. || pos2.z()>400.
2746  || pos1.Perp()>65. || pos2.Perp()>65.){
2747  break;
2748  }
2749  DVector3 B1,B2;
2750  double ds1=stepper1.Step(&pos1,&B1,-0.5);
2751  double ds2=stepper2.Step(&pos2,&B2,-0.5);
2752 
2753  // Compute the revised estimate for the doca
2754  diff=pos1-pos2;
2755  new_doca=diff.Mag();
2756 
2757  if(new_doca > doca){
2758  pos1=oldpos1;
2759  pos2=oldpos2;
2760  break;
2761  }
2762 
2763  // Propagate the covariance matrices along the trajectories
2764  if((track1_kd != NULL) && (track2_kd != NULL)){
2765  this->PropagateCovariance(ds1,q1,mass_sq1,mom1,oldpos1,B1,cov1);
2766  rt2->PropagateCovariance(ds2,q2,mass_sq2,mom2,oldpos2,B2,cov2);
2767  }
2768 
2769  // Store the current positions, doca and adjust flight times
2770  oldpos1=pos1;
2771  oldpos2=pos2;
2772  doca=new_doca;
2773 
2774  double one_over_p1_sq=1./mom1.Mag2();
2775  tflight1+=ds1*sqrt(1.+mass_sq1*one_over_p1_sq)/SPEED_OF_LIGHT;
2776 
2777  double one_over_p2_sq=1./mom2.Mag2();
2778  tflight2+=ds2*sqrt(1.+mass_sq2*one_over_p2_sq)/SPEED_OF_LIGHT;
2779 
2780  // New momenta
2781  stepper1.GetMomentum(mom1);
2782  stepper2.GetMomentum(mom2);
2783  }
2784  }
2785 
2786  // Use Brent's algorithm to find a better approximation for
2787  // the poca of the two tracks
2788  double ds=0.5;
2789  BrentsAlgorithm(pos1,mom1,pos2,mom2,ds,q2,doca);
2790 
2791  // "Vertex" is mid-point of line connecting the positions of closest
2792  // approach of the two tracks
2793  pos=0.5*(pos1+pos2);
2794 
2795  if((track1_kd != NULL) && (track2_kd != NULL)){
2796  if (DoFitVertex){
2797  // Use lagrange multiplier method to try to find a better common
2798  // vertex between the two tracks
2799  FitVertex(pos1,mom1,pos2,mom2,cov1,cov2,pos,vertex_chi2);
2800  }
2801  // Compute the variance on the doca
2802  diff=pos1-pos2;
2803  if (DoFitVertex) doca=diff.Mag(); // update if necessary
2804  double dx=diff.x();
2805  double dy=diff.y();
2806  double dz=diff.z();
2807  var_doca=(dx*dx*(cov1(kX,kX)+cov2(kX,kX))
2808  +dy*dy*(cov1(kY,kY)+cov2(kY,kY))
2809  +dz*dz*(cov1(kZ,kZ)+cov2(kZ,kZ))
2810  +2.*dx*dy*(cov1(kX,kY)+cov2(kX,kY))
2811  +2.*dx*dz*(cov1(kX,kZ)+cov2(kX,kZ))
2812  +2.*dy*dz*(cov1(kY,kZ)+cov2(kY,kZ)))
2813  /(doca*doca);
2814 
2815 
2816  // Adjust flight times
2817  double one_over_p1_sq=1./mom1.Mag2();
2818  tflight1+=ds*sqrt(1.+mass_sq1*one_over_p1_sq)/SPEED_OF_LIGHT;
2819 
2820  double one_over_p2_sq=1./mom2.Mag2();
2821  tflight2+=ds*sqrt(1.+mass_sq2*one_over_p2_sq)/SPEED_OF_LIGHT;
2822 
2823  *locCovarianceMatrix1 = cov1;
2824  track1_kd->setErrorMatrix(locCovarianceMatrix1);
2825  track1_kd->setMomentum(mom1);
2826  track1_kd->setPosition(pos1);
2827  track1_kd->setTime(track1_kd->time() + tflight1);
2828 
2829  *locCovarianceMatrix2 = cov2;
2830  track2_kd->setErrorMatrix(locCovarianceMatrix2);
2831  track2_kd->setMomentum(mom2);
2832  track2_kd->setPosition(pos2);
2833  track2_kd->setTime(track2_kd->time() + tflight2);
2834  }
2835  break;
2836  }
2837 
2838  // Propagate the covariance matrices along the trajectories
2839  if((track1_kd != NULL) && (track2_kd != NULL)){
2840  this->PropagateCovariance(this->swim_steps[i+1].s-swim_step1->s,q1,mass_sq1,swim_step1->mom,swim_step1->origin,swim_step1->B,cov1);
2841  rt2->PropagateCovariance(rt2->swim_steps[i+1].s-swim_step2->s,q2,mass_sq2,swim_step2->mom,swim_step2->origin,swim_step2->B,cov2);
2842  }
2843 
2844  // Store the current positions and doca
2845  tflight1=swim_step1->t;
2846  tflight2=swim_step2->t;
2847  doca=new_doca;
2848  }
2849 
2850  return NOERROR;
2851 }
2852 
2853 
2854 // Routine for finding the minimum of a function bracketed between two values
2855 // (see Numerical Recipes in C, pp. 404-405).
2856 #define ITMAX 20
2857 #define CGOLD 0.3819660
2858 #define EPS2 1.e-4
2859 #define ZEPS 1.0e-10
2860 #define SHFT(a,b,c,d) (a)=(b);(b)=(c);(c)=(d);
2861 #define SIGN(a,b) ((b)>=0.0?fabs(a):-fabs(a))
2863  DVector3 &pos2,DVector3 &mom2,
2864  double ds,double q2,
2865  double &doca) const{
2866  double d=0.,u=0.;
2867  double e=0.0; // will be distance moved on step before last
2868  double ax=0.;
2869  double bx=-ds;
2870  double cx=-2.*ds;
2871 
2872  double a=(ax<cx?ax:cx);
2873  double b=(ax>cx?ax:cx);
2874  double x=bx,w=bx,v=bx;
2875 
2876  // initialization
2877  double fw=doca;
2878  double fv=fw;
2879  double fx=fw;
2880  double u_old=x;
2881  DMagneticFieldStepper stepper1(this->bfield, this->q, &pos1, &mom1);
2882  DMagneticFieldStepper stepper2(this->bfield, q2, &pos2, &mom2);
2883 
2884  // main loop
2885  for (unsigned int iter=1;iter<=ITMAX;iter++){
2886  double xm=0.5*(a+b);
2887  double tol1=EPS2*fabs(x)+ZEPS;
2888  double tol2=2.0*tol1;
2889  if (fabs(x-xm)<=(tol2-0.5*(b-a))){
2890  doca=(pos1-pos2).Mag();
2891 
2892  // New momenta
2893  stepper1.GetMomentum(mom1);
2894  stepper2.GetMomentum(mom2);
2895 
2896  return NOERROR;
2897  }
2898  // trial parabolic fit
2899  if (fabs(e)>tol1){
2900  double x_minus_w=x-w;
2901  double x_minus_v=x-v;
2902  double r=x_minus_w*(fx-fv);
2903  double q=x_minus_v*(fx-fw);
2904  double p=x_minus_v*q-x_minus_w*r;
2905  q=2.0*(q-r);
2906  if (q>0.0) p=-p;
2907  q=fabs(q);
2908  double etemp=e;
2909  e=d;
2910  if (fabs(p)>=fabs(0.5*q*etemp) || p<=q*(a-x) || p>=q*(b-x))
2911  // fall back on the Golden Section technique
2912  d=CGOLD*(e=(x>=xm?a-x:b-x));
2913  else{
2914  // parabolic step
2915  d=p/q;
2916  u=x+d;
2917  if (u-a<tol2 || b-u <tol2)
2918  d=SIGN(tol1,xm-x);
2919  }
2920  } else{
2921  d=CGOLD*(e=(x>=xm?a-x:b-x));
2922  }
2923  u=(fabs(d)>=tol1 ? x+d: x+SIGN(tol1,d));
2924 
2925  // Function evaluation
2926  double du=u_old-u;
2927  stepper1.Step(&pos1,NULL,du);
2928  stepper2.Step(&pos2,NULL,du);
2929  DVector3 diff=pos1-pos2;
2930  double fu=diff.Mag();
2931  u_old=u;
2932 
2933  if (fu<=fx){
2934  if (u>=x) a=x; else b=x;
2935  SHFT(v,w,x,u);
2936  SHFT(fv,fw,fx,fu);
2937  }
2938  else {
2939  if (u<x) a=u; else b=u;
2940  if (fu<=fw || w==x){
2941  v=w;
2942  w=u;
2943  fv=fw;
2944  fw=fu;
2945  }
2946  else if (fu<=fv || v==x || v==w){
2947  v=u;
2948  fv=fu;
2949  }
2950  }
2951  }
2952 
2953  // We only get here if there is a convergence issue...
2954  doca=(pos1-pos2).Mag();
2955  stepper1.GetMomentum(mom1);
2956  stepper2.GetMomentum(mom2);
2957 
2958  return NOERROR;
2959 }
2960 
2961 
2962 // Use lagrange multiplier method to find better approximation for the vertex
2963 // position given the track momenta and positions and the corresponding error
2964 // matrices. See Frodesen, et al., Probability and Statistics in Particle
2965 // Physics, pp 298-320.
2966 // Assumes we can ignore the magnetic field.
2968  const DVector3 &pos2,const DVector3 &mom2,
2969  const TMatrixFSym &cov1,
2970  const TMatrixFSym &cov2,
2971  DVector3 &pos,double &vertex_chi2,double q1,double q2) const{
2972  // Vectors of measured quantities (eta0) and these quantities after adjustment
2973  // based on constraints (eta)
2974  TMatrixD eta0(12,1),eta(12,1);
2975  // Vectors of unmeasured quantities (common vertex {x,y,z})
2976  TMatrixD xi(3,1),xi_old(3,1);
2977  TMatrixD f(4,1); // vector representing constraint equations
2978  TMatrixD F_eta(4,12); // Matrix of derivatives of f with respect to eta
2979  TMatrixD F_xi(4,3); // Matrix of derivatives of f with respect to xi
2980 
2981  // Initialize to the guess to the vertex position
2982  xi(0,0)=pos.x();
2983  xi(1,0)=pos.y();
2984  xi(2,0)=pos.z();
2985 
2986  // B-field
2987  double Bz=fabs(bfield->GetBz(xi(0,0),xi(1,0),xi(2,0)));
2988  bool got_field=false;
2989  double a1=0.,a2=0.;
2990  if (Bz>0.01){
2991  a1=0.003*Bz*q1;
2992  a2=0.003*Bz*q2;
2993  got_field=true;
2994  }
2995 
2996  // Initialize to the measured values
2997  eta0(kX,0)=pos1.x();
2998  eta0(kY,0)=pos1.y();
2999  eta0(kZ,0)=pos1.z();
3000  eta0(kPx,0)=mom1.x();
3001  eta0(kPy,0)=mom1.y();
3002  eta0(kPz,0)=mom1.z();
3003  eta0(kX+6,0)=pos2.x();
3004  eta0(kY+6,0)=pos2.y();
3005  eta0(kZ+6,0)=pos2.z();
3006  eta0(kPx+6,0)=mom2.x();
3007  eta0(kPy+6,0)=mom2.y();
3008  eta0(kPz+6,0)=mom2.z();
3009 
3010  // Fill error matrix V from covariance matrices from the two tracks
3011  TMatrixD V(12,12);
3012  for (int m=0;m<6;m++){
3013  for (int n=0;n<6;n++){
3014  V(m,n)=cov1(n,m);
3015  V(m+6,n+6)=cov2(n,m);
3016  }
3017  }
3018 
3019  // Start with the actual measurement
3020  eta=eta0;
3021 
3022  // Define some matrices needed for some of the internal steps of the procedure
3023  TMatrixD LagMul(4,1); //Lagrange multipliers
3024  TMatrixD r(4,1); // r=f+F_eta*(eta0-eta)
3025  TMatrixD S(4,4); // S=F_eta*V*F_eta^T
3026 
3027  if (debug_level>1){
3028  cout << "Measured quantities:" <<endl;
3029  eta.Print();
3030  cout << "Covariance matrix for measured quantities:"<<endl;
3031  V.Print();
3032  cout << "Initial guess for common vertex position:"<<endl;
3033  xi.Print();
3034  }
3035  // Iterate to find the minimum chi^2
3036  double chi2=1.e6,chi2_old=1e6;
3037  for (int iter=0;iter<4;iter++){
3038  chi2_old=chi2;
3039 
3040  double dx1=xi(0,0)-eta(kX,0);
3041  double dy1=xi(1,0)-eta(kY,0);
3042  double dz1=xi(2,0)-eta(kZ,0);
3043  double dx2=xi(0,0)-eta(kX+6,0);
3044  double dy2=xi(1,0)-eta(kY+6,0);
3045  double dz2=xi(2,0)-eta(kZ+6,0);
3046  double px1=eta(kPx,0);
3047  double py1=eta(kPy,0);
3048  double pz1=eta(kPz,0);
3049  double px2=eta(kPx+6,0);
3050  double py2=eta(kPy+6,0);
3051  double pz2=eta(kPz+6,0);
3052 
3053  // Constraint equations
3054  if (got_field==false){
3055  // No field...
3056  f(0,0)=pz1*dx1-px1*dz1;
3057  f(1,0)=pz1*dy1-py1*dz1;
3058  f(2,0)=pz2*dx2-px2*dz2;
3059  f(3,0)=pz2*dy2-py2*dz2;
3060  }
3061  else{
3062  double twoks1=0.;
3063  if (fabs(pz1)>1e-8) twoks1=a1*dz1/pz1;
3064  double one_minus_cos_2ks1=1.-cos(twoks1);
3065  double sin_2ks1=sin(twoks1);
3066  f(0,0)=a1*dx1-px1*sin_2ks1+py1*one_minus_cos_2ks1;
3067  f(1,0)=a1*dy1-py1*sin_2ks1-px1*one_minus_cos_2ks1;
3068 
3069  double twoks2=0.;
3070  if (fabs(pz2)>1e-8) twoks2=a2*dz2/pz2;
3071  double one_minus_cos_2ks2=1.-cos(twoks2);
3072  double sin_2ks2=sin(twoks2);
3073  f(2,0)=a2*dx2-px2*sin_2ks2+py2*one_minus_cos_2ks2;
3074  f(3,0)=a2*dy2-py2*sin_2ks2-px2*one_minus_cos_2ks2;
3075 
3076  }
3077 
3078  if (iter>0){
3079  TMatrixD LagMul_T(TMatrixD::kTransposed,LagMul);
3080  chi2=(LagMul_T*S*LagMul)(0,0);
3081  LagMul_T*=2.;
3082  chi2+=(LagMul_T*f)(0,0);
3083 
3084  if (debug_level>0){
3085  cout << "iter " << iter << " : chi2=" << chi2 << endl;
3086  if (debug_level>1){
3087  cout << "Constraint equations:" << endl;
3088  f.Print();
3089  cout << "Unmeasured quantities:" << endl;
3090  xi.Print();
3091  cout << "Measured quantities:" << endl;
3092  eta.Print();
3093  }
3094  }
3095 
3096  //if (chi2>chi2_old) break;
3097  }
3098 
3099  if (got_field==false){
3100  // Compute the Jacobian matrix for eta
3101  F_eta(0,kX)=-pz1;
3102  F_eta(0,kZ)=+px1;
3103  F_eta(0,kPx)=-dz1;
3104  F_eta(0,kPz)=+dx1;
3105  F_eta(1,kY)=-pz1;
3106  F_eta(1,kZ)=+py1;
3107  F_eta(1,kPy)=-dz1;
3108  F_eta(1,kPz)=+dy1;
3109  F_eta(2,kX+6)=-pz2;
3110  F_eta(2,kZ+6)=+px2;
3111  F_eta(2,kPx+6)=-dz2;
3112  F_eta(2,kPz+6)=+dx2;
3113  F_eta(3,kY+6)=-pz2;
3114  F_eta(3,kZ+6)=+py2;
3115  F_eta(3,kPy+6)=-dz2;
3116  F_eta(3,kPz+6)=+dy2;
3117 
3118  // Compute the Jacobian matrix for xi
3119  F_xi(0,0)=pz1;
3120  F_xi(0,2)=-px1;
3121  F_xi(1,1)=pz1;
3122  F_xi(1,2)=-py1;
3123  F_xi(2,0)=pz2;
3124  F_xi(2,2)=-px2;
3125  F_xi(3,1)=pz2;
3126  F_xi(3,2)=-py2;
3127  }
3128  else{
3129  // Compute the Jacobian matrix for eta
3130  double twoks1=0.;
3131  if (fabs(pz1)>1e-8) twoks1=a1*dz1/pz1;
3132  double cos_2ks1=cos(twoks1);
3133  double one_minus_cos_2ks1=1.-cos_2ks1;
3134  double sin_2ks1=sin(twoks1);
3135  F_eta(0,kX)=-a1;
3136  F_eta(0,kZ)=(-a1/pz1)*(py1*sin_2ks1-px1*cos_2ks1);
3137  F_eta(0,kPx)=-sin_2ks1;
3138  F_eta(0,kPy)=one_minus_cos_2ks1;
3139  F_eta(0,kPz)=(a1*dz1/(pz1*pz1))*(px1*cos_2ks1-py1*sin_2ks1);
3140  F_eta(1,kY)=-a1;
3141  F_eta(1,kZ)=(a1/pz1)*(py1*cos_2ks1+px1*sin_2ks1);
3142  F_eta(1,kPy)=-sin_2ks1;
3143  F_eta(1,kPx)=-one_minus_cos_2ks1;
3144  F_eta(1,kPz)=(a1*dz1/(pz1*pz1))*(px1*cos_2ks1-py1*sin_2ks1);
3145  double twoks2=0.;
3146  if (fabs(pz2)>1e-8) twoks2=a2*dz2/pz2;
3147  double cos_2ks2=cos(twoks2);
3148  double one_minus_cos_2ks2=1.-cos_2ks2;
3149  double sin_2ks2=sin(twoks2);
3150  F_eta(2,kX+6)=-a2;
3151  F_eta(2,kZ+6)=(-a2/pz2)*(py2*sin_2ks2-px2*cos_2ks2);
3152  F_eta(2,kPx+6)=-sin_2ks2;
3153  F_eta(2,kPy+6)=one_minus_cos_2ks2;
3154  F_eta(2,kPz+6)=(a2*dz2/(pz2*pz2))*(px2*cos_2ks2-py2*sin_2ks2);
3155  F_eta(3,kY+6)=-a2;
3156  F_eta(3,kZ+6)=(a2/pz2)*(py2*cos_2ks2+px2*sin_2ks2);
3157  F_eta(3,kPy+6)=-sin_2ks2;
3158  F_eta(3,kPx+6)=-one_minus_cos_2ks2;
3159  F_eta(3,kPz+6)=(a2*dz2/(pz2*pz2))*(px2*cos_2ks2-py2*sin_2ks2);
3160 
3161  // Compute the Jacobian matrix for xi
3162  F_xi(0,0)=-F_eta(0,kX);
3163  F_xi(0,2)=-F_eta(0,kZ);
3164  F_xi(1,1)=-F_eta(1,kY);
3165  F_xi(1,2)=-F_eta(1,kZ);
3166  F_xi(2,0)=-F_eta(2,kX+6);
3167  F_xi(2,2)=-F_eta(2,kZ+6);
3168  F_xi(3,1)=-F_eta(3,kY+6);
3169  F_xi(3,2)=-F_eta(3,kZ+6);
3170 
3171  }
3172 
3173  if (debug_level>1){
3174  cout << "Jacobian matrices for eta and xi" <<endl;
3175  F_eta.Print();
3176  F_xi.Print();
3177  }
3178 
3179  TMatrixD F_eta_T(TMatrixD::kTransposed,F_eta);
3180  TMatrixD F_xi_T(TMatrixD::kTransposed,F_xi);
3181 
3182  // Adjust xi, then the Lagrange multipliers, then eta for this iteration
3183  xi_old=xi;
3184  S=F_eta*V*F_eta_T;
3185  TMatrixD Sinv(TMatrixD::kInverted,S);
3186  r=f+F_eta*(eta0-eta);
3187  TMatrixD Mtemp(TMatrixD::kInverted,F_xi_T*Sinv*F_xi);
3188  xi=xi_old-Mtemp*F_xi_T*Sinv*r;
3189  LagMul=Sinv*(r+F_xi*(xi-xi_old));
3190  eta=eta0-V*F_eta_T*LagMul;
3191  }
3192  pos.SetXYZ(xi_old(0,0),xi_old(1,0),xi_old(2,0));
3193  vertex_chi2=chi2_old;
3194 }
#define F(x, y, z)
void GetBField(DVector3 &B)
void setMomentum(const DVector3 &aMomentum)
void setTime(double locTime)
int InsertSteps(const swim_step_t *start_step, double delta_s, double step_size=0.02)
jerror_t BrentsAlgorithm(DVector3 &pos1, DVector3 &mom1, DVector3 &pos2, DVector3 &mom2, double ds, double q2, double &doca) const
double Straw_dx(const DCoordinateSystem *wire, double radius) const
bool GetCheckMaterialBoundaries(void) const
#define SPEED_OF_LIGHT
double GetBoundaryStepFraction(void) const
TMatrixD DMatrix
Definition: DMatrix.h:14
void FastSwim(const DVector3 &pos, const DVector3 &mom, double q, double smax=2000.0, double zmin=-100., double zmax=1000.0)
#define X0
void CopyWithShift(const DReferenceTrajectory *rt, DVector3 shift)
TVector2 DVector2
Definition: DVector2.h:9
TVector3 DVector3
Definition: DVector3.h:14
#define EPS2
const DMagneticFieldMap * bfield
Double_t x[NCHANNELS]
Definition: st_tw_resols.C:39
#define ITMAX
#define c
#define y
DVector3 GetLastDOCAPoint(void) const
swim_step_t * FindPlaneCrossing(const DVector3 &origin, DVector3 norm, int first_i=0, DetectorSystem_t detector=SYS_NULL) const
void Swim(const DVector3 &pos, const DVector3 &mom, double q=-1000.0, const TMatrixFSym *cov=NULL, double smax=2000.0, const DCoordinateSystem *wire=NULL)
double dPdx(double ptot, double KrhoZ_overA, double rhoZ_overA, double LogI) const
#define CGOLD
double zmax
DetectorSystem_t
Definition: GlueX.h:15
TMatrixDSym DMatrixDSym
Definition: DMatrixDSym.h:13
#define TWO_THIRD
#define ZEPS
void GetPosMom(DVector3 &pos, DVector3 &mom)
static char index(char c)
Definition: base64.cpp:115
jerror_t IntersectTracks(const DReferenceTrajectory *rt2, DKinematicData *track1_kd, DKinematicData *track2_kd, DVector3 &pos, double &doca, double &var_doca, double &vertex_chi2, bool DoFitVertex=false) const
#define X(str)
Definition: hddm-c.cpp:83
#define QuietNaN
DMagneticFieldStepper class.
swim_step_t * FindClosestSwimStep(const DCoordinateSystem *wire, int *istep_ptr=NULL) const
Definition: GlueX.h:19
TF1 * f
Definition: FitGains.C:21
double GetMaxStepSize(void) const
double time(void) const
jerror_t SetStartingParams(double q, const DVector3 *x, const DVector3 *p)
void setErrorMatrix(const shared_ptr< const TMatrixFSym > &aMatrix)
TEllipse * e
void FastSwimForHitSelection(const DVector3 &pos, const DVector3 &mom, double q)
double DistToRT(double x, double y, double z) const
const double alpha
#define H(x, y, z)
Definition: GlueX.h:20
double DistToRTBruteForce(const DCoordinateSystem *wire, double *s=NULL) const
double DistToRTwithTime(DVector3 hit, double *s=NULL, double *t=NULL, double *var_t=NULL, DetectorSystem_t detector=SYS_NULL) const
jerror_t FindPOCAtoPoint(const DVector3 &point, const DMatrixDSym *covpoint, DKinematicData *track_kd, double &doca, double &var_doca) const
Definition: GlueX.h:22
void GetMomentum(DVector3 &mom)
void Dump(double zmin=-1000.0, double zmax=1000.0)
#define ONE_THIRD
#define _DBG_
Definition: HDEVIO.h:12
DReferenceTrajectory::swim_step_t steps[256]
double GetMass(void) const
double FastStep(double stepsize=0.0)
double dPdx_from_A_Z_rho(double ptot, double A, double Z, double density) const
#define _DBG__
Definition: HDEVIO.h:13
TH1D * py[NCHANNELS]
void FitVertex(const DVector3 &pos1, const DVector3 &mom1, const DVector3 &pos2, const DVector3 &mom2, const TMatrixFSym &cov1, const TMatrixFSym &cov2, DVector3 &pos, double &vertex_chi2, double q1=1., double q2=1.) const
jerror_t GetIntersectionWithRadius(double R, DVector3 &mypos, double *s=NULL, double *t=NULL, DVector3 *dir=NULL) const
double sqrt(double)
double sin(double)
#define S(str)
Definition: hddm-c.cpp:84
DReferenceTrajectory & operator=(const DReferenceTrajectory &rt)
const swim_step_t * last_swim_step
last swim step used in DistToRT
double last_phi
last phi found in DistToRT
double GetMinStepSize(void) const
jerror_t FindPOCAtoLine(const DVector3 &origin, const DVector3 &dir, const DMatrixDSym *covpoint, DKinematicData *track_kd, DVector3 &commonpos, double &doca, double &var_doca) const
#define MIN_STEP_SIZE
double Step(DVector3 *newpos=NULL, DVector3 *B=NULL, double stepsize=0.0)
shared_ptr< const TMatrixFSym > errorMatrix(void) const
#define SHFT(a, b, c, d)
void setPosition(const DVector3 &aPosition)
TDirectory * dir
Definition: bcal_hist_eff.C:31
double zmin
const DRootGeom * RootGeom
jerror_t SetStepSize(double step)
TH2F * temp
union @6::@8 u
void SetStepSize(double step_size)
jerror_t GetIntersectionWithPlane(const DVector3 &origin, const DVector3 &norm, DVector3 &pos, double *s=NULL, double *t=NULL, double *var_t=NULL, DetectorSystem_t detector=SYS_NULL) const
#define SIGN(a, b)
#define I(x, y, z)
#define G(x, y, z)
void GetDirs(DVector3 &xdir, DVector3 &ydir, DVector3 &zdir)
#define EPS
jerror_t PropagateCovariance(double ds, double q, double mass_sq, const DVector3 &mom, const DVector3 &pos, const DVector3 &B, TMatrixFSym &C) const