MVAAlg.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////
2 // \fileMVAAlg.cxx
3 // m.haigh@warwick.ac.uk
4 ////////////////////////////////////////////////////////////////////
5 
8 #include "larcorealg/CoreUtils/quiet_Math_Functor.h" // remove the wrapper when ROOT header is fixed
13 
15 #include "canvas/Persistency/Common/FindManyP.h"
16 #include "canvas/Persistency/Common/FindOneP.h"
17 
18 #include "TPrincipal.h"
19 #include <Fit/Fitter.h>
20 
21 #include <cmath>
22 
24  : fCaloAlg(pset.get<fhicl::ParameterSet>("CalorimetryAlg")), fReader("")
25 {
26  fHitLabel = pset.get<std::string>("HitLabel");
27  fTrackLabel = pset.get<std::string>("TrackLabel");
28  fShowerLabel = pset.get<std::string>("ShowerLabel");
29  fSpacePointLabel = pset.get<std::string>("SpacePointLabel");
30  fTrackingLabel = pset.get<std::string>("TrackingLabel", "");
31 
32  fCheatVertex = pset.get<bool>("CheatVertex", false);
33 
34  fReader.AddVariable("evalRatio", &fResHolder.evalRatio);
35  fReader.AddVariable("coreHaloRatio", &fResHolder.coreHaloRatio);
36  fReader.AddVariable("concentration", &fResHolder.concentration);
37  fReader.AddVariable("conicalness", &fResHolder.conicalness);
38  fReader.AddVariable("dEdxStart", &fResHolder.dEdxStart);
39  fReader.AddVariable("dEdxEnd", &fResHolder.dEdxEnd);
40  fReader.AddVariable("dEdxEndRatio", &fResHolder.dEdxEndRatio);
41 
42  fMVAMethods = pset.get<std::vector<std::string>>("MVAMethods");
43  std::vector<std::string> weightFileBnames = pset.get<std::vector<std::string>>("WeightFiles");
44 
45  cet::search_path searchPath("FW_SEARCH_PATH");
46  for (auto fileIter = weightFileBnames.begin(); fileIter != weightFileBnames.end(); ++fileIter) {
47  std::string fileWithPath;
48  if (!searchPath.find_file(*fileIter, fileWithPath)) {
49  fWeightFiles.clear();
50  fMVAMethods.clear();
51  throw cet::exception("MVAPID") << "Unable to find weight file " << *fileIter
52  << " in search path " << searchPath.to_string() << std::endl;
53  }
54  fWeightFiles.push_back(fileWithPath);
55  }
56 
57  if (fMVAMethods.size() != fWeightFiles.size()) {
58  std::cerr << "Mismatch in number of MVA methods and weight files!" << std::endl;
59  exit(1);
60  }
61 
62  for (unsigned int iMethod = 0; iMethod != fMVAMethods.size(); ++iMethod) {
63  fReader.BookMVA(fMVAMethods[iMethod], fWeightFiles[iMethod]);
64  }
65 }
66 
67 int
69 {
70  const double fiducialDist = 5.0;
71 
72  if (pos.X() > (fDetMinX + fiducialDist) && pos.X() < (fDetMaxX - fiducialDist) &&
73  pos.Y() > (fDetMinY + fiducialDist) && pos.Y() < (fDetMaxY - fiducialDist) &&
74  pos.Z() > (fDetMinZ + fiducialDist) && pos.Z() < (fDetMaxZ - fiducialDist))
75  return 1;
76  else
77  return 0;
78 }
79 
80 void
82 {
83 
85 
86  fDetMinX = 999999.9;
87  fDetMaxX = -999999.9;
88  fDetMinY = 999999.9;
89  fDetMaxY = -999999.9;
90  fDetMinZ = 999999.9;
91  fDetMaxZ = -999999.9;
92 
93  for (unsigned int t = 0; t < geom->TotalNTPC(); t++) {
94  if (geom->TPC(t).MinX() < fDetMinX) fDetMinX = geom->TPC(t).MinX();
95  if (geom->TPC(t).MaxX() > fDetMaxX) fDetMaxX = geom->TPC(t).MaxX();
96  if (geom->TPC(t).MinY() < fDetMinY) fDetMinY = geom->TPC(t).MinY();
97  if (geom->TPC(t).MaxY() > fDetMaxY) fDetMaxY = geom->TPC(t).MaxY();
98  if (geom->TPC(t).MinZ() < fDetMinZ) fDetMinZ = geom->TPC(t).MinZ();
99  if (geom->TPC(t).MaxZ() > fDetMaxZ) fDetMaxZ = geom->TPC(t).MaxZ();
100  }
101 }
102 
103 void
105 {
106 
108 
109  fNormToWiresY.clear();
110  fNormToWiresZ.clear();
111 
112  int planeKey;
113 
114  //Get normals to wires for each plane in the detector
115  //This assumes equal numbers of TPCs in each cryostat and equal numbers of planes in each TPC
116  for (geo::PlaneGeo const& plane : geom->IteratePlanes()) {
117  std::string id = std::string(plane.ID());
118  int pcryo = id.find("C");
119  int ptpc = id.find("T");
120  int pplane = id.find("P");
121  std::string scryo = id.substr(pcryo + 2, 2);
122  std::string stpc = id.substr(ptpc + 2, 2);
123  std::string splane = id.substr(pplane + 2, 2);
124  int icryo = std::stoi(scryo);
125  int itpc = std::stoi(stpc);
126  int iplane = std::stoi(splane);
127  planeKey = icryo * geom->NTPC(0) * geom->Nplanes(0, 0) + itpc * geom->Nplanes(0, 0) +
128  iplane; //single index for all planes in detector
129  fNormToWiresY.insert(
130  std::make_pair(planeKey, -plane.Wire(0).Direction().Z())); //y component of normal
131  fNormToWiresZ.insert(
132  std::make_pair(planeKey, plane.Wire(0).Direction().Y())); //z component of normal
133  }
134 }
135 
136 void
138  std::vector<anab::MVAPIDResult>& result,
141 {
142  auto const clockData =
144  auto const detProp =
146  this->PrepareEvent(evt, clockData);
147 
148  for (auto trackIter = fTracks.begin(); trackIter != fTracks.end(); ++trackIter) {
149  mvapid::MVAAlg::SortedObj sortedObj;
150 
151  std::vector<double> eVals, eVecs;
152  int isStoppingReco;
153  this->RunPCA(fTracksToHits[*trackIter], eVals, eVecs);
154  double evalRatio;
155  if (eVals[0] < 0.0001)
156  evalRatio = 0.0;
157  else
158  evalRatio = std::sqrt(eVals[1] * eVals[1] + eVals[2] * eVals[2]) / eVals[0];
159  this->FitAndSortTrack(*trackIter, isStoppingReco, sortedObj);
160  double coreHaloRatio, concentration, conicalness;
161  this->_Var_Shape(sortedObj, coreHaloRatio, concentration, conicalness);
162  double dEdxStart = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0., 0.05);
163  double dEdxEnd = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.9, 1.0);
164  double dEdxPenultimate = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.8, 0.9);
165 
166  fResHolder.isTrack = 1;
167  fResHolder.isStoppingReco = isStoppingReco;
168  fResHolder.nSpacePoints = sortedObj.hitMap.size();
169  fResHolder.trackID = (*trackIter)->ID();
170  fResHolder.evalRatio = evalRatio;
171  fResHolder.concentration = concentration;
172  fResHolder.coreHaloRatio = coreHaloRatio;
173  fResHolder.conicalness = conicalness;
174  fResHolder.dEdxStart = dEdxStart;
175  fResHolder.dEdxEnd = dEdxEnd;
176  if (dEdxPenultimate < 0.1)
177  fResHolder.dEdxEndRatio = 1.0;
178  else
179  fResHolder.dEdxEndRatio = dEdxEnd / dEdxPenultimate;
180  fResHolder.length = sortedObj.length;
181 
182  for (auto methodIter = fMVAMethods.begin(); methodIter != fMVAMethods.end(); ++methodIter) {
183  fResHolder.mvaOutput[*methodIter] = fReader.EvaluateMVA(*methodIter);
184  }
185  result.push_back(fResHolder);
186  util::CreateAssn(evt, result, *trackIter, trackAssns);
187  }
188 
189  for (auto showerIter = fShowers.begin(); showerIter != fShowers.end(); ++showerIter) {
190  mvapid::MVAAlg::SortedObj sortedObj;
191 
192  std::vector<double> eVals, eVecs;
193  int isStoppingReco;
194 
195  this->RunPCA(fShowersToHits[*showerIter], eVals, eVecs);
196 
197  double evalRatio;
198  if (eVals[0] < 0.0001)
199  evalRatio = 0.0;
200  else
201  evalRatio = std::sqrt(eVals[1] * eVals[1] + eVals[2] * eVals[2]) / eVals[0];
202 
203  this->SortShower(*showerIter, isStoppingReco, sortedObj);
204 
205  double coreHaloRatio, concentration, conicalness;
206  this->_Var_Shape(sortedObj, coreHaloRatio, concentration, conicalness);
207  double dEdxStart = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0., 0.05);
208  double dEdxEnd = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.9, 1.0);
209  double dEdxPenultimate = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.8, 0.9);
210 
211  fResHolder.isTrack = 0;
212  fResHolder.isStoppingReco = isStoppingReco;
213  fResHolder.nSpacePoints = sortedObj.hitMap.size();
215  (*showerIter)->ID() + 1000; //For the moment label showers by adding 1000 to ID
216 
217  fResHolder.evalRatio = evalRatio;
218  fResHolder.concentration = concentration;
219  fResHolder.coreHaloRatio = coreHaloRatio;
220  fResHolder.conicalness = conicalness;
221  fResHolder.dEdxStart = dEdxStart;
222  fResHolder.dEdxEnd = dEdxEnd;
223  if (dEdxPenultimate < 0.1)
224  fResHolder.dEdxEndRatio = 1.0;
225  else
226  fResHolder.dEdxEndRatio = dEdxEnd / dEdxPenultimate;
227  fResHolder.length = sortedObj.length;
228 
229  for (auto methodIter = fMVAMethods.begin(); methodIter != fMVAMethods.end(); ++methodIter) {
230  fResHolder.mvaOutput[*methodIter] = fReader.EvaluateMVA(*methodIter);
231  }
232  result.push_back(fResHolder);
233  util::CreateAssn(evt, result, *showerIter, showerAssns);
234  }
235 }
236 
237 void
239 {
240 
241  fHits.clear();
242  fSpacePoints.clear();
243  fTracks.clear();
244  fShowers.clear();
245  fSpacePointsToHits.clear();
246  fHitsToSpacePoints.clear();
247  fTracksToHits.clear();
248  fTracksToSpacePoints.clear();
249  fShowersToHits.clear();
250  fShowersToSpacePoints.clear();
251 
252  fEventT0 = trigger_offset(clockData);
253 
255  evt.getByLabel(fHitLabel, hitsHandle);
256 
257  for (unsigned int iHit = 0; iHit < hitsHandle->size(); ++iHit) {
258  const art::Ptr<recob::Hit> hit(hitsHandle, iHit);
259  fHits.push_back(hit);
260  }
261 
263  evt.getByLabel(fTrackLabel, tracksHandle);
264 
265  for (unsigned int iTrack = 0; iTrack < tracksHandle->size(); ++iTrack) {
266  const art::Ptr<recob::Track> track(tracksHandle, iTrack);
267  fTracks.push_back(track);
268  }
269 
271  evt.getByLabel(fShowerLabel, showersHandle);
272 
273  for (unsigned int iShower = 0; iShower < showersHandle->size(); ++iShower) {
274  const art::Ptr<recob::Shower> shower(showersHandle, iShower);
275  fShowers.push_back(shower);
276  }
277 
279  evt.getByLabel(fSpacePointLabel, spHandle);
280 
281  for (unsigned int iSP = 0; iSP < spHandle->size(); ++iSP) {
282  const art::Ptr<recob::SpacePoint> spacePoint(spHandle, iSP);
283  fSpacePoints.push_back(spacePoint);
284  }
285 
286  art::FindManyP<recob::Hit> findTracksToHits(fTracks, evt, fTrackLabel);
287  art::FindManyP<recob::Hit> findShowersToHits(fShowers, evt, fShowerLabel);
288  art::FindOneP<recob::Hit> findSPToHits(fSpacePoints, evt, fSpacePointLabel);
289 
290  for (unsigned int iSP = 0; iSP < fSpacePoints.size(); ++iSP) {
291  const art::Ptr<recob::SpacePoint> spacePoint = fSpacePoints.at(iSP);
292 
293  const art::Ptr<recob::Hit> hit = findSPToHits.at(iSP);
294  fSpacePointsToHits[spacePoint] = hit;
295  fHitsToSpacePoints[hit] = spacePoint;
296  }
297 
298  for (unsigned int iTrack = 0; iTrack < fTracks.size(); ++iTrack) {
299  const art::Ptr<recob::Track> track = fTracks.at(iTrack);
300 
301  const std::vector<art::Ptr<recob::Hit>> trackHits = findTracksToHits.at(iTrack);
302 
303  for (unsigned int iHit = 0; iHit < trackHits.size(); ++iHit) {
304  const art::Ptr<recob::Hit> hit = trackHits.at(iHit);
305  fTracksToHits[track].push_back(hit);
306  if (fHitsToSpacePoints.count(hit)) {
307  fTracksToSpacePoints[track].push_back(fHitsToSpacePoints.at(hit));
308  }
309  }
310  }
311 
312  for (unsigned int iShower = 0; iShower < fShowers.size(); ++iShower) {
313  const art::Ptr<recob::Shower> shower = fShowers.at(iShower);
314  const std::vector<art::Ptr<recob::Hit>> showerHits = findShowersToHits.at(iShower);
315 
316  for (unsigned int iHit = 0; iHit < showerHits.size(); ++iHit) {
317  const art::Ptr<recob::Hit> hit = showerHits.at(iHit);
318  fShowersToHits[shower].push_back(hit);
319  if (fHitsToSpacePoints.count(hit)) {
320  fShowersToSpacePoints[shower].push_back(fHitsToSpacePoints.at(hit));
321  }
322  }
323  }
324 
325  if (fCheatVertex) {
327  evt.getByLabel(fTrackingLabel, partHandle);
328 
329  if (partHandle->size() == 0 || partHandle->at(0).TrackId() != 1) {
330  std::cout << "Error, ID of first track in largeant list is not 0" << std::endl;
331  exit(1);
332  }
333  fVertex4Vect = partHandle->at(0).Position();
334  }
335 }
336 
337 void
339  int& isStoppingReco,
340  mvapid::MVAAlg::SortedObj& sortedTrack)
341 {
342 
343  sortedTrack.hitMap.clear();
344  TVector3 trackPoint, trackDir;
345  this->LinFit(track, trackPoint, trackDir);
346 
347  TVector3 nearestPointStart, nearestPointEnd;
348 
349  //For single-particle events can opt to cheat vertex from start of primary trajectory.
350  //Ok since in real events it should be possible to identify the true vertex.
351  if (fCheatVertex) {
352  if ((track->End<TVector3>() - fVertex4Vect.Vect()).Mag() >
353  (track->Vertex<TVector3>() - fVertex4Vect.Vect()).Mag()) {
354  nearestPointStart =
355  trackPoint +
356  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
357  nearestPointEnd = trackPoint + trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) /
358  trackDir.Mag2());
359  isStoppingReco = this->IsInActiveVol(track->End<TVector3>());
360  }
361  else {
362  nearestPointStart =
363  trackPoint +
364  trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) / trackDir.Mag2());
365  nearestPointEnd =
366  trackPoint +
367  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
368  isStoppingReco = this->IsInActiveVol(track->Vertex<TVector3>());
369  trackDir *= -1.;
370  }
371  }
372  else {
373  if (track->End<TVector3>().Z() >=
374  track->Vertex<TVector3>().Z()) { //Otherwise assume particle is forward-going for now...
375  nearestPointStart =
376  trackPoint +
377  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
378  nearestPointEnd = trackPoint + trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) /
379  trackDir.Mag2());
380  isStoppingReco = this->IsInActiveVol(track->End<TVector3>());
381  }
382  else {
383  nearestPointStart =
384  trackPoint +
385  trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) / trackDir.Mag2());
386  nearestPointEnd =
387  trackPoint +
388  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
389  isStoppingReco = this->IsInActiveVol(track->Vertex<TVector3>());
390  }
391 
392  if (trackDir.Z() <= 0) {
393  trackDir.SetX(-trackDir.X());
394  trackDir.SetY(-trackDir.Y());
395  trackDir.SetZ(-trackDir.Z());
396  }
397  }
398 
399  sortedTrack.start = nearestPointStart;
400  sortedTrack.end = nearestPointEnd;
401  sortedTrack.dir = trackDir;
402  sortedTrack.length = (nearestPointEnd - nearestPointStart).Mag();
403 
404  std::vector<art::Ptr<recob::Hit>> hits = fTracksToHits[track];
405 
406  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
407 
408  if (!fHitsToSpacePoints.count(*hitIter)) continue;
410 
411  TVector3 nearestPoint =
412  trackPoint + trackDir * (trackDir.Dot(TVector3(sp->XYZ()) - trackPoint) / trackDir.Mag2());
413  double lengthAlongTrack = (nearestPointStart - nearestPoint).Mag();
414  sortedTrack.hitMap.insert(std::pair<double, art::Ptr<recob::Hit>>(lengthAlongTrack, *hitIter));
415  }
416 }
417 
418 //void mvapid::MVAAlg::SortShower(art::Ptr<recob::Shower> shower,TVector3 dir,int& isStoppingReco,
419 // mvapid::MVAAlg::SortedObj& sortedShower){
420 void
422  int& isStoppingReco,
423  mvapid::MVAAlg::SortedObj& sortedShower)
424 {
425  sortedShower.hitMap.clear();
426 
427  std::vector<art::Ptr<recob::Hit>> hits = fShowersToHits[shower];
428 
429  TVector3 showerEnd(0, 0, 0);
430  double furthestHitFromStart = -999.9;
431  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
432 
433  if (!fHitsToSpacePoints.count(*hitIter)) continue;
435  if ((TVector3(sp->XYZ()) - shower->ShowerStart()).Mag() > furthestHitFromStart) {
436  showerEnd = TVector3(sp->XYZ());
437  furthestHitFromStart = (TVector3(sp->XYZ()) - shower->ShowerStart()).Mag();
438  }
439  }
440 
441  TVector3 showerPoint, showerDir;
442  this->LinFitShower(shower, showerPoint, showerDir);
443 
444  TVector3 nearestPointStart, nearestPointEnd;
445 
446  //Ensure that shower is fitted in correct direction (assuming for now that particle moves in +z direction)
447 
448  if (fCheatVertex) {
449  if ((showerEnd - fVertex4Vect.Vect()).Mag() >
450  (shower->ShowerStart() - fVertex4Vect.Vect()).Mag()) {
451  nearestPointStart =
452  showerPoint +
453  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
454  nearestPointEnd =
455  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
456  isStoppingReco = this->IsInActiveVol(showerEnd);
457  }
458  else {
459  nearestPointStart =
460  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
461  nearestPointEnd =
462  showerPoint +
463  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
464  isStoppingReco = this->IsInActiveVol(shower->ShowerStart());
465  showerDir *= -1.;
466  }
467  }
468  else {
469  if (showerEnd.Z() >= shower->ShowerStart().Z()) {
470  nearestPointStart =
471  showerPoint +
472  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
473  nearestPointEnd =
474  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
475  isStoppingReco = this->IsInActiveVol(showerEnd);
476  }
477  else {
478  nearestPointStart =
479  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
480  nearestPointEnd =
481  showerPoint +
482  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
483  isStoppingReco = this->IsInActiveVol(shower->ShowerStart());
484  }
485 
486  if (showerDir.Z() <= 0) {
487  showerDir.SetX(-showerDir.X());
488  showerDir.SetY(-showerDir.Y());
489  showerDir.SetZ(-showerDir.Z());
490  }
491  }
492 
493  sortedShower.start = nearestPointStart;
494  sortedShower.end = nearestPointEnd;
495  sortedShower.dir = showerDir;
496  sortedShower.length = (nearestPointEnd - nearestPointStart).Mag();
497 
498  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
499 
500  if (!fHitsToSpacePoints.count(*hitIter)) continue;
502 
503  TVector3 nearestPoint =
504  showerPoint +
505  showerDir * (showerDir.Dot(TVector3(sp->XYZ()) - showerPoint) / showerDir.Mag2());
506  double lengthAlongShower = (nearestPointStart - nearestPoint).Mag();
507  sortedShower.hitMap.insert(
508  std::pair<double, art::Ptr<recob::Hit>>(lengthAlongShower, *hitIter));
509  }
510 }
511 void
513  std::vector<double>& eVals,
514  std::vector<double>& eVecs)
515 {
516  TPrincipal* principal = new TPrincipal(3, "D");
517 
518  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
519 
520  if (fHitsToSpacePoints.count(*hitIter)) {
521  principal->AddRow(fHitsToSpacePoints.at(*hitIter)->XYZ());
522  }
523  }
524 
525  // PERFORM PCA
526  principal->MakePrincipals();
527  // GET EIGENVALUES AND EIGENVECTORS
528  for (unsigned int i = 0; i < 3; ++i) {
529  eVals.push_back(principal->GetEigenValues()->GetMatrixArray()[i]);
530  }
531 
532  for (unsigned int i = 0; i < 9; ++i) {
533  eVecs.push_back(principal->GetEigenVectors()->GetMatrixArray()[i]);
534  }
535 }
536 void
538  double& coreHaloRatio,
539  double& concentration,
540  double& conicalness)
541 {
542 
543  static const unsigned int conMinHits = 3;
544  static const double minCharge = 0.1;
545  static const double conFracRange = 0.2;
546  static const double MoliereRadius = 10.1;
547  static const double MoliereRadiusFraction = 0.2;
548 
549  double totalCharge = 0;
550  double totalChargeStart = 0;
551  double totalChargeEnd = 0;
552 
553  double chargeCore = 0;
554  double chargeHalo = 0;
555  double chargeCon = 0;
556  unsigned int nHits = 0;
557 
558  //stuff for conicalness
559  double chargeConStart = 0;
560  double chargeConEnd = 0;
561  unsigned int nHitsConStart = 0;
562  unsigned int nHitsConEnd = 0;
563 
564  for (auto hitIter = track.hitMap.begin(); hitIter != track.hitMap.end(); ++hitIter) {
565  if (fHitsToSpacePoints.count(hitIter->second)) {
566  art::Ptr<recob::SpacePoint> sp = fHitsToSpacePoints.at(hitIter->second);
567 
568  double distFromTrackFit = ((TVector3(sp->XYZ()) - track.start).Cross(track.dir)).Mag();
569 
570  ++nHits;
571 
572  if (distFromTrackFit < MoliereRadiusFraction * MoliereRadius)
573  chargeCore += hitIter->second->Integral();
574  else
575  chargeHalo += hitIter->second->Integral();
576 
577  totalCharge += hitIter->second->Integral();
578 
579  chargeCon += hitIter->second->Integral() / std::max(1.E-2, distFromTrackFit);
580  if (hitIter->first / track.length < conFracRange) {
581  chargeConStart += distFromTrackFit * distFromTrackFit * hitIter->second->Integral();
582  ++nHitsConStart;
583  totalChargeStart += hitIter->second->Integral();
584  }
585  else if (1. - hitIter->first / track.length < conFracRange) {
586  chargeConEnd += distFromTrackFit * distFromTrackFit * hitIter->second->Integral();
587  ++nHitsConEnd;
588  totalChargeEnd += hitIter->second->Integral();
589  }
590  }
591  }
592 
593  coreHaloRatio = chargeHalo / TMath::Max(1.0E-3, chargeCore);
594  coreHaloRatio = TMath::Min(100.0, coreHaloRatio);
595  concentration = chargeCon / totalCharge;
596  if (nHitsConStart >= conMinHits && nHitsConEnd >= conMinHits && totalChargeEnd > minCharge &&
597  sqrt(chargeConStart) > minCharge && totalChargeStart > minCharge) {
598  conicalness = (sqrt(chargeConEnd) / totalChargeEnd) / (sqrt(chargeConStart) / totalChargeStart);
599  }
600  else {
601  conicalness = 1.;
602  }
603 }
604 
605 double
607  const detinfo::DetectorPropertiesData& det_prop,
608  const mvapid::MVAAlg::SortedObj& track,
609  double start,
610  double end)
611 {
612 
613  double trackLength = (track.end - track.start).Mag();
614  return CalcSegmentdEdxDist(clock_data, det_prop, track, start * trackLength, end * trackLength);
615 }
616 
617 double
619  const detinfo::DetectorPropertiesData& det_prop,
620  const mvapid::MVAAlg::SortedObj& track,
621  double distAtEnd)
622 {
623 
624  double trackLength = (track.end - track.start).Mag();
625  return CalcSegmentdEdxDist(clock_data, det_prop, track, trackLength - distAtEnd, trackLength);
626 }
627 
628 double
630  const detinfo::DetectorPropertiesData& det_prop,
631  const mvapid::MVAAlg::SortedObj& track,
632  double start,
633  double end)
634 {
636 
637  double totaldEdx = 0;
638  unsigned int nHits = 0;
639 
640  //Loop over hits again to calculate average dE/dx and shape variables
641  for (auto hitIter = track.hitMap.begin(); hitIter != track.hitMap.end(); ++hitIter) {
642 
643  if (hitIter->first < start) continue;
644  if (hitIter->first >= end) break;
645 
646  art::Ptr<recob::Hit> hit = hitIter->second;
647 
648  //Pitch to use in dEdx calculation
649  double yzPitch =
650  geom->WirePitch(hit->WireID().Plane,
651  hit->WireID().TPC); //pitch not taking into account angle of track or shower
652  double xComponent, pitch3D;
653 
654  TVector3 dir = track.dir;
655 
656  //This assumes equal numbers of TPCs in each cryostat and equal numbers of planes in each TPC
657  int planeKey = hit->WireID().Cryostat * geom->NTPC(0) * geom->Nplanes(0, 0) +
658  hit->WireID().TPC * geom->Nplanes(0, 0) + hit->WireID().Plane;
659 
660  if (fNormToWiresY.count(planeKey) && fNormToWiresZ.count(planeKey)) {
661  TVector3 normToWires(0.0, fNormToWiresY.at(planeKey), fNormToWiresZ.at(planeKey));
662  yzPitch =
663  geom->WirePitch(hit->WireID().Plane, hit->WireID().TPC) / fabs(dir.Dot(normToWires));
664  }
665 
666  xComponent = yzPitch * dir[0] / sqrt(dir[1] * dir[1] + dir[2] * dir[2]);
667  pitch3D = sqrt(xComponent * xComponent + yzPitch * yzPitch);
668 
669  double dEdx = fCaloAlg.dEdx_AREA(clock_data, det_prop, *hit, pitch3D, fEventT0);
670  if (dEdx < 50.) {
671  ++nHits;
672  totaldEdx += dEdx;
673  }
674  }
675 
676  return nHits ? totaldEdx / nHits : 0;
677 }
678 
679 int
680 mvapid::MVAAlg::LinFit(const art::Ptr<recob::Track> track, TVector3& trackPoint, TVector3& trackDir)
681 {
682 
683  const std::vector<art::Ptr<recob::SpacePoint>>& sp = fTracksToSpacePoints.at(track);
684 
685  TGraph2D grFit(1);
686  unsigned int iPt = 0;
687  for (auto spIter = sp.begin(); spIter != sp.end(); ++spIter) {
688  TVector3 point = (*spIter)->XYZ();
689  grFit.SetPoint(iPt++, point.X(), point.Y(), point.Z());
690  }
691 
692  //Lift from the ROOT line3Dfit.C tutorial
693  ROOT::Fit::Fitter fitter;
694  // make the functor object
695  mvapid::MVAAlg::SumDistance2 sdist(&grFit);
696 
697  ROOT::Math::Functor fcn(sdist, 6);
698 
699  //Initial fit parameters from track start and end...
700  TVector3 trackStart = track->Vertex<TVector3>();
701  TVector3 trackEnd = track->End<TVector3>();
702  trackDir = (trackEnd - trackStart).Unit();
703 
704  TVector3 x0 = trackStart - trackDir;
705  TVector3 u = trackDir;
706 
707  double pStart[6] = {x0.X(), u.X(), x0.Y(), u.Y(), x0.Z(), u.Z()};
708 
709  fitter.SetFCN(fcn, pStart);
710 
711  bool ok = fitter.FitFCN();
712  if (!ok) {
713  trackPoint.SetXYZ(x0.X(), x0.Y(), x0.Z());
714  trackDir.SetXYZ(u.X(), u.Y(), u.Z());
715  trackDir = trackDir.Unit();
716  return 1;
717  }
718  else {
719  const ROOT::Fit::FitResult& result = fitter.Result();
720  const double* parFit = result.GetParams();
721  trackPoint.SetXYZ(parFit[0], parFit[2], parFit[4]);
722  trackDir.SetXYZ(parFit[1], parFit[3], parFit[5]);
723  trackDir = trackDir.Unit();
724  return 0;
725  }
726 }
727 
728 int
730  TVector3& showerPoint,
731  TVector3& showerDir)
732 {
733 
734  const std::vector<art::Ptr<recob::SpacePoint>>& sp = fShowersToSpacePoints.at(shower);
735 
736  TGraph2D grFit(1);
737  unsigned int iPt = 0;
738  for (auto spIter = sp.begin(); spIter != sp.end(); ++spIter) {
739  TVector3 point = (*spIter)->XYZ();
740  grFit.SetPoint(iPt++, point.X(), point.Y(), point.Z());
741  }
742 
743  //Lift from the ROOT line3Dfit.C tutorial
744  ROOT::Fit::Fitter fitter;
745  // make the functor object
746  mvapid::MVAAlg::SumDistance2 sdist(&grFit);
747 
748  ROOT::Math::Functor fcn(sdist, 6);
749 
750  //Initial fit parameters from shower start and end...
751  TVector3 showerStart = shower->ShowerStart();
752  showerDir = shower->Direction().Unit();
753 
754  TVector3 x0 = showerStart - showerDir;
755  TVector3 u = showerDir;
756 
757  double pStart[6] = {x0.X(), u.X(), x0.Y(), u.Y(), x0.Z(), u.Z()};
758 
759  fitter.SetFCN(fcn, pStart);
760 
761  bool ok = fitter.FitFCN();
762  if (!ok) {
763  showerPoint.SetXYZ(x0.X(), x0.Y(), x0.Z());
764  showerDir.SetXYZ(u.X(), u.Y(), u.Z());
765  showerDir = showerDir.Unit();
766  return 1;
767  }
768  else {
769  const ROOT::Fit::FitResult& result = fitter.Result();
770  const double* parFit = result.GetParams();
771  showerPoint.SetXYZ(parFit[0], parFit[2], parFit[4]);
772  showerDir.SetXYZ(parFit[1], parFit[3], parFit[5]);
773  showerDir = showerDir.Unit();
774  return 0;
775  }
776 }
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
std::vector< art::Ptr< recob::Shower > > fShowers
Definition: MVAAlg.h:147
const TVector3 & ShowerStart() const
Definition: Shower.h:192
anab::MVAPIDResult fResHolder
Definition: MVAAlg.h:159
double CalcSegmentdEdxDistAtEnd(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const mvapid::MVAAlg::SortedObj &track, double distAtEnd)
Definition: MVAAlg.cxx:618
const calo::CalorimetryAlg fCaloAlg
Definition: MVAAlg.h:131
unsigned int TotalNTPC() const
Returns the total number of TPCs in the detector.
int LinFit(const art::Ptr< recob::Track > track, TVector3 &trackPoint, TVector3 &trackDir)
Definition: MVAAlg.cxx:680
void RunPCA(std::vector< art::Ptr< recob::Hit >> &hits, std::vector< double > &eVals, std::vector< double > &eVecs)
Definition: MVAAlg.cxx:512
static QCString result
IteratorBox< plane_iterator,&GeometryCore::begin_plane,&GeometryCore::end_plane > IteratePlanes() const
Enables ranged-for loops on all planes of the detector.
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:156
double fDetMaxX
Definition: MVAAlg.h:135
std::vector< std::string > fMVAMethods
Definition: MVAAlg.h:163
std::string string
Definition: nybbler.cc:12
geo::WireID WireID() const
Definition: Hit.h:233
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:88
struct vector vector
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
int IsInActiveVol(const TVector3 &pos)
Definition: MVAAlg.cxx:68
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::Hit > > > fTracksToHits
Definition: MVAAlg.h:151
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:91
std::string fSpacePointLabel
Definition: MVAAlg.h:143
string dir
void PrepareEvent(const art::Event &event, const detinfo::DetectorClocksData &clockData)
Definition: MVAAlg.cxx:238
MVAAlg(fhicl::ParameterSet const &pset)
Definition: MVAAlg.cxx:23
TMVA::Reader fReader
Definition: MVAAlg.h:161
Particle class.
bool fCheatVertex
Definition: MVAAlg.h:166
void RunPID(art::Event &evt, std::vector< anab::MVAPIDResult > &result, art::Assns< recob::Track, anab::MVAPIDResult, void > &trackAssns, art::Assns< recob::Shower, anab::MVAPIDResult, void > &showerAssns)
Definition: MVAAlg.cxx:137
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
art framework interface to geometry description
std::vector< art::Ptr< recob::Hit > > fHits
Definition: MVAAlg.h:149
int LinFitShower(const art::Ptr< recob::Shower > shower, TVector3 &showerPoint, TVector3 &showerDir)
Definition: MVAAlg.cxx:729
void SortShower(art::Ptr< recob::Shower > shower, int &isStoppingReco, mvapid::MVAAlg::SortedObj &sortedShower)
Definition: MVAAlg.cxx:421
void GetDetectorEdges()
Definition: MVAAlg.cxx:81
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
double dEdx(float dqdx, float Efield)
Definition: doAna.cpp:21
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
double fDetMinY
Definition: MVAAlg.h:135
std::map< double, const art::Ptr< recob::Hit > > hitMap
Definition: MVAAlg.h:44
std::map< int, double > fNormToWiresY
Definition: MVAAlg.h:137
double fDetMinZ
Definition: MVAAlg.h:135
T get(std::string const &key) const
Definition: ParameterSet.h:271
void GetWireNormals()
Definition: MVAAlg.cxx:104
std::vector< art::Ptr< recob::SpacePoint > > fSpacePoints
Definition: MVAAlg.h:148
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:82
Point_t const & Vertex() const
Definition: Track.h:124
double MinZ() const
Returns the world z coordinate of the start of the box.
double fDetMinX
Definition: MVAAlg.h:135
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
TLorentzVector fVertex4Vect
Definition: MVAAlg.h:168
const TVector3 & Direction() const
Definition: Shower.h:189
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::Hit > > > fShowersToHits
Definition: MVAAlg.h:153
static int max(int a, int b)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
Detector simulation of raw signals on wires.
double MaxY() const
Returns the world y coordinate of the end of the box.
std::string fTrackingLabel
Definition: MVAAlg.h:144
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
double CalcSegmentdEdxFrac(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
Definition: MVAAlg.cxx:606
std::string fTrackLabel
Definition: MVAAlg.h:140
std::string fHitLabel
Definition: MVAAlg.h:142
Contains all timing reference information for the detector.
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::SpacePoint > > > fTracksToSpacePoints
Definition: MVAAlg.h:152
double MaxZ() const
Returns the world z coordinate of the end of the box.
double dEdx_AREA(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, recob::Hit const &hit, double pitch, double T0=0) const
void _Var_Shape(const SortedObj &track, double &coreHaloRatio, double &concentration, double &conicalness)
Definition: MVAAlg.cxx:537
std::vector< art::Ptr< recob::Track > > fTracks
Definition: MVAAlg.h:146
int trigger_offset(DetectorClocksData const &data)
const Double32_t * XYZ() const
Definition: SpacePoint.h:76
Point_t const & End() const
Definition: Track.h:125
std::string fShowerLabel
Definition: MVAAlg.h:141
double fDetMaxZ
Definition: MVAAlg.h:135
TPCGeo const & TPC(unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified TPC.
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::SpacePoint > > > fShowersToSpacePoints
Definition: MVAAlg.h:155
TCEvent evt
Definition: DataStructs.cxx:7
auto const & get(AssnsNode< L, R, D > const &r)
Definition: AssnsNode.h:115
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
std::map< std::string, double > mvaOutput
Definition: MVAPIDResult.h:27
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > fSpacePointsToHits
Definition: MVAAlg.h:157
double fDetMaxY
Definition: MVAAlg.h:135
void FitAndSortTrack(art::Ptr< recob::Track > track, int &isStoppingReco, SortedObj &sortedObj)
Definition: MVAAlg.cxx:338
unsigned int trackID
Definition: MVAPIDResult.h:21
std::map< int, double > fNormToWiresZ
Definition: MVAAlg.h:138
double MinY() const
Returns the world y coordinate of the start of the box.
double fEventT0
Definition: MVAAlg.h:133
double CalcSegmentdEdxDist(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
Definition: MVAAlg.cxx:629
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
std::vector< std::string > fWeightFiles
Definition: MVAAlg.h:164