All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ClusterParamsAlg.cxx
Go to the documentation of this file.
1 #include "ClusterParamsAlg.h"
2 
3 // LArSoft includes
4 #include "lardata/Utilities/SimpleFits.h" // LinearFit<>
5 #include "lardataalg/Utilities/StatCollector.h" // StatCollector<>
6 
7 //-----Math-------
8 #include <algorithm>
9 #include <cmath>
10 #include <iostream>
11 #include <vector>
12 
13 #include "TCanvas.h"
14 #include "TH1.h"
15 #include "TLegend.h"
16 #include "TMath.h"
17 #include "TPrincipal.h"
18 #include "TStopwatch.h"
19 
23 
24 #include "cetlib/pow.h"
25 
26 namespace {
27  constexpr double PI{3.14159265};
28 }
29 
30 namespace cluster {
31 
33  {
34  fMinNHits = 10;
35  enableFANN = false;
36  verbose = true;
37  Initialize();
38  }
39 
40  ClusterParamsAlg::ClusterParamsAlg(const std::vector<util::PxHit>& inhitlist)
41  {
42  fMinNHits = 10;
43  enableFANN = false;
44  verbose = true;
45  SetHits(inhitlist);
46  }
47 
48  int
49  ClusterParamsAlg::SetHits(const std::vector<util::PxHit>& inhitlist)
50  {
51  Initialize();
52 
53  // Make default values
54  // Is done by the struct
55  if (!(inhitlist.size())) {
56  throw CRUException("Provided empty hit list!");
57  return -1;
58  }
59 
60  fHitVector.reserve(inhitlist.size());
61  for (auto h : inhitlist)
62  fHitVector.push_back(h);
63 
64  fPlane = fHitVector[0].plane;
65 
66  if (fHitVector.size() < fMinNHits) {
67  if (verbose)
68  std::cout << " the hitlist is too small. Continuing to run may result in crash!!! "
69  << std::endl;
70  return -1;
71  }
72 
73  return fHitVector.size();
74  }
75 
76  void
78  {
79  fPlane = p;
80  for (auto& h : fHitVector)
81  h.plane = p;
84  }
85 
86  void
87  ClusterParamsAlg::GetFANNVector(std::vector<float>& data)
88  {
89  unsigned int length = 13;
90  if (data.size() != length) data.resize(length);
91  data[0] = fParams.opening_angle / PI;
93  data[2] = fParams.closing_angle / PI;
95  data[4] = fParams.eigenvalue_principal;
96  data[5] = fParams.eigenvalue_secondary;
97  data[6] = fParams.width / fParams.length;
100  data[9] = fParams.N_Hits_HC / fParams.N_Hits;
101  data[10] = fParams.modified_hit_density;
102  data[11] = fParams.RMS_charge / fParams.mean_charge;
103  data[12] = log(fParams.sum_charge / fParams.length);
104  return;
105  }
106 
107  void
109  {
110  std::vector<float> data;
111  GetFANNVector(data);
112  if (verbose) {
113  std::cout << "Printing FANN input vector:\n"
114  << " Opening Angle (normalized) ... : " << data[0] << "\n"
115  << " Opening Angle charge weight .. : " << data[1] << "\n"
116  << " Closing Angle (normalized) ... : " << data[2] << "\n"
117  << " Closing Angle charge weight .. : " << data[3] << "\n"
118  << " Principal Eigenvalue ......... : " << data[4] << "\n"
119  << " Secondary Eigenvalue ......... : " << data[5] << "\n"
120  << " Width / Length ............... : " << data[6] << "\n"
121  << " Hit Density / M.H.D. ......... : " << data[7] << "\n"
122  << " Percent MultiHit Wires ....... : " << data[8] << "\n"
123  << " Percent High Charge Hits ..... : " << data[9] << "\n"
124  << " Modified Hit Density ......... : " << data[10] << "\n"
125  << " Charge RMS / Mean Charge ...... : " << data[11] << "\n"
126  << " log(Sum Charge / Length) ...... : " << data[12] << "\n";
127  }
128  }
129 
130  void
132  {
133  fTimeRecord_ProcName.clear();
134  fTimeRecord_ProcTime.clear();
135 
136  TStopwatch localWatch;
137  localWatch.Start();
138 
139  //--- Initilize attributes values ---//
140  fFinishedGetAverages = false;
141  fFinishedGetRoughAxis = false;
142  fFinishedGetProfileInfo = false;
144  fFinishedRefineDirection = false;
145  fFinishedGetFinalSlope = false;
147  fFinishedTrackShowerSep = false;
148  fFinishedGetEndCharges = false;
149 
150  fRough2DSlope = -999.999; // slope
151  fRough2DIntercept = -999.999; // slope
152 
153  fRoughBeginPoint.w = -999.999;
154  fRoughEndPoint.w = -999.999;
155 
156  fRoughBeginPoint.t = -999.999;
157  fRoughEndPoint.t = -999.999;
158 
159  fProfileIntegralForward = 999.999;
160  fProfileIntegralBackward = 999.999;
161  fProfileMaximumBin = -999;
162 
163  fChargeCutoffThreshold.clear();
164  fChargeCutoffThreshold.reserve(3);
165  fChargeCutoffThreshold.push_back(500);
166  fChargeCutoffThreshold.push_back(500);
167  fChargeCutoffThreshold.push_back(1000);
168 
169  fHitVector.clear();
170 
171  fParams.Clear();
172 
173  fTimeRecord_ProcName.push_back("Initialize");
174  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
175  }
176 
177  void
179  {
180  enableFANN = true;
181  }
182 
183  void
185  bool override_DoGetAverages,
186  bool override_DoGetRoughAxis,
187  bool override_DoGetProfileInfo,
188  bool override_DoStartPointsAndDirection,
189  bool override_DoGetFinalSlope,
190  bool override_DoTrackShowerSep,
191  bool override_DoEndCharge)
192  {
193  GetAverages(override_DoGetAverages);
194  GetRoughAxis(override_DoGetRoughAxis);
195  GetProfileInfo(gser, override_DoGetProfileInfo);
196  RefineStartPointAndDirection(gser, override_DoStartPointsAndDirection);
197  GetFinalSlope(gser, override_DoGetFinalSlope);
198  GetEndCharges(gser, override_DoEndCharge);
199  TrackShowerSeparation(override_DoTrackShowerSep);
200  }
201 
202  void
204  {
205  if (!override) { //Override being set, we skip all this logic.
206  //OK, no override. Stop if we're already finshed.
207  if (fFinishedGetAverages) return;
208  }
209 
210  TStopwatch localWatch;
211  localWatch.Start();
212 
213  TPrincipal fPrincipal(2, "D");
214 
215  fParams.N_Hits = fHitVector.size();
216 
217  std::map<double, int> wireMap;
218 
219  lar::util::StatCollector<double> charge, sumADC;
220 
221  int uniquewires = 0;
222  int multi_hit_wires = 0;
223  for (auto& hit : fHitVector) {
224  double data[2];
225  data[0] = hit.w;
226  data[1] = hit.t;
227  fPrincipal.AddRow(data);
228  fParams.charge_wgt_x += hit.w * hit.charge;
229  fParams.charge_wgt_y += hit.t * hit.charge;
230  charge.add(hit.charge);
231  sumADC.add(hit.sumADC);
232 
233  if (wireMap[hit.w] == 0) { uniquewires++; }
234  if (wireMap[hit.w] == 1) { multi_hit_wires++; }
235  wireMap[hit.w]++;
236  }
237 
238  fParams.sum_charge = charge.Sum();
239  fParams.mean_charge = charge.Average();
240  fParams.rms_charge = charge.RMS();
241 
242  fParams.sum_ADC = sumADC.Sum();
243  fParams.mean_ADC = sumADC.Average();
244  fParams.rms_ADC = sumADC.RMS();
245 
246  fParams.N_Wires = uniquewires;
247  fParams.multi_hit_wires = multi_hit_wires;
248 
249  if (fPrincipal.GetMeanValues()->GetNrows() < 2) { throw cluster::CRUException(); }
250 
251  fParams.mean_x = (*fPrincipal.GetMeanValues())[0];
252  fParams.mean_y = (*fPrincipal.GetMeanValues())[1];
254 
255  if (fParams.sum_charge != 0.) {
258  }
259  else { // "SNAFU"; use the mean
262  }
263 
264  fPrincipal.MakePrincipals();
265 
266  fParams.eigenvalue_principal = (*fPrincipal.GetEigenValues())[0];
267  fParams.eigenvalue_secondary = (*fPrincipal.GetEigenValues())[1];
268 
269  fFinishedGetAverages = true;
270 
271  fTimeRecord_ProcName.push_back("GetAverages");
272  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
273  }
274 
275  // Also does the high hitlist
276  void
278  {
279  if (!override) { //Override being set, we skip all this logic.
280  //OK, no override. Stop if we're already finshed.
281  if (fFinishedGetRoughAxis) return;
282  //Try to run the previous function if not yet done.
283  if (!fFinishedGetAverages) GetAverages(true);
284  }
285  else {
286  //Try to run the previous function if not yet done.
287  if (!fFinishedGetAverages) GetAverages(true);
288  }
289 
290  TStopwatch localWatch;
291  localWatch.Start();
292 
293  double rmsx = 0.0;
294  double rmsy = 0.0;
295  double rmsq = 0.0;
296  //using the charge weighted coordinates find the axis from slope
297  double ncw = 0;
298  double sumtime = 0; //from sum averages
299  double sumwire = 0; //from sum averages
300  double sumwiretime = 0; //sum over (wire*time)
301  double sumwirewire = 0; //sum over (wire*wire)
302  //next loop over all hits again
303 
304  for (auto& hit : fHitVector) {
305  // First, abuse this loop to calculate rms in x and y
306  rmsx += pow(fParams.mean_x - hit.w, 2) / fParams.N_Hits;
307  rmsy += pow(fParams.mean_y - hit.t, 2) / fParams.N_Hits;
308  rmsq += pow(fParams.mean_charge - hit.charge, 2) / fParams.N_Hits;
309  //if charge is above avg_charge
310 
311  if (hit.charge > fParams.mean_charge) {
312  ncw += 1;
313  sumwire += hit.w;
314  sumtime += hit.t;
315  sumwiretime += hit.w * hit.t;
316  sumwirewire += pow(hit.w, 2);
317  } //for high charge
318  } //For hh loop
319 
320  fParams.rms_x = sqrt(rmsx);
321  fParams.rms_y = sqrt(rmsy);
322  fParams.RMS_charge = sqrt(rmsq);
323 
324  fParams.N_Hits_HC = ncw;
325  //Looking for the slope and intercept of the line above avg_charge hits
326 
327  if ((ncw * sumwirewire - sumwire * sumwire) > 0.00001)
328  fRough2DSlope = (ncw * sumwiretime - sumwire * sumtime) /
329  (ncw * sumwirewire - sumwire * sumwire); //slope for cw
330  else
331  fRough2DSlope = 999;
334  fParams.charge_wgt_y - fRough2DSlope * (fParams.charge_wgt_x) : //intercept for cw
335  0.;
336 
337  //Getthe 2D_angle
338  fParams.cluster_angle_2d = atan(fRough2DSlope) * 180 / PI;
339 
340  fFinishedGetRoughAxis = true;
341 
342  fTimeRecord_ProcName.push_back("GetRoughAxis");
343  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
344  }
345 
346  void
348  {
349  if (!override) { //Override being set, we skip all this logic.
350  //OK, no override. Stop if we're already finshed.
351  if (fFinishedGetProfileInfo) return;
352  //Try to run the previous function if not yet done.
354  }
355  else {
357  }
358 
359  TStopwatch localWatch;
360  localWatch.Start();
361 
362  bool drawOrtHistos = false;
363 
364  //these variables need to be initialized to other values?
365  if (fRough2DSlope == -999.999 || fRough2DIntercept == -999.999) GetRoughAxis(true);
366 
367  //get slope of lines orthogonal to those found crossing the shower.
368  double inv_2d_slope = 0;
369  if (fabs(fRough2DSlope) > 0.001)
370  inv_2d_slope = -1. / fRough2DSlope;
371  else
372  inv_2d_slope = -999999.;
373 
374  double InterHigh = -999999;
375  double InterLow = 999999;
376  double WireHigh = -999999;
377  double WireLow = 999999;
378  fInterHigh_side = -999999;
379  fInterLow_side = 999999;
380  double min_ortdist(999), max_ortdist(-999); // needed to calculate width
381  //loop over all hits. Create coarse and fine profiles
382  // of the charge weight to help refine the start/end
383  // points and have a first guess of direction
384 
385  for (auto const& hit : fHitVector) {
386 
387  //calculate intercepts along
388  double intercept = hit.t - inv_2d_slope * (double)(hit.w);
389  double side_intercept = hit.t - fRough2DSlope * (double)(hit.w);
390 
391  util::PxPoint OnlinePoint;
392  // get coordinates of point on axis.
393  gser.GetPointOnLine(fRough2DSlope, fRough2DIntercept, &hit, OnlinePoint);
394 
395  double ortdist = gser.Get2DDistance(&OnlinePoint, &hit);
396  if (ortdist < min_ortdist) min_ortdist = ortdist;
397  if (ortdist > max_ortdist) max_ortdist = ortdist;
398 
399  if (inv_2d_slope != -999999) //almost all cases
400  {
401  if (intercept > InterHigh) { InterHigh = intercept; }
402 
403  if (intercept < InterLow) { InterLow = intercept; }
404  }
405  else //slope is practically horizontal. Care only about wires.
406  {
407  if (hit.w > WireHigh) { WireHigh = hit.w; }
408  if (hit.w < WireLow) { WireLow = hit.w; }
409  }
410 
411  if (side_intercept > fInterHigh_side) { fInterHigh_side = side_intercept; }
412 
413  if (side_intercept < fInterLow_side) { fInterLow_side = side_intercept; }
414  }
415  // end of first HitIter loop, at this point we should
416  // have the extreme intercepts
417 
418  /////////////////////////////////////////////
419  // Second loop. Fill profiles.
420  /////////////////////////////////////////////
421 
422  // get projected points at the beginning and end of the axis.
423 
424  util::PxPoint HighOnlinePoint, LowOnlinePoint, BeginOnlinePoint, EndOnlinePoint;
425 
426  if (inv_2d_slope != -999999) //almost all cases
427  {
428  gser.GetPointOnLineWSlopes(fRough2DSlope, fRough2DIntercept, InterHigh, HighOnlinePoint);
429  gser.GetPointOnLineWSlopes(fRough2DSlope, fRough2DIntercept, InterLow, LowOnlinePoint);
430  }
431  else //need better treatment of horizontal events.
432  {
433  util::PxPoint ptemphigh(fPlane, WireHigh, (fInterHigh_side + fInterLow_side) / 2);
434  util::PxPoint ptemplow(fPlane, WireLow, (fInterHigh_side + fInterLow_side) / 2);
435  gser.GetPointOnLine(fRough2DSlope, fRough2DIntercept, &ptemphigh, HighOnlinePoint);
436  gser.GetPointOnLine(fRough2DSlope, fRough2DIntercept, &ptemplow, LowOnlinePoint);
437  }
438  if (verbose) {
439  std::cout << " Low online point: " << LowOnlinePoint.w << ", " << LowOnlinePoint.t
440  << " High: " << HighOnlinePoint.w << ", " << HighOnlinePoint.t << std::endl;
441  //define BeginOnlinePoint as the one with lower wire number (for now), adjust intercepts accordingly
442  }
443  if (HighOnlinePoint.w >= LowOnlinePoint.w) {
444  BeginOnlinePoint = LowOnlinePoint;
445  fBeginIntercept = InterLow;
446  EndOnlinePoint = HighOnlinePoint;
447  fEndIntercept = InterHigh;
448  }
449  else {
450  BeginOnlinePoint = HighOnlinePoint;
451  fBeginIntercept = InterHigh;
452  EndOnlinePoint = LowOnlinePoint;
453  fEndIntercept = InterLow;
454  }
455 
456  fProjectedLength = gser.Get2DDistance(&BeginOnlinePoint, &EndOnlinePoint);
457 
458  /////////////////// THe binning is now set here
459  fCoarseNbins = 2;
460 
462  if (fProfileNbins < 10) fProfileNbins = 10;
463 
464  fChargeProfile.clear();
465  fCoarseChargeProfile.clear();
466  fChargeProfile.resize(fProfileNbins, 0);
468 
469  ////////////////////////// end of new binning
470  // Some fitting variables to make a histogram:
471 
472  // TODO this is nonsense for small clusters
473  std::vector<double> ort_profile;
474  const int NBINS = 100;
475  ort_profile.resize(NBINS);
476 
477  std::vector<double> ort_dist_vect;
478  ort_dist_vect.reserve(fHitVector.size());
479 
480  double current_maximum = 0;
481  for (auto& hit : fHitVector) {
482 
483  util::PxPoint OnlinePoint;
484  // get coordinates of point on axis.
485  gser.GetPointOnLine(fRough2DSlope, &BeginOnlinePoint, &hit, OnlinePoint);
486  double ortdist = gser.Get2DDistance(&OnlinePoint, &hit);
487 
488  double linedist = gser.Get2DDistance(&OnlinePoint, &BeginOnlinePoint);
489  ort_dist_vect.push_back(ortdist);
490  int ortbin;
491  if (ortdist == 0)
492  ortbin = 0;
493  else if (max_ortdist == min_ortdist)
494  ortbin = 0;
495  else
496  ortbin = (int)(ortdist - min_ortdist) / (max_ortdist - min_ortdist) * (NBINS - 1);
497 
498  ort_profile.at(ortbin) += hit.charge;
499 
500  //////////////////////////////////////////////////////////////////////
501  //calculate the weight along the axis, this formula is based on rough guessology.
502  // there is no physics motivation behind the particular numbers, A.S.
503  // \todo: switch to something that is motivated and easier to
504  // spell than guessology. C.A.
505  ///////////////////////////////////////////////////////////////////////
506  double weight = (ortdist < 1.) ? 10 * (hit.charge) : 5 * (hit.charge) / ortdist;
507 
508  int fine_bin =
509  fProjectedLength ? (int)(linedist / fProjectedLength * (fProfileNbins - 1)) : 0;
510  int coarse_bin =
511  fProjectedLength ? (int)(linedist / fProjectedLength * (fCoarseNbins - 1)) : 0;
512 
513  if (fine_bin < fProfileNbins) //only fill if bin number is in range
514  {
515  fChargeProfile.at(fine_bin) += weight;
516 
517  //find maximum bin on the fly:
518  if (fChargeProfile.at(fine_bin) > current_maximum && fine_bin != 0 &&
519  fine_bin != fProfileNbins - 1) {
520  current_maximum = fChargeProfile.at(fine_bin);
521  fProfileMaximumBin = fine_bin;
522  }
523  }
524 
525  if (coarse_bin < fCoarseNbins) //only fill if bin number is in range
526  fCoarseChargeProfile.at(coarse_bin) += weight;
527 
528  } // end second loop on hits. Now should have filled profile vectors.
529 
530  if (verbose) std::cout << "end second loop " << std::endl;
531 
532  double integral = 0;
533  int ix = 0;
534  for (ix = 0; ix < NBINS; ix++) {
535  integral += ort_profile.at(ix);
536  if (integral >= 0.95 * fParams.sum_charge) { break; }
537  }
538 
539  fParams.width =
540  2 * (double)ix / (double)(NBINS - 1) *
541  (double)(max_ortdist -
542  min_ortdist); // multiply by two because ortdist is folding in both sides.
543 
544  if (verbose) std::cout << " after width " << std::endl;
545 
546  if (drawOrtHistos) {
547  TH1F* ortDistHist = new TH1F(
548  "ortDistHist", "Orthogonal Distance to axis;Distance (cm)", 100, min_ortdist, max_ortdist);
549  TH1F* ortDistHistCharge =
550  new TH1F("ortDistHistCharge",
551  "Orthogonal Distance to axis (Charge Weighted);Distance (cm)",
552  100,
553  min_ortdist,
554  max_ortdist);
555  TH1F* ortDistHistAndrzej = new TH1F("ortDistHistAndrzej",
556  "Orthogonal Distance Andrzej weighting",
557  100,
558  min_ortdist,
559  max_ortdist);
560 
561  for (int index = 0; index < fParams.N_Hits; index++) {
562  ortDistHist->Fill(ort_dist_vect.at(index));
563  ortDistHistCharge->Fill(ort_dist_vect.at(index), fHitVector.at(index).charge);
564  double weight = (ort_dist_vect.at(index) < 1.) ?
565  10 * (fHitVector.at(index).charge) :
566  (5 * (fHitVector.at(index).charge) / ort_dist_vect.at(index));
567  ortDistHistAndrzej->Fill(ort_dist_vect.at(index), weight);
568  }
569  ortDistHist->Scale(1.0 / ortDistHist->Integral());
570  ortDistHistCharge->Scale(1.0 / ortDistHistCharge->Integral());
571  ortDistHistAndrzej->Scale(1.0 / ortDistHistAndrzej->Integral());
572 
573  TCanvas* ortCanv = new TCanvas("ortCanv", "ortCanv", 600, 600);
574  ortCanv->cd();
575  ortDistHistAndrzej->SetLineColor(3);
576  ortDistHistAndrzej->Draw("");
577  ortDistHistCharge->Draw("same");
578  ortDistHist->SetLineColor(2);
579  ortDistHist->Draw("same");
580 
581  TLegend* leg = new TLegend(.4, .6, .8, .9);
582  leg->SetHeader("Charge distribution");
583  leg->AddEntry(ortDistHist, "Unweighted Hits");
584  leg->AddEntry(ortDistHistCharge, "Charge Weighted Hits");
585  leg->AddEntry(ortDistHistAndrzej, "Charge Weighted by Guessology");
586  leg->Draw();
587 
588  ortCanv->Update();
589  }
590 
593 
594  //calculate the forward and backward integrals counting int the maximum bin.
595 
596  for (int ibin = 0; ibin < fProfileNbins; ibin++) {
599  }
600 
601  // now, we have the forward and backward integral so start
602  // stepping forward and backward to find the trunk of the
603  // shower. This is done making sure that the running
604  // integral until less than 1% is left and the bin is above
605  // a set theshold many empty bins.
606 
607  //forward loop
608  double running_integral = fProfileIntegralForward;
609  int startbin, endbin;
610  for (startbin = fProfileMaximumBin; startbin > 1 && startbin < (int)(fChargeProfile.size());
611  startbin--) {
612  running_integral -= fChargeProfile.at(startbin);
613  if (fChargeProfile.at(startbin) < fChargeCutoffThreshold.at(fPlane) &&
614  running_integral / fProfileIntegralForward < 0.01)
615  break;
616  else if (fChargeProfile.at(startbin) < fChargeCutoffThreshold.at(fPlane) &&
617  startbin - 1 > 0 &&
618  fChargeProfile.at(startbin - 1) < fChargeCutoffThreshold.at(fPlane) &&
619  startbin - 2 > 0 &&
620  fChargeProfile.at(startbin - 2) < fChargeCutoffThreshold.at(fPlane))
621  break;
622  }
623 
624  //backward loop
625  running_integral = fProfileIntegralBackward;
626  for (endbin = fProfileMaximumBin; endbin > 0 && endbin < fProfileNbins - 1; endbin++) {
627  running_integral -= fChargeProfile.at(endbin);
628  if (fChargeProfile.at(endbin) < fChargeCutoffThreshold.at(fPlane) &&
629  running_integral / fProfileIntegralBackward < 0.01)
630  break;
631  else if (fChargeProfile.at(endbin) < fChargeCutoffThreshold.at(fPlane) &&
632  endbin + 1 < fProfileNbins - 1 && endbin + 2 < fProfileNbins - 1 &&
633  fChargeProfile.at(endbin + 1) < fChargeCutoffThreshold.at(fPlane) &&
634  fChargeProfile.at(endbin + 2) < fChargeCutoffThreshold.at(fPlane))
635  break;
636  }
637 
638  // now have profile start and endpoints. Now translate to wire/time.
639  // Will use wire/time that are on the rough axis.
640  // fProjectedLength is the length from fInterLow to interhigh a
641  // long the rough_2d_axis on bin distance is:
642  // util::PxPoint OnlinePoint;
643 
644  if (inv_2d_slope != -999999) //almost all cases
645  {
646  double ort_intercept_begin =
647  fBeginIntercept + (fEndIntercept - fBeginIntercept) / fProfileNbins * startbin;
649  fRough2DSlope, fRough2DIntercept, ort_intercept_begin, fRoughBeginPoint);
650 
651  double ort_intercept_end =
652  fBeginIntercept + (fEndIntercept - fBeginIntercept) / fProfileNbins * endbin;
654 
656  fRough2DSlope, fRough2DIntercept, ort_intercept_end, fRoughEndPoint);
658  }
659  else {
660  double wire_begin = WireLow + (WireHigh - WireLow) / fProfileNbins * startbin;
661 
662  util::PxPoint ptemphigh(fPlane, wire_begin, (fInterHigh_side + fInterLow_side) / 2);
665 
666  double wire_end = WireLow + (WireHigh - WireLow) / fProfileNbins * endbin;
667 
668  util::PxPoint ptemplow(fPlane, wire_end, (fInterHigh_side + fInterLow_side) / 2);
671  }
672 
673  if (verbose)
674  std::cout << " rough start points " << fRoughBeginPoint.w << ", " << fRoughBeginPoint.t
675  << " end: " << fRoughEndPoint.w << " " << fRoughEndPoint.t << std::endl;
676 
677  // ok. now have fRoughBeginPoint and fRoughEndPoint. No decision about direction has been made yet.
680 
682 
683  fTimeRecord_ProcName.push_back("GetProfileInfo");
684  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
685  }
686 
687  /**
688  * Calculates the following variables:
689  * length
690  * width
691  * hit_density_1D
692  * hit_density_2D
693  * direction
694  * @param override [description]
695  */
696  void
698  {
699  TStopwatch localWatch;
700  localWatch.Start();
701 
702  // need to define physical direction with openind angles and pass that to Ryan's line finder.
703 
704  // Direction decision has been moved entirely to Refine Direction()
705  // opening and closing angles are already set
706  // they are associated with the start and endpoints.
707  // If refine direction opted to switch the shower direction,
708  // it also switched opening and closing angles.
709 
710  /////////////////////////////////
711  // fParams.direction= ...
712 
713  // Direction is decided by RefineDirection()
714  util::PxPoint startHit, endHit;
715  startHit = fRoughBeginPoint;
716  endHit = fRoughEndPoint;
717 
718  ///////////////////////////////
719  //Ryan's Shower Strip finder work here.
720  //First we need to define the strip width that we want
721  double d = 0.6; //this is the width of the strip.... this needs to be tuned to something.
722  //===============================================================================================================
723  // Will need to feed in the set of hits that we want.
724  // const std::vector<util::PxHit*> whole;
725  std::vector<const util::PxHit*> subhit;
726  double linearlimit = 8;
727  double ortlimit = 12;
728  double lineslopetest(fRough2DSlope);
729  util::PxHit averageHit;
730  //also are we sure this is actually doing what it is supposed to???
731  //are we sure this works?
732  //is anybody even listening? Were they eaten by a grue?
733  gser.SelectLocalHitlist(
734  fHitVector, subhit, startHit, linearlimit, ortlimit, lineslopetest, averageHit);
735 
736  if (!(subhit.size()) || subhit.size() < 3) {
737  if (verbose)
738  std::cout << "Subhit list is empty or too small. Using rough start/end points..."
739  << std::endl;
742 
744 
745  fTimeRecord_ProcName.push_back("RefineStartPoints");
746  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
747 
748  return;
749  }
750 
751  double avgwire = averageHit.w;
752  double avgtime = averageHit.t;
753  //vertex in tilda-space pair(x-til,y-til)
754  std::vector<std::pair<double, double>> vertil;
755  //vector of the sum of the distance of a vector to every vertex in tilda-space
756  std::vector<double> vs;
757  // $$This needs to be corrected//this is the good hits that are between strip
758  std::vector<const util::PxHit*> ghits;
759  ghits.reserve(subhit.size());
760  int n = 0;
761  double fardistcurrent = 0;
762  util::PxHit startpoint;
763  double gwiretime = 0;
764  double gwire = 0;
765  double gtime = 0;
766  double gwirewire = 0;
767  //KAZU!!! I NEED this too//this will need to come from somewhere...
768  //This is some hit that is from the way far side of the entire shower cluster...
769  util::PxPoint farhit = fRoughEndPoint;
770 
771  //=============building the local vector===================
772  //this is clumsy... but just want to get something to work for now.
773  //THIS is REAL Horrible and CLUMSY... We can make things into helper functions soon.
774  std::vector<util::PxHit> returnhits;
775  std::vector<double> radiusofhit;
776  std::vector<int> closehits;
777  //==============================================================================
778  //Now we need to do the transformation into "tilda-space"
779  for (unsigned int a = 0; a < subhit.size(); a++) {
780  for (unsigned int b = a + 1; b < subhit.size(); b++) {
781  if (subhit.at(a)->w != subhit.at(b)->w) {
782  double xtil = ((subhit.at(a)->t - avgtime) - (subhit.at(b)->t - avgtime));
783  xtil /= ((subhit.at(a)->w - avgwire) - (subhit.at(b)->w - avgwire));
784  double ytil = (subhit.at(a)->w - avgwire) * xtil - (subhit.at(a)->t - avgtime);
785  //push back the tilda vertex point on the pair
786  std::pair<double, double> tv(xtil, ytil);
787  vertil.push_back(tv);
788  } //if Wires are not the same
789  } //for over b
790  } //for over a
791  // look at the distance from a tilda-vertex to all other tilda-verticies
792  for (unsigned int z = 0; z < vertil.size(); z++) {
793  double vd = 0; //the distance for vertex... just needs to be something 0
794  for (unsigned int c = 0; c < vertil.size(); c++)
795  vd += std::hypot(vertil.at(z).first - vertil.at(c).first,
796  vertil.at(z).second - vertil.at(c).second);
797  vs.push_back(vd);
798  vd = 0.0;
799  } //for z loop
800 
801  if (vs.size() == 0) //al hits on same wire?!
802  {
803  if (verbose)
804  std::cout << "vertil list is empty. all subhits are on the same wire?" << std::endl;
807 
809 
810  fTimeRecord_ProcName.push_back("RefineStartPoints");
811  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
812  return;
813  }
814  //need to find the min of the distance of vertex in tilda-space
815  //this will get me the area where things are most linear
816  int minvs = std::min_element(vs.begin(), vs.end()) - vs.begin();
817  // now use the min position to find the vertex in tilda-space
818  //now need to look a which hits are between the tilda lines from the points
819  //in the tilda space everything in wire time is based on the new origin which is at the average wire/time
820  double tilwire = vertil.at(minvs).first; //slope
821  double tiltimet = -vertil.at(minvs).second +
822  d * sqrt(1 + pow(tilwire, 2)); //negative cept is accounted for top strip
823  double tiltimeb = -vertil.at(minvs).second -
824  d * sqrt(1 + pow(tilwire, 2)); //negative cept is accounted for bottom strip
825  // look over the subhit list and ask for which are inside of the strip
826  for (unsigned int a = 0; a < subhit.size(); a++) {
827  double dtstrip =
828  (-tilwire * (subhit.at(a)->w - avgwire) + (subhit.at(a)->t - avgtime) - tiltimet) /
829  sqrt(tilwire * tilwire + 1);
830  double dbstrip =
831  (-tilwire * (subhit.at(a)->w - avgwire) + (subhit.at(a)->t - avgtime) - tiltimeb) /
832  sqrt(tilwire * tilwire + 1);
833 
834  if ((dtstrip < 0.0 && dbstrip > 0.0) || (dbstrip < 0.0 && dtstrip > 0.0)) {
835  ghits.push_back(subhit.at(a));
836  } //if between strip
837  } //for a loop
838  //=======================Do stuff with the good hits============================
839 
840  //of these small set of hits just fit a simple line.
841  //then we will need to see which of these hits is farthest away
842 
843  for (unsigned int g = 0; g < ghits.size(); g++) {
844  // should call the helper funtion to do the fit
845  //but for now since there is no helper function I will just write it here......again
846  n += 1;
847  gwiretime += ghits.at(g)->w * ghits.at(g)->t;
848  gwire += ghits.at(g)->w;
849  gtime += ghits.at(g)->t;
850  gwirewire += ghits.at(g)->w * ghits.at(g)->w;
851  // now work on calculating the distance in wire time space from the far point
852  //farhit needs to be a hit that is given to me
853  double fardist = std::hypot(ghits.at(g)->w - farhit.w, ghits.at(g)->t - farhit.t);
854  //come back to this... there is a better way to do this probably in the loop
855  //there should also be a check that the hit that is farthest away has subsequent hits after it on a few wires
856  if (fardist > fardistcurrent) {
857  fardistcurrent = fardist;
858  //if fardist... this is the point to use for the start point
859  startpoint.t = ghits.at(g)->t;
860  startpoint.w = ghits.at(g)->w;
861  startpoint.plane = ghits.at(g)->plane;
862  startpoint.charge = ghits.at(g)->charge;
863  }
864  } //for ghits loop
865  //This can be the new start point
866 
867  //should this be here? Id argue this might be done once outside.
868  fParams.length = gser.Get2DDistance((util::PxPoint*)&startpoint, &endHit);
869  if (fParams.length)
871  else
873 
874  if (fParams.length && fParams.width)
876  else
878 
879  fParams.start_point = startpoint;
880 
882 
883  fTimeRecord_ProcName.push_back("RefineStartPoints");
884  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
885  }
886 
887  ///////////////////////////////////////////////////////
888  ////
889  /////////////////////////////////////////////////////////////
890 
891  void
893  {
894  if (!override) { //Override being set, we skip all this logic.
895  //OK, no override. Stop if we're already finshed.
896  if (fFinishedGetFinalSlope) return;
897  //Try to run the previous function if not yet done.
899  }
900  else {
901  //Try to run the previous function if not yet done.
903  }
904 
905  TStopwatch localWatch;
906  localWatch.Start();
907  /**
908  * Calculates the following variables:
909  * angle_2d
910  * modified_hit_density
911  */
912 
913  const int NBINS = 720;
914  std::vector<int> fh_omega_single(NBINS, 0); //720,-180., 180.
915 
916  double current_maximum = 0;
917  double curr_max_bin = -1;
918 
919  for (auto hit : fHitVector) {
920 
921  double omx = gser.Get2Dangle((util::PxPoint*)&hit,
922  &fParams.start_point); // in rad and assuming cm/cm space.
923  int nbin = (omx + TMath::Pi()) * (NBINS - 1) / (2 * TMath::Pi());
924  if (nbin >= NBINS) nbin = NBINS - 1;
925  if (nbin < 0) nbin = 0;
926  fh_omega_single[nbin] += hit.charge;
927 
928  if (fh_omega_single[nbin] > current_maximum) {
929  current_maximum = fh_omega_single[nbin];
930  curr_max_bin = nbin;
931  }
932  }
933  // FIXME: using two different definitions of PI in the same calculation?
934  // 2022-04-18 CHG
935  fParams.angle_2d = (curr_max_bin / 720 * (2 * TMath::Pi())) - TMath::Pi();
936  fParams.angle_2d *= 180 / PI;
937  if (verbose) std::cout << " Final 2D angle: " << fParams.angle_2d << " degrees " << std::endl;
938 
939  double mod_angle =
940  (fabs(fParams.angle_2d) <= 90) ?
941  fabs(fParams.angle_2d) :
942  180 - fabs(fParams.angle_2d); //want to transfer angle to radians and from 0 to 90.
943 
944  if (
945  cos(
946  mod_angle * PI /
947  180.)) { //values here based on fit of hit_density_1D vs. mod_angle defined as above (in ArgoNeuT).
948  if (mod_angle <= 75.)
950  fParams.hit_density_1D / (2.45 * cos(mod_angle * PI / 180.) + 0.2595);
951  else
953  fParams.hit_density_1D / (1.036 * mod_angle * PI / 180. - 0.2561);
954 
955  //calculate modified mean_charge:
956  fParams.modmeancharge = fParams.mean_charge / (1264 - 780 * cos(mod_angle * PI / 180.));
957  }
958  else
960 
961  /////////////////// testing
962  // do not use for now.
963  bool drawProfileHisto = false;
964  if (drawProfileHisto) {
965 
967  double corr_factor = 1;
968  if (
969  cos(
970  mod_angle * PI /
971  180.)) { //values here based on fit of hit_density_1D vs. mod_angle defined as above (in ArgoNeuT).
972 
973  if (mod_angle <= 75.)
974  corr_factor *= (2.45 * cos(mod_angle * PI / 180.) + 0.2595);
975  else
976  corr_factor *= (1.036 * mod_angle * PI / 180. - 0.2561);
977  }
978 
979  fProfileNbins =
980  (int)(fProfileNbins / 2 * corr_factor + 0.5); // +0.5 to round to nearest sensible value
981 
982  if (verbose) std::cout << " number of final profile bins " << fProfileNbins << std::endl;
983 
984  fChargeProfile.clear();
985  fCoarseChargeProfile.clear();
986  fChargeProfile.resize(fProfileNbins, 0);
988 
989  fChargeProfileNew.clear();
990  fChargeProfileNew.resize(fProfileNbins, 0);
991 
992  util::PxPoint BeginOnlinePoint;
994 
995  for (auto hit : fHitVector) {
996 
997  util::PxPoint OnlinePoint;
998  // get coordinates of point on axis.
999  gser.GetPointOnLine(fRough2DSlope, &BeginOnlinePoint, &hit, OnlinePoint);
1000 
1001  double linedist = gser.Get2DDistance(&OnlinePoint, &BeginOnlinePoint);
1002  int fine_bin = (int)(linedist / fProjectedLength * (fProfileNbins - 1));
1003 
1004  if (fine_bin < fProfileNbins) //only fill if bin number is in range
1005  {
1006  fChargeProfileNew.at(fine_bin) += hit.charge;
1007  }
1008  }
1009 
1010  TH1F* charge_histo = new TH1F("charge dist", "charge dist", fProfileNbins, 0, fProfileNbins);
1011 
1012  TH1F* charge_diff = new TH1F("charge diff", "charge diff", fProfileNbins, 0, fProfileNbins);
1013 
1014  for (int ix = 0; ix < fProfileNbins - 1; ix++) {
1015  charge_histo->SetBinContent(ix, fChargeProfileNew[ix]);
1016 
1017  if (ix > 2 && ix < fProfileNbins - 3) {
1018  double diff =
1019  (1. / 12. * fChargeProfileNew[ix - 2] - 2. / 3. * fChargeProfileNew[ix - 1] +
1020  2. / 3. * fChargeProfileNew[ix + 1] - 1. / 12. * fChargeProfileNew[ix + 2]) /
1021  (double)fProfileNbins;
1022  charge_diff->SetBinContent(ix, diff);
1023  }
1024  }
1025 
1026  TCanvas* chCanv = new TCanvas("chCanv", "chCanv", 600, 600);
1027  chCanv->cd();
1028  charge_histo->SetLineColor(3);
1029  charge_histo->Draw("");
1030  charge_diff->SetLineColor(2);
1031  charge_diff->Draw("same");
1032  chCanv->Update();
1033  }
1034 
1035  /// end testing
1036 
1037  fTimeRecord_ProcName.push_back("GetFinalSlope");
1038  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1039 
1040  fFinishedGetFinalSlope = true;
1041  return;
1042  }
1043 
1044  /**
1045  * This function calculates the opening/closing angles
1046  * It also makes a decision on whether or not to flip the direction
1047  * the cluster. Then the start point is refined later.
1048  * @param override [description]
1049  */
1050  void
1052  {
1053  //
1054  // We don't use "override"? Should we remove? 05/01/14
1055  //
1056  if (!override) override = true;
1057 
1058  TStopwatch localWatch;
1059  localWatch.Start();
1060 
1061  // if(!override) { //Override being set, we skip all this logic.
1062  // //OK, no override. Stop if we're already finshed.
1063  // if (fFinishedRefineDirection) return;
1064  // //Try to run the previous function if not yet done.
1065  // if (!fFinishedGetProfileInfo) GetProfileInfo(true);
1066  // } else {
1067  // //Try to run the previous function if not yet done.
1068  // if (!fFinishedGetProfileInfo) GetProfileInfo(true);
1069  // }
1070 
1071  // double wire_2_cm = gser.WireToCm();
1072  // double time_2_cm = gser.TimeToCm();
1073 
1074  // Removing the conversion since these points are set above in cm-cm space
1075  // -- Corey
1076 
1077  util::PxPoint this_startPoint, this_endPoint;
1078 
1080  this_startPoint = fParams.start_point;
1081  this_endPoint = fParams.end_point;
1082  }
1083  else {
1084  this_startPoint = fRoughBeginPoint;
1085  this_endPoint = fRoughEndPoint;
1086  }
1087  if (verbose) {
1088  std::cout << "Angle: Start point: (" << this_startPoint.w << ", " << this_startPoint.t
1089  << ")\n";
1090  std::cout << "Angle: End point : (" << this_endPoint.w << ", " << this_endPoint.t << ")\n";
1091  }
1092  double endStartDiff_x = (this_endPoint.w - this_startPoint.w);
1093  double endStartDiff_y = (this_endPoint.t - this_startPoint.t);
1094  double rms_forward = 0;
1095  double rms_backward = 0;
1096  double mean_forward = 0;
1097  double mean_backward = 0;
1098  //double weight_total = 0;
1099  double hit_counter_forward = 0;
1100  double hit_counter_backward = 0;
1101 
1102  if (verbose && endStartDiff_y == 0 && endStartDiff_x == 0) {
1103  std::cerr << "Error: end point and start point are the same!\n";
1104 
1105  fTimeRecord_ProcName.push_back("RefineDirection");
1106  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1107  return;
1108  }
1109 
1110  double percentage = 0.90;
1111  double percentage_HC = 0.90 * fParams.N_Hits_HC / fParams.N_Hits;
1112  const int NBINS = 200;
1113  const double wgt = 1.0 / fParams.N_Hits;
1114 
1115  // Containers for the angle histograms
1116  std::vector<float> opening_angle_bin(NBINS, 0.0);
1117  std::vector<float> closing_angle_bin(NBINS, 0.0);
1118  std::vector<float> opening_angle_highcharge_bin(NBINS, 0.0);
1119  std::vector<float> closing_angle_highcharge_bin(NBINS, 0.0);
1120  std::vector<float> opening_angle_chargeWgt_bin(NBINS, 0.0);
1121  std::vector<float> closing_angle_chargeWgt_bin(NBINS, 0.0);
1122  //hard coding this for now, should use SetRefineDirectionQMin function
1123  fQMinRefDir = 25;
1124 
1125  int index = -1;
1126  for (auto& hit : fHitVector) {
1127  index++;
1128  //skip this hit if below minimum cutoff param
1129  if (hit.charge < fQMinRefDir) continue;
1130 
1131  //weight_total = hit.charge;
1132 
1133  // Compute forward mean
1134  double hitStartDiff_x = (hit.w - this_startPoint.w);
1135  double hitStartDiff_y = (hit.t - this_startPoint.t);
1136 
1137  if (hitStartDiff_x == 0 && hitStartDiff_y == 0) continue;
1138 
1139  double cosangle_start = (endStartDiff_x * hitStartDiff_x + endStartDiff_y * hitStartDiff_y);
1140  cosangle_start /= (pow(pow(endStartDiff_x, 2) + pow(endStartDiff_y, 2), 0.5) *
1141  pow(pow(hitStartDiff_x, 2) + pow(hitStartDiff_y, 2), 0.5));
1142 
1143  if (cosangle_start > 0) {
1144  // Only take into account for hits that are in "front"
1145  //no weighted average, works better as flat average w/ min charge cut
1146  mean_forward += cosangle_start;
1147  rms_forward += pow(cosangle_start, 2);
1148  hit_counter_forward++;
1149  }
1150 
1151  // Compute backward mean
1152  double hitEndDiff_x = (hit.w - this_endPoint.w);
1153  double hitEndDiff_y = (hit.t - this_endPoint.t);
1154  if (hitEndDiff_x == 0 && hitEndDiff_y == 0) continue;
1155 
1156  double cosangle_end = (endStartDiff_x * hitEndDiff_x + endStartDiff_y * hitEndDiff_y) * (-1.);
1157  cosangle_end /= (pow(pow(endStartDiff_x, 2) + pow(endStartDiff_y, 2), 0.5) *
1158  pow(pow(hitEndDiff_x, 2) + pow(hitEndDiff_y, 2), 0.5));
1159 
1160  if (cosangle_end > 0) {
1161  //no weighted average, works better as flat average w/ min charge cut
1162  mean_backward += cosangle_end;
1163  rms_backward += pow(cosangle_end, 2);
1164  hit_counter_backward++;
1165  }
1166 
1167  int N_bins_OPEN = (NBINS - 1) * acos(cosangle_start) / PI;
1168  int N_bins_CLOSE = (NBINS - 1) * acos(cosangle_end) / PI;
1169  if (N_bins_OPEN < 0) N_bins_OPEN = 0;
1170  if (N_bins_CLOSE < 0) N_bins_CLOSE = 0;
1171 
1172  opening_angle_chargeWgt_bin[N_bins_OPEN] += hit.charge / fParams.sum_charge;
1173  closing_angle_chargeWgt_bin[N_bins_CLOSE] += hit.charge / fParams.sum_charge;
1174  opening_angle_bin[N_bins_OPEN] += wgt;
1175  closing_angle_bin[N_bins_CLOSE] += wgt;
1176 
1177  //Also fill bins for high charge hits
1178  if (hit.charge > fParams.mean_charge) {
1179  opening_angle_highcharge_bin[N_bins_OPEN] += wgt;
1180  closing_angle_highcharge_bin[N_bins_CLOSE] += wgt;
1181  }
1182  }
1183 
1184  //initialize iterators for the angles
1185  int iBin(0), jBin(0), kBin(0), lBin(0), mBin(0), nBin(0);
1186  double percentage_OPEN(0.0), percentage_CLOSE(0.0), percentage_OPEN_HC(0.0),
1187  percentage_CLOSE_HC(0.0), //The last 2 are for High Charge (HC)
1188  percentage_OPEN_CHARGEWGT(0.0), percentage_CLOSE_CHARGEWGT(0.0);
1189 
1190  for (iBin = 0; percentage_OPEN <= percentage && iBin < NBINS; iBin++) {
1191  percentage_OPEN += opening_angle_bin[iBin];
1192  }
1193 
1194  for (jBin = 0; percentage_CLOSE <= percentage && jBin < NBINS; jBin++) {
1195  percentage_CLOSE += closing_angle_bin[jBin];
1196  }
1197 
1198  for (kBin = 0; percentage_OPEN_HC <= percentage_HC && kBin < NBINS; kBin++) {
1199  percentage_OPEN_HC += opening_angle_highcharge_bin[kBin];
1200  }
1201 
1202  for (lBin = 0; percentage_CLOSE_HC <= percentage_HC && lBin < NBINS; lBin++) {
1203  percentage_CLOSE_HC += closing_angle_highcharge_bin[lBin];
1204  }
1205 
1206  for (mBin = 0; percentage_OPEN_CHARGEWGT <= percentage && mBin < NBINS; mBin++) {
1207  percentage_OPEN_CHARGEWGT += opening_angle_chargeWgt_bin[mBin];
1208  }
1209 
1210  for (nBin = 0; percentage_CLOSE_CHARGEWGT <= percentage && nBin < NBINS; nBin++) {
1211  percentage_CLOSE_CHARGEWGT += closing_angle_chargeWgt_bin[nBin];
1212  }
1213 
1214  double opening_angle = iBin * PI / NBINS;
1215  double closing_angle = jBin * PI / NBINS;
1216  double opening_angle_highcharge = kBin * PI / NBINS;
1217  double closing_angle_highcharge = lBin * PI / NBINS;
1218  double opening_angle_charge_wgt = mBin * PI / NBINS;
1219  double closing_angle_charge_wgt = nBin * PI / NBINS;
1220 
1221  double value_1 = closing_angle / opening_angle - 1;
1222  double value_2 = (fCoarseChargeProfile[0] / fCoarseChargeProfile[1]);
1223  if (value_2 < 100.0 && value_2 > 0.01)
1224  value_2 = log(value_2);
1225  else
1226  value_2 = 0.0;
1227  double value_3 = closing_angle_charge_wgt / opening_angle_charge_wgt - 1;
1228 
1229  // Using a sigmoid function to determine flipping.
1230  // I am going to weigh all of the values above (1, 2, 3) equally.
1231  // But, since value 2 in particular, I'm going to cut things off at 5.
1232 
1233  if (value_1 > 3) value_1 = 3.0;
1234  if (value_1 < -3) value_1 = -3.0;
1235  if (value_2 > 3) value_2 = 3.0;
1236  if (value_2 < -3) value_2 = -3.0;
1237  if (value_3 > 3) value_3 = 3.0;
1238  if (value_3 < -3) value_3 = -3.0;
1239 
1240  double Exp = exp(value_1 + value_2 + value_3);
1241  double sigmoid = (Exp - 1) / (Exp + 1);
1242 
1243  bool flip = false;
1244  if (sigmoid < 0) flip = true;
1245  if (flip) {
1246  if (verbose) std::cout << "Flipping!" << std::endl;
1247  std::swap(opening_angle, closing_angle);
1248  std::swap(opening_angle_highcharge, closing_angle_highcharge);
1249  std::swap(opening_angle_charge_wgt, closing_angle_charge_wgt);
1250  std::swap(fParams.start_point, fParams.end_point);
1251  std::swap(fRoughBeginPoint, fRoughEndPoint);
1252  }
1253  else if (verbose) {
1254  std::cout << "Not Flipping!\n";
1255  }
1256 
1257  fParams.opening_angle = opening_angle;
1258  fParams.opening_angle_charge_wgt = opening_angle_charge_wgt;
1259  fParams.closing_angle = closing_angle;
1260  fParams.closing_angle_charge_wgt = closing_angle_charge_wgt;
1261 
1262  fFinishedRefineDirection = true;
1263 
1264  fTimeRecord_ProcName.push_back("RefineDirection");
1265  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1266  } //end RefineDirection
1267 
1268  void
1270  {
1271  TStopwatch localWatch;
1272  localWatch.Start();
1273 
1274  if (fHitVector.size()) {
1275  std::vector<const util::PxHit*> container_polygon;
1276  gser.SelectPolygonHitList(fHitVector, container_polygon);
1277  //now making Polygon Object
1278  std::pair<float, float> tmpvertex;
1279  //make Polygon Object as in mac/PolyOverlap.cc
1280  std::vector<std::pair<float, float>> vertices;
1281  for (unsigned int i = 0; i < container_polygon.size(); i++) {
1282  tmpvertex = std::make_pair(container_polygon.at(i)->w, container_polygon.at(i)->t);
1283  vertices.push_back(tmpvertex);
1284  }
1285  fParams.PolyObject = Polygon2D(vertices);
1286  }
1287 
1288  fTimeRecord_ProcName.push_back("FillPolygon");
1289  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1290  }
1291 
1292  ///////////////////////////////////////////////////////////////////////////////////
1293 
1294  void
1296  {
1297  // This function is meant to pick the direction.
1298  // It refines both the start and end point, and then asks
1299  // if it should flip.
1300 
1301  TStopwatch localWatch;
1302  localWatch.Start();
1303 
1304  if (verbose) std::cout << " here!!! " << std::endl;
1305 
1306  if (!override) { //Override being set, we skip all this logic.
1307  //OK, no override. Stop if we're already finshed.
1309  fTimeRecord_ProcName.push_back("RefineStartPointAndDirection");
1310  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1311  return;
1312  }
1313  //Try to run the previous function if not yet done.
1314  if (!fFinishedGetProfileInfo) GetProfileInfo(gser, true);
1315  }
1316  else {
1317  //Try to run the previous function if not yet done.
1318  if (!fFinishedGetProfileInfo) GetProfileInfo(gser, true);
1319  }
1320  if (verbose) {
1321  std::cout << "REFINING .... " << std::endl;
1322  std::cout << " Rough start and end point: " << std::endl;
1323  std::cout << " s: (" << fParams.start_point.w << ", " << fParams.start_point.t << ")"
1324  << std::endl;
1325  std::cout << " e: (" << fParams.end_point.w << ", " << fParams.end_point.t << ")"
1326  << std::endl;
1327  }
1328  RefineStartPoints(gser);
1329  if (verbose) {
1330  std::cout << " Once Refined start and end point: " << std::endl;
1331  std::cout << " s: (" << fParams.start_point.w << ", " << fParams.start_point.t << ")"
1332  << std::endl;
1333  std::cout << " e: (" << fParams.end_point.w << ", " << fParams.end_point.t << ")"
1334  << std::endl;
1335  std::swap(fParams.start_point, fParams.end_point);
1336  std::swap(fRoughBeginPoint, fRoughEndPoint);
1337  }
1338  RefineStartPoints(gser);
1339  if (verbose) {
1340  std::cout << " Twice Refined start and end point: " << std::endl;
1341  std::cout << " s: (" << fParams.start_point.w << ", " << fParams.start_point.t << ")"
1342  << std::endl;
1343  std::cout << " e: (" << fParams.end_point.w << ", " << fParams.end_point.t << ")"
1344  << std::endl;
1345  std::swap(fParams.start_point, fParams.end_point);
1346  std::swap(fRoughBeginPoint, fRoughEndPoint);
1347  }
1348  RefineDirection();
1349  if (verbose) {
1350  std::cout << " Final start and end point: " << std::endl;
1351  std::cout << " s: (" << fParams.start_point.w << ", " << fParams.start_point.t << ")"
1352  << std::endl;
1353  std::cout << " e: (" << fParams.end_point.w << ", " << fParams.end_point.t << ")"
1354  << std::endl;
1355  }
1357 
1358  fTimeRecord_ProcName.push_back("RefineStartPointAndDirection");
1359  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1360  }
1361 
1362  void
1364  {
1365  if (!override) return;
1366  if (verbose) std::cout << " ---- Trying T/S sep. ------ \n";
1367  if (enableFANN) {
1368  if (verbose) std::cout << " ---- Doing T/S sep. ------- \n";
1369  std::vector<float> FeatureVector, outputVector;
1370  GetFANNVector(FeatureVector);
1371  }
1372  else {
1373  if (verbose) std::cout << " ---- Failed T/S sep. ------ \n";
1374  }
1375  }
1376 
1377  //----------------------------------------------------------------------------
1378  void
1379  ClusterParamsAlg::GetEndCharges(util::GeometryUtilities const& gser, bool override_ /* = false */)
1380  {
1381  if (!override_) { //Override being set, we skip all this logic.
1382  //OK, no override. Stop if we're already finshed.
1383  if (fFinishedGetEndCharges) return;
1384  }
1385 
1386  TStopwatch localWatch;
1387  localWatch.Start();
1388 
1390  fParams.end_charge = EndCharge(gser);
1391 
1392  fFinishedGetEndCharges = true;
1393 
1394  fTimeRecord_ProcName.push_back("GetEndCharges");
1395  fTimeRecord_ProcTime.push_back(localWatch.RealTime());
1396  } // ClusterParamsAlg::GetEndCharges()
1397 
1398  //----------------------------------------------------------------------------
1399  double
1400  ClusterParamsAlg::LinearIntegral(double m, double q, double x1, double x2)
1401  {
1402  return m / 2. * (cet::square(x2) - cet::square(x1)) + q * (x2 - x1);
1403  }
1404 
1405  //----------------------------------------------------------------------------
1406  double
1408  double from_length,
1409  double to_length,
1410  unsigned int fit_first_bin,
1411  unsigned int fit_end_bin)
1412  {
1413  // first compute the information on the charge profile
1414  GetProfileInfo(gser);
1415 
1416  // first things first
1417  if (fit_first_bin > fit_end_bin) std::swap(fit_first_bin, fit_end_bin);
1418 
1419  // how many bins can we use?
1420  const unsigned int nbins =
1421  std::min(fit_end_bin - fit_first_bin, (unsigned int)fChargeProfile.size());
1422  if (nbins == 0) return 0;
1423 
1424  // move the specified range within reason
1425  if (fit_end_bin > fChargeProfile.size()) {
1426  fit_end_bin = fChargeProfile.size();
1427  fit_first_bin = fit_end_bin - nbins;
1428  }
1429 
1430  // if we have to use only one bin, so be it
1431  if (nbins < 2) return fChargeProfile[fit_first_bin];
1432 
1433  // fit the charge profile vs. bin number;
1434  // we assume that the first bin (#0) starts from the very beginning of the
1435  // cluster, that is at axis coordinate 0.
1437  for (unsigned int iBin = fit_first_bin; iBin < fit_end_bin; ++iBin) {
1438  // should we use a Poisson error instead of no error?
1439  fitter.add((double)iBin, fChargeProfile[iBin]);
1440  } // for
1441 
1442  // get the linear fit parameters; [0] intercept [1] slope
1444  try {
1445  fit = fitter.FitParameters();
1446  }
1447  catch (std::range_error const&) { // oops...
1448  // this is actually unexpected, since the only reason for the fit to fail
1449  // would be a singular fit matrix, that should be pretty much impossible
1450  // given that the bin coordinates are well behaved and there are at least
1451  // two of them
1452  std::cerr << "IntegrateFitCharge(): linear fit failed!" << std::endl;
1453  return 0.;
1454  }
1455 
1456  // now determine the bins corresponding to the length to integrate;
1457  // note that length can be pathologic (0, negative...); not our problem!
1458  const double length_to_bin = (double(fChargeProfile.size() - 1) / fProjectedLength);
1459  const double from_bin = from_length * length_to_bin, to_bin = to_length * length_to_bin;
1460 
1461  // we want the integral between from_bin and to_bin now:
1462  return LinearIntegral(fit[1] /* m */, fit[0] /* q */, from_bin, to_bin);
1463  } // ClusterParamsAlg::IntegrateFitCharge()
1464 
1465  //----------------------------------------------------------------------------
1466  double
1468  float length /* = 1. */,
1469  unsigned int nbins /* = 10 */)
1470  {
1471  switch (fHitVector.size()) {
1472  case 0: return 0.;
1473  case 1:
1474  return fHitVector.back().charge;
1475  // the "default" is the rest of the function
1476  } // switch
1477 
1478  // this is the range of the fit:
1479  const unsigned int fit_first_bin = 0, fit_last_bin = nbins;
1480 
1481  return IntegrateFitCharge(gser, 0., length, fit_first_bin, fit_last_bin);
1482  } // ClusterParamsAlg::StartCharge()
1483 
1484  //----------------------------------------------------------------------------
1485  double
1487  float length /* = 1. */,
1488  unsigned int nbins /* = 10 */)
1489  {
1490  switch (fHitVector.size()) {
1491  case 0: return 0.;
1492  case 1:
1493  return fHitVector.back().charge;
1494  // the "default" is the rest of the function
1495  } // switch
1496 
1497  // need the available number of bins and the axis length
1498  GetProfileInfo(gser);
1499  const unsigned int MaxBins = fChargeProfile.size();
1500 
1501  // this is the range of the fit:
1502  const unsigned int fit_first_bin = MaxBins > nbins ? MaxBins - nbins : 0,
1503  fit_last_bin = MaxBins;
1504 
1505  // now determine the integration range, in bin units;
1506  // get to the end, and go length backward;
1507  // note that length can be pathologic (0, negative...); not our problem!
1508  const double from = fProjectedLength - length, to = fProjectedLength;
1509 
1510  return IntegrateFitCharge(gser, from, to, fit_first_bin, fit_last_bin);
1511  } // ClusterParamsAlg::EndCharge()
1512 
1513  //------------------------------------------------------------------------------
1514  float
1516  {
1517  if (fHitVector.size() < 2) return 0.0F;
1518 
1519  // compute all the averages
1520  GetAverages();
1521 
1522  // return the relevant information
1523  return std::isnormal(fParams.N_Wires) ? fParams.multi_hit_wires / fParams.N_Wires : 0.;
1524  } // ClusterParamsAlg::MultipleHitWires()
1525 
1526  //------------------------------------------------------------------------------
1527  float
1529  {
1530  if (fHitVector.size() < 2) return 0.0F;
1531 
1532  // compute all the averages
1533  GetAverages();
1534  RefineStartPoints(gser); // fParams.length
1535 
1536  // return the relevant information
1537  return std::isnormal(fParams.length) ? fParams.multi_hit_wires / fParams.length : 0.;
1538  } // ClusterParamsAlg::MultipleHitDensity()
1539 
1540 } //end namespace
double rms_ADC
RMS (standard deviation of sample) of ADC counts of hits in ADC.
Definition: ClusterParams.h:31
process_name opflash particleana ie ie ie z
void GetRoughAxis(bool override=false)
2D polygon object
void GetFinalSlope(util::GeometryUtilities const &gser, bool override=false)
double closing_angle_charge_wgt
Same for charge_wgt.
Definition: ClusterParams.h:43
std::vector< double > fCoarseChargeProfile
std::vector< std::string > fTimeRecord_ProcName
double charge_wgt_y
Mean of hits along y, charge weighted.
Definition: ClusterParams.h:37
void FillPolygon(util::GeometryUtilities const &gser)
Double_t Get2DDistance(Double_t wire1, Double_t time1, Double_t wire2, Double_t time2) const
BEGIN_PROLOG could also be cerr
process_name cluster
Definition: cheaterreco.fcl:51
Double_t Get2Dangle(Double_t deltawire, Double_t deltatime) const
double mean_x
Mean of hits along x, peaks only.
Definition: ClusterParams.h:32
double start_charge
Charge at the start of the cluster.
Definition: ClusterParams.h:44
double IntegrateFitCharge(util::GeometryUtilities const &gser, double from_length, double to_length, unsigned int fit_first_bin, unsigned int fit_end_bin)
Integrates the charge between two positions in the cluster axis.
pdgs p
Definition: selectors.fcl:22
double rms_charge
RMS (standard deviation of sample) of charge of hits in ADC.
Definition: ClusterParams.h:28
void GetProfileInfo(util::GeometryUtilities const &gser, bool override=false)
Polygon2D PolyObject
Polygon Object...see Polygon2D.hh.
Definition: ClusterParams.h:21
std::vector< util::PxHit > fHitVector
double eigenvalue_principal
the principal eigenvalue from PCA
Definition: ClusterParams.h:46
double mean_ADC
Mean (average) of ADC counts of hits, in ADC.
Definition: ClusterParams.h:30
process_name hit
Definition: cheaterreco.fcl:51
BEGIN_PROLOG g
double mean_y
Mean of hits along y, peaks only.
Definition: ClusterParams.h:33
Classes gathering simple statistics.
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
Weight_t RMS() const
Returns the root mean square.
constexpr double PI
Classes performing simple fits.
void RefineStartPoints(util::GeometryUtilities const &gser)
Performs a linear regression of data.
Definition: SimpleFits.h:849
Weight_t Average() const
Returns the value average.
Class def header for exception classes in ClusterRecoUtil package.
process_name gaushit a
while getopts h
double eigenvalue_secondary
the secondary eigenvalue from PCA
Definition: ClusterParams.h:47
float MultipleHitDensity(util::GeometryUtilities const &gser)
Returns the number of multiple hits per wire.
double sum_charge
Sum charge of hits in ADC.
Definition: ClusterParams.h:26
bool isfinite(Vector const &v)
Returns whether all components of the vector are finite.
cluster::cluster_params fParams
Weight_t Sum() const
Returns the weighted sum of the values.
void GetEndCharges(util::GeometryUtilities const &gser, bool override_=false)
double t
Definition: PxUtils.h:11
void TrackShowerSeparation(bool override=false)
std::vector< double > fTimeRecord_ProcTime
util::PxPoint start_point
start point
Definition: ClusterParams.h:23
double angle_2d
Angle of axis in wire/hit view.
Definition: ClusterParams.h:39
double EndCharge(util::GeometryUtilities const &gser, float length=1., unsigned int nbins=10)
Returns the expected charge at the end of the cluster.
std::array< Data_t, NParams > FitParameters_t
type of set of fit parameters
Definition: SimpleFits.h:866
void GetAverages(bool override=false)
double w
Definition: PxUtils.h:10
std::vector< double > fChargeProfile
void RefineDirection(bool override=false)
Int_t GetPointOnLine(Double_t slope, Double_t intercept, Double_t wire1, Double_t time1, Double_t &wireout, Double_t &timeout) const
void RefineStartPointAndDirection(util::GeometryUtilities const &gser, bool override=false)
std::vector< double > fChargeProfileNew
double sum_ADC
Sum charge of ADC counts of hits, in ADC.
Definition: ClusterParams.h:29
util::PxPoint end_point
end point
Definition: ClusterParams.h:24
size_t fMinNHits
Cut value for # hits: below this value clusters are not evaluated.
double mean_charge
Mean (average) charge of hits in ADC.
Definition: ClusterParams.h:27
int SetHits(const std::vector< util::PxHit > &)
double charge_wgt_x
Mean of hits along x, charge weighted.
Definition: ClusterParams.h:36
double opening_angle_charge_wgt
Same for charge_wgt.
Definition: ClusterParams.h:41
double charge
area charge
Definition: PxUtils.h:39
virtual FitParameters_t FitParameters() const override
Computes and returns all the parameters of the fit result.
void FillParams(util::GeometryUtilities const &gser, bool override_DoGetAverages=false, bool override_DoGetRoughAxis=false, bool override_DoGetProfileInfo=false, bool override_DoRefineStartPointsAndDirection=false, bool override_DoGetFinalSlope=false, bool override_DoTrackShowerSep=false, bool override_DoEndCharge=false)
double opening_angle
Width of angular distubtion wrt vertx.
Definition: ClusterParams.h:40
void SelectPolygonHitList(const std::vector< util::PxHit > &hitlist, std::vector< const util::PxHit * > &hitlistlocal) const
Int_t GetPointOnLineWSlopes(Double_t slope, Double_t intercept, Double_t ort_intercept, Double_t &wireout, Double_t &timeout) const
double cluster_angle_2d
Linear best fit to high-charge hits in the cluster.
Definition: ClusterParams.h:38
double rms_x
rms of hits along x (wires)
Definition: ClusterParams.h:34
double rms_y
rms of hits along y, (time)
Definition: ClusterParams.h:35
bool add(Data_t x, Data_t y, Data_t sy=Data_t(1.0))
Definition: SimpleFits.h:366
static double LinearIntegral(double m, double q, double x1, double x2)
Returns the integral of f(x) = mx + q defined in [x1, x2].
double closing_angle
Width of angular distubtion wrt endpoint.
Definition: ClusterParams.h:42
double StartCharge(util::GeometryUtilities const &gser, float length=1., unsigned int nbins=10)
Returns the expected charge at the beginning of the cluster.
BEGIN_PROLOG could also be cout
unsigned int plane
Definition: PxUtils.h:12
void SelectLocalHitlist(const std::vector< util::PxHit > &hitlist, std::vector< const util::PxHit * > &hitlistlocal, util::PxPoint &startHit, Double_t &linearlimit, Double_t &ortlimit, Double_t &lineslopetest) const
std::vector< double > fChargeCutoffThreshold
float MultipleHitWires()
Returns the number of multiple hits per wire.
void add(Data_t value, Weight_t weight=Weight_t(1.0))
Adds one entry with specified value and weight.
void GetFANNVector(std::vector< float > &data)
double end_charge
Charge at the (other) end of the cluster.
Definition: ClusterParams.h:45