LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
cmtool::CFAlgoShowerCompat Class Reference

#include "CFAlgoShowerCompat.h"

Inheritance diagram for cmtool::CFAlgoShowerCompat:
cmtool::CFloatAlgoBase cmtool::CMAlgoBase

Public Member Functions

 CFAlgoShowerCompat ()
 Default constructor. More...
 
virtual ~CFAlgoShowerCompat ()
 Default destructor. More...
 
virtual float Float (const std::vector< const cluster::ClusterParamsAlg * > &clusters)
 
virtual void Report ()
 
virtual void Reset ()
 Function to reset the algorithm instance, called together with manager's Reset() More...
 
void PrintClusterInfo (const cluster::ClusterParamsAlg &c)
 
void WriteHaxFile ()
 
virtual void EventBegin (const std::vector< cluster::ClusterParamsAlg > &clusters)
 
virtual void EventEnd ()
 
virtual void IterationBegin (const std::vector< cluster::ClusterParamsAlg > &clusters)
 
virtual void IterationEnd ()
 
void SetAnaFile (TFile *fout)
 Setter function for an output plot TFile pointer. More...
 
virtual void SetVerbose (bool doit=true)
 Setter function for verbosity. More...
 

Protected Attributes

TFile * _fout
 TFile pointer to an output file. More...
 
bool _verbose
 Boolean to choose verbose mode. Turned on if CMergeManager/CMatchManager's verbosity level is >= kPerMerging. More...
 

Private Attributes

TTree * _ana_tree
 
double _o_ang_avg
 
double _o_ang_rms
 
double _o_ang_wt_avg
 
double _o_ang_wt_rms
 
double _max_trackness
 
double _max_len_over_width
 
double _min_oa_over_len
 
double _max_poly_perim_over_A
 
double _min_modhitdens
 
TFile * _fout_hax
 

Detailed Description

User implementation for CFloatAlgoBase class doxygen documentation!

Definition at line 31 of file CFAlgoShowerCompat.h.

Constructor & Destructor Documentation

cmtool::CFAlgoShowerCompat::CFAlgoShowerCompat ( )

Default constructor.

Definition at line 10 of file CFAlgoShowerCompat.cxx.

References _ana_tree, _fout_hax, _max_len_over_width, _max_poly_perim_over_A, _max_trackness, _min_modhitdens, _min_oa_over_len, _o_ang_avg, _o_ang_rms, _o_ang_wt_avg, and _o_ang_wt_rms.

10  : CFloatAlgoBase()
11  //-------------------------------------------------------
12  {
13 
14  _fout_hax = 0;
15  _ana_tree = 0;
16 
17  if(!_fout_hax)
18  _fout_hax = new TFile("fout_hax.root","RECREATE");
19 
20  if(!_ana_tree){
21  _ana_tree = new TTree("ana_tree","ana_tree");
22  _ana_tree->Branch("o_ang_avg",&_o_ang_avg,"o_ang_avg/D");
23  _ana_tree->Branch("o_ang_rms",&_o_ang_rms,"o_ang_rms/D");
24  _ana_tree->Branch("o_ang_wt_avg",&_o_ang_wt_avg,"o_ang_wt_avg/D");
25  _ana_tree->Branch("o_ang_wt_rms",&_o_ang_wt_rms,"o_ang_wt_rms/D");
26  _ana_tree->Branch("max_trackness",&_max_trackness,"max_trackness/D");
27  _ana_tree->Branch("max_len_over_width",&_max_len_over_width,"max_len_over_width/D");
28  _ana_tree->Branch("min_oa_over_len",&_min_oa_over_len,"min_oa_over_len/D");
29  _ana_tree->Branch("max_poly_perim_over_A",&_max_poly_perim_over_A,"max_poly_perim_over_A/D");
30  _ana_tree->Branch("min_modhitdens",&_min_modhitdens,"min_modhitdens/D");
31  }
32 
33  }
CFloatAlgoBase()
Default constructor.
virtual cmtool::CFAlgoShowerCompat::~CFAlgoShowerCompat ( )
inlinevirtual

Default destructor.

Definition at line 39 of file CFAlgoShowerCompat.h.

References Float(), PrintClusterInfo(), Report(), and Reset().

39 {};

Member Function Documentation

virtual void cmtool::CMAlgoBase::EventBegin ( const std::vector< cluster::ClusterParamsAlg > &  clusters)
inlinevirtualinherited

Optional function: called at the beginning of 1st iteration. This is called per event.

Reimplemented in cmtool::CFAlgoArray, cmtool::CPAlgoArray, cmtool::CBAlgoArray, and cmtool::CBAlgoPolyShortestDist.

Definition at line 45 of file CMAlgoBase.h.

Referenced by cmtool::CMergeManager::EventBegin().

46  { if(clusters.size()) return; }
virtual void cmtool::CMAlgoBase::EventEnd ( )
inlinevirtualinherited

Optional function: called at the end of event ... after the last merging iteration is over.

Reimplemented in cmtool::CFAlgoArray, cmtool::CPAlgoArray, and cmtool::CBAlgoArray.

Definition at line 51 of file CMAlgoBase.h.

Referenced by cmtool::CMatchManager::EventEnd(), and cmtool::CMergeManager::EventEnd().

52  {return;}
float cmtool::CFAlgoShowerCompat::Float ( const std::vector< const cluster::ClusterParamsAlg * > &  clusters)
virtual

Core function: given a set of CPANs, return a float which indicates the compatibility the cluster combination.

Reimplemented from cmtool::CFloatAlgoBase.

Definition at line 43 of file CFAlgoShowerCompat.cxx.

References _ana_tree, _max_len_over_width, _max_poly_perim_over_A, _max_trackness, _min_modhitdens, _min_oa_over_len, _o_ang_avg, _o_ang_rms, _o_ang_wt_avg, and _o_ang_wt_rms.

Referenced by ~CFAlgoShowerCompat().

45  {
46  _o_ang_avg = 0;
47  _o_ang_rms = 0;
48  _o_ang_wt_avg = 0;
49  _o_ang_wt_rms = 0;
50  _max_trackness = -9999.;
51  _max_len_over_width = -9999.;
52  _min_oa_over_len = 9999.;
53  _max_poly_perim_over_A = -9999.;
54  _min_modhitdens = 9999.;
55 
56  //Take the smallest and largest opening angles of clusters in this
57  //permutation and average them
58  double min_OA = 99999.;
59  double max_OA = -99999.;
60  double min_OA_wt = 99999.;
61  double max_OA_wt = -99999.;
62  for(auto const& c : clusters){
63  // PrintClusterInfo(*c);
64  double this_OA = c->GetParams().opening_angle;
65  if(this_OA > max_OA) max_OA = this_OA;
66  if(this_OA < min_OA) min_OA = this_OA;
67  double this_OA_wt = c->GetParams().opening_angle_charge_wgt;
68  if(this_OA_wt > max_OA) max_OA_wt = this_OA_wt;
69  if(this_OA_wt < min_OA) min_OA_wt = this_OA_wt;
70  double this_trackness = c->GetParams().trackness;
71  if(this_trackness > _max_trackness) _max_trackness = this_trackness;
72  double this_L_over_W = c->GetParams().length / c->GetParams().width;
73  if(this_L_over_W > _max_len_over_width) _max_len_over_width = this_L_over_W;
74  double this_OA_over_L = this_OA/c->GetParams().length;
75  if(this_OA_over_L < _min_oa_over_len) _min_oa_over_len = this_OA_over_L;
76  double this_poly_perim_over_A = c->GetParams().PolyObject.Perimeter()/c->GetParams().PolyObject.Area();
77  if(this_poly_perim_over_A > _max_poly_perim_over_A) _max_poly_perim_over_A=this_poly_perim_over_A;
78  double this_modhitdens = c->GetParams().modified_hit_density;
79  if(this_modhitdens < _min_modhitdens) _min_modhitdens = this_modhitdens;
80  }
81 
82  _o_ang_avg = (min_OA + max_OA)/2;
83  _o_ang_rms = pow( (pow(min_OA,2)+pow(max_OA,2))/2 , 0.5);
84  _o_ang_wt_avg = (min_OA_wt + max_OA_wt)/2;
85  _o_ang_wt_rms = pow( (pow(min_OA_wt,2)+pow(max_OA_wt,2))/2 , 0.5);
86 
87  _ana_tree->Fill();
88 
89  bool accept_match = true;
90  //Reject match if it is very track-like
91  if(_min_oa_over_len < 0.0007) accept_match = false;
92  if(_o_ang_avg*_o_ang_rms < 0.01) accept_match = false;
93  if(_max_len_over_width > 20) accept_match = false;
94 
95  return accept_match ? 1 : -1;
96  }
virtual void cmtool::CMAlgoBase::IterationBegin ( const std::vector< cluster::ClusterParamsAlg > &  clusters)
inlinevirtualinherited

Optional function: called at the beggining of each iteration over all pairs of clusters. This provides all clusters' information in case the algorithm need them. Note this is called per iteration which may be more than once per event.

Reimplemented in cmtool::CFAlgoArray, cmtool::CPAlgoArray, and cmtool::CBAlgoArray.

Definition at line 59 of file CMAlgoBase.h.

Referenced by cmtool::CMatchManager::EventBegin(), cmtool::CMatchManager::IterationBegin(), and cmtool::CMergeManager::IterationBegin().

60  { if(clusters.size()) return;}
virtual void cmtool::CMAlgoBase::IterationEnd ( )
inlinevirtualinherited

Optional function: called at the end of each iteration over all pairs of clusters.

Reimplemented in cmtool::CFAlgoArray, cmtool::CPAlgoArray, and cmtool::CBAlgoArray.

Definition at line 65 of file CMAlgoBase.h.

Referenced by cmtool::CMatchManager::IterationEnd(), and cmtool::CMergeManager::IterationEnd().

66  {return; }
void cmtool::CFAlgoShowerCompat::PrintClusterInfo ( const cluster::ClusterParamsAlg c)

Definition at line 105 of file CFAlgoShowerCompat.cxx.

References cluster::ClusterParamsAlg::GetParams(), and cluster::cluster_params::opening_angle.

Referenced by ~CFAlgoShowerCompat().

105  {
106  std::cout<<" This cluster's info is as follows:"<<std::endl;
107  std::cout<<" Opening Angle: "<<c.GetParams().opening_angle<<std::endl;
108  // std::cout<<" Opening Angle Charge Weight: "<<c.GetParams().opening_angle_charge_wgt<<std::endl;
109 
110  }
const cluster_params & GetParams() const
double opening_angle
Width of angular distubtion wrt vertx.
Definition: ClusterParams.h:41
void cmtool::CFAlgoShowerCompat::Report ( )
virtual

Optional function: called after each iterative approach if a manager class is run with verbosity level <= kPerIteration. Maybe useful for debugging.

Reimplemented from cmtool::CMAlgoBase.

Definition at line 99 of file CFAlgoShowerCompat.cxx.

Referenced by ~CFAlgoShowerCompat().

101  {
102 
103  }
void cmtool::CFAlgoShowerCompat::Reset ( )
virtual

Function to reset the algorithm instance, called together with manager's Reset()

Reimplemented from cmtool::CMAlgoBase.

Definition at line 36 of file CFAlgoShowerCompat.cxx.

Referenced by ~CFAlgoShowerCompat().

38  {
39 
40  }
void cmtool::CMAlgoBase::SetAnaFile ( TFile *  fout)
inlineinherited

Setter function for an output plot TFile pointer.

Definition at line 77 of file CMAlgoBase.h.

References cmtool::CMAlgoBase::_fout.

Referenced by cmtool::CMergeManager::EventBegin().

77 { _fout = fout; }
TFile * _fout
TFile pointer to an output file.
Definition: CMAlgoBase.h:85
virtual void cmtool::CMAlgoBase::SetVerbose ( bool  doit = true)
inlinevirtualinherited
void cmtool::CFAlgoShowerCompat::WriteHaxFile ( )
inline

Definition at line 64 of file CFAlgoShowerCompat.h.

References _ana_tree, and _fout_hax.

65  {
66  _fout_hax->cd();
67  _ana_tree->Write();
68  _fout_hax->Close();
69  };

Member Data Documentation

TTree* cmtool::CFAlgoShowerCompat::_ana_tree
private

Optional function: called at the beginning of 1st iteration. This is called per event. Optional function: called at the end of event ... after the last merging iteration is over. Optional function: called at the beggining of each iterative loop. This provides all clusters' information in case the algorithm need them. Note this is called per iteration which may be more than once per event. Optional function: called at the end of each iterative loop.

Definition at line 69 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), Float(), and WriteHaxFile().

TFile* cmtool::CMAlgoBase::_fout
protectedinherited

TFile pointer to an output file.

Definition at line 85 of file CMAlgoBase.h.

Referenced by cmtool::CMAlgoBase::CMAlgoBase(), and cmtool::CMAlgoBase::SetAnaFile().

TFile* cmtool::CFAlgoShowerCompat::_fout_hax
private

Definition at line 107 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and WriteHaxFile().

double cmtool::CFAlgoShowerCompat::_max_len_over_width
private

Definition at line 102 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_max_poly_perim_over_A
private

Definition at line 104 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_max_trackness
private

Definition at line 101 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_min_modhitdens
private

Definition at line 105 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_min_oa_over_len
private

Definition at line 103 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_o_ang_avg
private

Definition at line 97 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_o_ang_rms
private

Definition at line 98 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_o_ang_wt_avg
private

Definition at line 99 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().

double cmtool::CFAlgoShowerCompat::_o_ang_wt_rms
private

Definition at line 100 of file CFAlgoShowerCompat.h.

Referenced by CFAlgoShowerCompat(), and Float().


The documentation for this class was generated from the following files: