All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CFAlgoShowerCompat.cxx
Go to the documentation of this file.
1 #include "TFile.h"
2 #include "TTree.h"
3 
5 
6 namespace cmtool {
7 
8  //-------------------------------------------------------
10  //-------------------------------------------------------
11  {
12 
13  _fout_hax = 0;
14  _ana_tree = 0;
15 
16  if (!_fout_hax) _fout_hax = new TFile("fout_hax.root", "RECREATE");
17 
18  if (!_ana_tree) {
19  _ana_tree = new TTree("ana_tree", "ana_tree");
20  _ana_tree->Branch("o_ang_avg", &_o_ang_avg, "o_ang_avg/D");
21  _ana_tree->Branch("o_ang_rms", &_o_ang_rms, "o_ang_rms/D");
22  _ana_tree->Branch("o_ang_wt_avg", &_o_ang_wt_avg, "o_ang_wt_avg/D");
23  _ana_tree->Branch("o_ang_wt_rms", &_o_ang_wt_rms, "o_ang_wt_rms/D");
24  _ana_tree->Branch("max_trackness", &_max_trackness, "max_trackness/D");
25  _ana_tree->Branch("max_len_over_width", &_max_len_over_width, "max_len_over_width/D");
26  _ana_tree->Branch("min_oa_over_len", &_min_oa_over_len, "min_oa_over_len/D");
27  _ana_tree->Branch(
28  "max_poly_perim_over_A", &_max_poly_perim_over_A, "max_poly_perim_over_A/D");
29  _ana_tree->Branch("min_modhitdens", &_min_modhitdens, "min_modhitdens/D");
30  }
31  }
32 
33  //-----------------------------
34  void
36  //-----------------------------
37  {}
38 
39  //----------------------------------------------------------------------------------------------
40  float
42  const std::vector<const cluster::ClusterParamsAlg*>& clusters)
43  //----------------------------------------------------------------------------------------------
44  {
45  _o_ang_avg = 0;
46  _o_ang_rms = 0;
47  _o_ang_wt_avg = 0;
48  _o_ang_wt_rms = 0;
49  _max_trackness = -9999.;
50  _max_len_over_width = -9999.;
51  _min_oa_over_len = 9999.;
52  _max_poly_perim_over_A = -9999.;
53  _min_modhitdens = 9999.;
54 
55  //Take the smallest and largest opening angles of clusters in this
56  //permutation and average them
57  double min_OA = 99999.;
58  double max_OA = -99999.;
59  double min_OA_wt = 99999.;
60  double max_OA_wt = -99999.;
61  for (auto const& c : clusters) {
62  // PrintClusterInfo(*c);
63  double this_OA = c->GetParams().opening_angle;
64  if (this_OA > max_OA) max_OA = this_OA;
65  if (this_OA < min_OA) min_OA = this_OA;
66  double this_OA_wt = c->GetParams().opening_angle_charge_wgt;
67  if (this_OA_wt > max_OA) max_OA_wt = this_OA_wt;
68  if (this_OA_wt < min_OA) min_OA_wt = this_OA_wt;
69  double this_trackness = c->GetParams().trackness;
70  if (this_trackness > _max_trackness) _max_trackness = this_trackness;
71  double this_L_over_W = c->GetParams().length / c->GetParams().width;
72  if (this_L_over_W > _max_len_over_width) _max_len_over_width = this_L_over_W;
73  double this_OA_over_L = this_OA / c->GetParams().length;
74  if (this_OA_over_L < _min_oa_over_len) _min_oa_over_len = this_OA_over_L;
75  double this_poly_perim_over_A =
76  c->GetParams().PolyObject.Perimeter() / c->GetParams().PolyObject.Area();
77  if (this_poly_perim_over_A > _max_poly_perim_over_A)
78  _max_poly_perim_over_A = this_poly_perim_over_A;
79  double this_modhitdens = c->GetParams().modified_hit_density;
80  if (this_modhitdens < _min_modhitdens) _min_modhitdens = this_modhitdens;
81  }
82 
83  _o_ang_avg = (min_OA + max_OA) / 2;
84  _o_ang_rms = pow((pow(min_OA, 2) + pow(max_OA, 2)) / 2, 0.5);
85  _o_ang_wt_avg = (min_OA_wt + max_OA_wt) / 2;
86  _o_ang_wt_rms = pow((pow(min_OA_wt, 2) + pow(max_OA_wt, 2)) / 2, 0.5);
87 
88  _ana_tree->Fill();
89 
90  bool accept_match = true;
91  //Reject match if it is very track-like
92  if (_min_oa_over_len < 0.0007) accept_match = false;
93  if (_o_ang_avg * _o_ang_rms < 0.01) accept_match = false;
94  if (_max_len_over_width > 20) accept_match = false;
95 
96  return accept_match ? 1 : -1;
97  }
98 
99  //------------------------------
100  void
102  //------------------------------
103  {}
104 
105  void
107  {
108  std::cout << " This cluster's info is as follows:" << std::endl;
109  std::cout << " Opening Angle: " << c.GetParams().opening_angle << std::endl;
110  // std::cout<<" Opening Angle Charge Weight: "<<c.GetParams().opening_angle_charge_wgt<<std::endl;
111  }
112 }
This algo only matches clusters if they are not track-like. This is implemented in an algo because it...
float Float(util::GeometryUtilities const &, const std::vector< const cluster::ClusterParamsAlg * > &clusters) override
const cluster_params & GetParams() const
void Reset() override
Function to reset the algorithm instance, called together with manager&#39;s Reset()
void PrintClusterInfo(const cluster::ClusterParamsAlg &c)
double opening_angle
Width of angular distubtion wrt vertx.
Definition: ClusterParams.h:40
CFAlgoShowerCompat()
Default constructor.
BEGIN_PROLOG could also be cout