TrajectoryMCSFitter.cxx
Go to the documentation of this file.
1 #include "TrajectoryMCSFitter.h"
5 #include "TMatrixDSym.h"
6 #include "TMatrixDSymEigen.h"
7 
8 using namespace std;
9 using namespace trkf;
10 using namespace recob::tracking;
11 
12 recob::MCSFitResult TrajectoryMCSFitter::fitMcs(const recob::TrackTrajectory& traj, int pid) const {
13  //
14  // Break the trajectory in segments of length approximately equal to segLen_
15  //
16  vector<size_t> breakpoints;
17  vector<float> segradlengths;
18  vector<float> cumseglens;
19  breakTrajInSegments(traj, breakpoints, segradlengths, cumseglens);
20  //
21  // Fit segment directions, and get 3D angles between them
22  //
23  if (segradlengths.size()<2) return recob::MCSFitResult();
24  vector<float> dtheta;
25  Vector_t pcdir0;
26  Vector_t pcdir1;
27  for (unsigned int p = 0; p<segradlengths.size(); p++) {
28  linearRegression(traj, breakpoints[p], breakpoints[p+1], pcdir1);
29  if (p>0) {
30  if (segradlengths[p]<-100. || segradlengths[p-1]<-100.) {
31  dtheta.push_back(-999.);
32  } else {
33  const double cosval = pcdir0.X()*pcdir1.X()+pcdir0.Y()*pcdir1.Y()+pcdir0.Z()*pcdir1.Z();
34  //assert(std::abs(cosval)<=1);
35  //units are mrad
36  double dt = 1000.*acos(cosval);//should we try to use expansion for small angles?
37  dtheta.push_back(dt);
38  }
39  }
40  pcdir0 = pcdir1;
41  }
42  //
43  // Perform likelihood scan in forward and backward directions
44  //
45  vector<float> cumLenFwd;
46  vector<float> cumLenBwd;
47  for (unsigned int i = 0; i<cumseglens.size()-2; i++) {
48  cumLenFwd.push_back(cumseglens[i]);
49  cumLenBwd.push_back(cumseglens.back()-cumseglens[i+2]);
50  }
51  double detAngResol = DetectorAngularResolution(std::abs(traj.StartDirection().Z()));
52  const ScanResult fwdResult = doLikelihoodScan(dtheta, segradlengths, cumLenFwd, true , pid, detAngResol);
53  const ScanResult bwdResult = doLikelihoodScan(dtheta, segradlengths, cumLenBwd, false, pid, detAngResol);
54  //
55  return recob::MCSFitResult(pid,
56  fwdResult.p,fwdResult.pUnc,fwdResult.logL,
57  bwdResult.p,bwdResult.pUnc,bwdResult.logL,
58  segradlengths,dtheta);
59 }
60 
61 void TrajectoryMCSFitter::breakTrajInSegments(const recob::TrackTrajectory& traj, vector<size_t>& breakpoints, vector<float>& segradlengths, vector<float>& cumseglens) const {
62  //
64  auto const* _SCE = (applySCEcorr_ ? lar::providerFrom<spacecharge::SpaceChargeService>() : NULL);
65  //
66  const double trajlen = traj.Length();
67  const double thisSegLen = (trajlen>(segLen_*minNSegs_) ? segLen_ : trajlen/double(minNSegs_) );
68  // std::cout << "track with length=" << trajlen << " broken in nseg=" << std::max(minNSegs_,int(trajlen/segLen_)) << " of length=" << thisSegLen << " where segLen_=" << segLen_ << std::endl;
69  //
70  constexpr double lar_radl_inv = 1./14.0;
71  cumseglens.push_back(0.);//first segment has zero cumulative length from previous segments
72  double thislen = 0.;
73  double totlen = 0.;
74  auto nextValid=traj.FirstValidPoint();
75  breakpoints.push_back(nextValid);
76  auto pos0 = traj.LocationAtPoint(nextValid);
77  if (applySCEcorr_) {
78  const double Position[3] = {pos0.X(), pos0.Y(), pos0.Z()};
79  geo::TPCID tpcid = geom->FindTPCAtPosition(Position);
80  geo::Vector_t pos0_offset = geo::Vector_t(0., 0., 0.);
81  if(tpcid.isValid) {
82  geo::Point_t p0(pos0.X(), pos0.Y(), pos0.Z());
83  pos0_offset = _SCE->GetCalPosOffsets(p0, tpcid.TPC);
84  }
85  pos0.SetX(pos0.X() - pos0_offset.X());
86  pos0.SetY(pos0.Y() + pos0_offset.Y());
87  pos0.SetZ(pos0.Z() + pos0_offset.Z());
88  }
89  auto dir0 = traj.DirectionAtPoint(nextValid);
90  nextValid = traj.NextValidPoint(nextValid+1);
91  int npoints = 0;
92  while (nextValid!=recob::TrackTrajectory::InvalidIndex) {
93  if (npoints==0) dir0 = traj.DirectionAtPoint(nextValid);
94  auto pos1 = traj.LocationAtPoint(nextValid);
95  if (applySCEcorr_) {
96  const double Position[3] = {pos1.X(), pos1.Y(), pos1.Z()};
97  geo::TPCID tpcid = geom->FindTPCAtPosition(Position);
98  geo::Vector_t pos1_offset = geo::Vector_t(0., 0., 0.);
99  if(tpcid.isValid) {
100  geo::Point_t p1(pos1.X(), pos1.Y(), pos1.Z());
101  pos1_offset = _SCE->GetCalPosOffsets(p1, tpcid.TPC);
102  }
103  pos1.SetX(pos1.X() - pos1_offset.X());
104  pos1.SetY(pos1.Y() + pos1_offset.Y());
105  pos1.SetZ(pos1.Z() + pos1_offset.Z());
106  }
107  //increments along the initial direction of the segment
108  auto step = (pos1-pos0).R();
109  thislen += dir0.Dot(pos1-pos0);
110  totlen += step;
111  pos0=pos1;
112  // //fixme: testing alternative approaches here
113  // //test1: increments following scatters
114  // auto step = (pos1-pos0).R();
115  // thislen += step;
116  // totlen += step;
117  // pos0=pos1;
118  // //test2: end-start distance along the initial direction of the segment
119  // thislen = dir0.Dot(pos1-pos0);
120  // totlen = (pos1-pos0).R();
121  //
122  npoints++;
123  if (thislen>=(thisSegLen-segLenTolerance_)) {
124  breakpoints.push_back(nextValid);
125  if (npoints>=minHitsPerSegment_) segradlengths.push_back(thislen*lar_radl_inv);
126  else segradlengths.push_back(-999.);
127  cumseglens.push_back(totlen);
128  thislen = 0.;
129  npoints = 0;
130  }
131  nextValid = traj.NextValidPoint(nextValid+1);
132  }
133  //then add last segment
134  if (thislen>0.) {
135  breakpoints.push_back(traj.LastValidPoint()+1);
136  segradlengths.push_back(thislen*lar_radl_inv);
137  cumseglens.push_back(cumseglens.back()+thislen);
138  }
139  return;
140 }
141 
142 const TrajectoryMCSFitter::ScanResult TrajectoryMCSFitter::doLikelihoodScan(std::vector<float>& dtheta, std::vector<float>& seg_nradlengths, std::vector<float>& cumLen,
143  bool fwdFit, int pid, float pmin, float pmax, float pstep, float detAngResol) const {
144  int best_idx = -1;
145  float best_logL = std::numeric_limits<float>::max();
146  float best_p = -1.0;
147  std::vector<float> vlogL;
148  for (float p_test = pmin; p_test <= pmax; p_test+=pstep) {
149  float logL = mcsLikelihood(p_test, detAngResol, dtheta, seg_nradlengths, cumLen, fwdFit, pid);
150  if (logL < best_logL) {
151  best_p = p_test;
152  best_logL = logL;
153  best_idx = vlogL.size();
154  }
155  vlogL.push_back(logL);
156  }
157  //
158  //uncertainty from left side scan
159  float lunc = -1.;
160  if (best_idx>0) {
161  for (int j=best_idx-1;j>=0;j--) {
162  float dLL = vlogL[j]-vlogL[best_idx];
163  if ( dLL>=0.5 ) {
164  lunc = (best_idx-j)*pstep;
165  break;
166  }
167  }
168  }
169  //uncertainty from right side scan
170  float runc = -1.;
171  if (best_idx<int(vlogL.size()-1)) {
172  for (unsigned int j=best_idx+1;j<vlogL.size();j++) {
173  float dLL = vlogL[j]-vlogL[best_idx];
174  if ( dLL>=0.5 ) {
175  runc = (j-best_idx)*pstep;
176  break;
177  }
178  }
179  }
180  return ScanResult(best_p, std::max(lunc,runc), best_logL);
181 }
182 
183 const TrajectoryMCSFitter::ScanResult TrajectoryMCSFitter::doLikelihoodScan(std::vector<float>& dtheta, std::vector<float>& seg_nradlengths, std::vector<float>& cumLen,
184  bool fwdFit, int pid, float detAngResol) const {
185 
186  //do a first, coarse scan
187  const ScanResult& coarseRes = doLikelihoodScan(dtheta, seg_nradlengths, cumLen, fwdFit, pid, pMin_, pMax_, pStepCoarse_, detAngResol);
188 
189  float pmax = std::min(coarseRes.p+fineScanWindow_,pMax_);
190  float pmin = std::max(coarseRes.p-fineScanWindow_,pMin_);
191  if (coarseRes.pUnc < (std::numeric_limits<float>::max()-1.)) {
192  pmax = std::min(coarseRes.p+2*coarseRes.pUnc,pMax_);
193  pmin = std::max(coarseRes.p-2*coarseRes.pUnc,pMin_);
194  }
195 
196  //do the fine grained scan in a smaller region
197  const ScanResult& refineRes = doLikelihoodScan(dtheta, seg_nradlengths, cumLen, fwdFit, pid, pmin, pmax, pStep_, detAngResol);
198 
199  return refineRes;
200 }
201 
202 void TrajectoryMCSFitter::linearRegression(const recob::TrackTrajectory& traj, const size_t firstPoint, const size_t lastPoint, Vector_t& pcdir) const {
203  //
205  auto const* _SCE = (applySCEcorr_ ? lar::providerFrom<spacecharge::SpaceChargeService>() : NULL);
206  //
207  int npoints = 0;
208  geo::vect::MiddlePointAccumulator middlePointCalc;
209  size_t nextValid = firstPoint;
210  //fixme explore a max number of points to use for linear regression
211  //while (nextValid<std::min(firstPoint+10,lastPoint)) {
212  while (nextValid<lastPoint) {
213  auto tempP = traj.LocationAtPoint(nextValid);
214  if (applySCEcorr_) {
215  const double Position[3] = {tempP.X(), tempP.Y(), tempP.Z()};
216  geo::TPCID tpcid = geom->FindTPCAtPosition(Position);
217  geo::Vector_t tempP_offset = geo::Vector_t(0., 0., 0.);
218  if(tpcid.isValid) {
219  geo::Point_t ptemp(tempP.X(), tempP.Y(), tempP.Z());
220  tempP_offset = _SCE->GetCalPosOffsets(ptemp, tpcid.TPC);
221  }
222  tempP.SetX(tempP.X() - tempP_offset.X());
223  tempP.SetY(tempP.Y() + tempP_offset.Y());
224  tempP.SetZ(tempP.Z() + tempP_offset.Z());
225  }
226  middlePointCalc.add(tempP);
227  //middlePointCalc.add(traj.LocationAtPoint(nextValid));
228  nextValid = traj.NextValidPoint(nextValid+1);
229  npoints++;
230  }
231  const auto avgpos = middlePointCalc.middlePoint();
232  const double norm = 1./double(npoints);
233  //
234  //assert(npoints>0);
235  //
236  TMatrixDSym m(3);
237  nextValid = firstPoint;
238  while (nextValid<lastPoint) {
239  auto p = traj.LocationAtPoint(nextValid);
240  if (applySCEcorr_) {
241  const double Position[3] = {p.X(), p.Y(), p.Z()};
242  geo::TPCID tpcid = geom->FindTPCAtPosition(Position);
243  geo::Vector_t p_offset = geo::Vector_t(0., 0., 0.);
244  if(tpcid.isValid) {
245  geo::Point_t point(p.X(), p.Y(), p.Z());
246  p_offset = _SCE->GetCalPosOffsets(point, tpcid.TPC);
247  }
248  p.SetX(p.X() - p_offset.X());
249  p.SetY(p.Y() + p_offset.Y());
250  p.SetZ(p.Z() + p_offset.Z());
251  }
252  const double xxw0 = p.X()-avgpos.X();
253  const double yyw0 = p.Y()-avgpos.Y();
254  const double zzw0 = p.Z()-avgpos.Z();
255  m(0, 0) += xxw0*xxw0*norm;
256  m(0, 1) += xxw0*yyw0*norm;
257  m(0, 2) += xxw0*zzw0*norm;
258  m(1, 0) += yyw0*xxw0*norm;
259  m(1, 1) += yyw0*yyw0*norm;
260  m(1, 2) += yyw0*zzw0*norm;
261  m(2, 0) += zzw0*xxw0*norm;
262  m(2, 1) += zzw0*yyw0*norm;
263  m(2, 2) += zzw0*zzw0*norm;
264  nextValid = traj.NextValidPoint(nextValid+1);
265  }
266  //
267  const TMatrixDSymEigen me(m);
268  const auto& eigenval = me.GetEigenValues();
269  const auto& eigenvec = me.GetEigenVectors();
270  //
271  int maxevalidx = 0;
272  double maxeval = eigenval(0);
273  for (int i=1; i<3; ++i) {
274  if (eigenval(i)>maxeval) {
275  maxevalidx = i;
276  maxeval = eigenval(i);
277  }
278  }
279  //
280  pcdir = Vector_t(eigenvec(0, maxevalidx),eigenvec(1, maxevalidx),eigenvec(2, maxevalidx));
281  if (traj.DirectionAtPoint(firstPoint).Dot(pcdir)<0.) pcdir*=-1.;
282  //
283 }
284 
285 double TrajectoryMCSFitter::mcsLikelihood(double p, double theta0x, std::vector<float>& dthetaij, std::vector<float>& seg_nradl, std::vector<float>& cumLen, bool fwd, int pid) const {
286  //
287  const int beg = (fwd ? 0 : (dthetaij.size()-1));
288  const int end = (fwd ? dthetaij.size() : -1);
289  const int incr = (fwd ? +1 : -1);
290  //
291  // bool print = false;
292  //
293  const double m = mass(pid);
294  const double m2 = m*m;
295  const double Etot = sqrt(p*p + m2);//Initial energy
296  double Eij2 = 0.;
297  //
298  double const fixedterm = 0.5 * std::log( 2.0 * M_PI );
299  double result = 0;
300  for (int i = beg; i != end; i+=incr ) {
301  if (dthetaij[i]<0) {
302  //cout << "skip segment with too few points" << endl;
303  continue;
304  }
305  //
306  const double Eij = GetE(Etot,cumLen[i],m);
307  Eij2 = Eij*Eij;
308  //
309  if ( Eij2 <= m2 ) {
311  break;
312  }
313  const double pij = sqrt(Eij2 - m2);//momentum at this segment
314  const double beta = sqrt( 1. - ((m2)/(pij*pij + m2)) );
315  constexpr double HighlandSecondTerm = 0.038;
316  const double tH0 = ( HighlandFirstTerm(pij) / (pij*beta) ) * ( 1.0 + HighlandSecondTerm * std::log( seg_nradl[i] ) ) * sqrt( seg_nradl[i] );
317  const double rms = sqrt( 2.0*( tH0 * tH0 + theta0x * theta0x ) );
318  if (rms==0.0) {
319  std::cout << " Error : RMS cannot be zero ! " << std::endl;
321  }
322  const double arg = dthetaij[i]/rms;
323  result += ( std::log( rms ) + 0.5 * arg * arg + fixedterm);
324  }
325  return result;
326 }
327 
328 double TrajectoryMCSFitter::energyLossLandau(const double mass2,const double e2, const double x) const {
329  //
330  // eq. (33.11) in http://pdg.lbl.gov/2016/reviews/rpp2016-rev-passage-particles-matter.pdf (except density correction is ignored)
331  //
332  if (x<=0.) return 0.;
333  constexpr double Iinv2 = 1./(188.E-6*188.E-6);
334  constexpr double matConst = 1.4*18./40.;//density*Z/A
335  constexpr double me = 0.511;
336  constexpr double kappa = 0.307075;
337  constexpr double j = 0.200;
338  //
339  const double beta2 = (e2-mass2)/e2;
340  const double gamma2 = 1./(1.0 - beta2);
341  const double epsilon = 0.5*kappa*x*matConst/beta2;
342  //
343  return 0.001*epsilon*( log(2.*me*beta2*gamma2*epsilon*Iinv2) + j - beta2 );
344 }
345 //
346 double TrajectoryMCSFitter::energyLossBetheBloch(const double mass,const double e2) const {
347  // stolen, mostly, from GFMaterialEffects.
348  constexpr double Iinv = 1./188.E-6;
349  constexpr double matConst = 1.4*18./40.;//density*Z/A
350  constexpr double me = 0.511;
351  constexpr double kappa = 0.307075;
352  //
353  const double beta2 = (e2-mass*mass)/e2;
354  const double gamma2 = 1./(1.0 - beta2);
355  const double massRatio = me/mass;
356  const double argument = (2.*me*gamma2*beta2*Iinv) * std::sqrt(1+2*std::sqrt(gamma2)*massRatio + massRatio*massRatio);
357  //
358  double dedx = kappa*matConst/beta2;
359  //
360  if (mass==0.0) return(0.0);
361  if (argument <= exp(beta2)) {
362  dedx = 0.;
363  } else{
364  dedx *= (log(argument)-beta2)*1.E-3; // Bethe-Bloch, converted to GeV/cm
365  if (dedx<0.) dedx = 0.;
366  }
367  return dedx;
368 }
369 //
370 double TrajectoryMCSFitter::GetE(const double initial_E, const double length_travelled, const double m) const {
371  //
372  if (eLossMode_==1) {
373  // ELoss mode: MIP (constant)
374  constexpr double kcal = 0.002105;
375  return (initial_E - kcal*length_travelled);//energy at this segment
376  }
377  //
378  // Non constant energy loss distribution
379  const double step_size = length_travelled / nElossSteps_;
380  //
381  double current_E = initial_E;
382  const double m2 = m*m;
383  //
384  for (auto i = 0; i < nElossSteps_; ++i) {
385  if (eLossMode_==2) {
386  double dedx = energyLossBetheBloch(m,current_E);
387  current_E -= (dedx * step_size);
388  } else {
389  // MPV of Landau energy loss distribution
390  current_E -= energyLossLandau(m2,current_E*current_E,step_size);
391  }
392  if ( current_E <= m ) {
393  // std::cout<<"WARNING: current_E less than mu mass. it is "<<current_E<<std::endl;
394  return 0.;
395  }
396  }
397  return current_E;
398 }
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
Utilities to extend the interface of geometry vectors.
double beta(double KE, const simb::MCParticle *part)
static QCString result
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:40
T DirectionAtPoint(unsigned int p) const
Direction at point p. Use e.g. as:
recob::tracking::Vector_t Vector_t
Definition: TrackState.h:19
Helper class to compute the middle point in a point set.
size_t LastValidPoint() const
Returns the index of the last valid point in the trajectory.
bool isValid
Whether this ID points to a valid element.
Definition: geo_types.h:211
STL namespace.
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:164
void add(Point const &p)
Accumulates a point.
art framework interface to geometry description
geo::TPCID FindTPCAtPosition(double const worldLoc[3]) const
Returns the ID of the TPC at specified location.
T abs(T value)
double Length(size_t startAt=0) const
Returns the approximate length of the trajectory.
static constexpr double m2
Definition: Units.h:72
A trajectory in space reconstructed from hits.
T LocationAtPoint(unsigned int p) const
Position at point p. Use e.g. as:
p
Definition: test.py:223
static int max(int a, int b)
The data type to uniquely identify a TPC.
Definition: geo_types.h:386
auto norm(Vector const &v)
Return norm of the specified vector.
Class storing the result of the Maximum Likelihood fit of Multiple Coulomb Scattering angles between ...
Definition: MCSFitResult.h:19
#define M_PI
Definition: includeROOT.h:54
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
size_t FirstValidPoint() const
Returns the index of the first valid point in the trajectory.
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
static constexpr size_t InvalidIndex
Value returned on failed index queries.
size_t NextValidPoint(size_t index) const
Returns the index of the next valid point in the trajectory.
list x
Definition: train.py:276
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
Vector_t StartDirection() const
Returns the direction of the trajectory at the first point.
recob::tracking::Vector_t Vector_t
QTextStream & endl(QTextStream &s)