LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
PIDAAlg.cxx
Go to the documentation of this file.
1 
10 #include <iostream>
11 #include <sstream>
12 #include <cmath>
13 
14 #include "PIDAAlg.h"
15 
17  fExponentConstant = p.get<float>("ExponentConstant",0.42);
18  fMinResRange = p.get<float>("MinResRange",0);
19  fMaxResRange = p.get<float>("MaxResRange",30);
20  fMaxPIDAValue = p.get<float>("MaxPIDAValue",50);
21 
22  fKDEEvalMaxSigma = p.get<float>("KDEEvalMaxSigma",3);
23  fKDEEvalStepSize = p.get<float>("KDEEvalStepSize",0.01);
24  fKDEBandwidths = p.get< std::vector<float> >("KDEBandwidths");
25 
26  fnormalDist = util::NormalDistribution(fKDEEvalMaxSigma,fKDEEvalStepSize);
27 
28  fPIDAHistNbins = p.get<unsigned int>("PIDAHistNbins",100);
29  fPIDAHistMin = p.get<float>("PIDAHistMin",0.0);
30  fPIDAHistMax = p.get<float>("PIDAHistMax",50.0);
31 
33 }
34 
40  fpida_values.clear();
41  fpida_errors.clear();
42 
43  fpida_kde_mp = std::vector<float>(fKDEBandwidths.size(),fPIDA_BOGUS);
44  fpida_kde_fwhm = std::vector<float>(fKDEBandwidths.size(),fPIDA_BOGUS);
45  fpida_kde_b = std::vector<float>(fKDEBandwidths.size(),fPIDA_BOGUS);
46  fkde_distribution = std::vector< std::vector<float> >(fKDEBandwidths.size());
47  fkde_dist_min = std::vector<float>(fKDEBandwidths.size(),fPIDA_BOGUS);
48  fkde_dist_max = std::vector<float>(fKDEBandwidths.size(),fPIDA_BOGUS);
49 }
50 
51 void pid::PIDAAlg::SetPIDATree(TTree *tree, TH1F* hist_vals, std::vector<TH1F*> hist_kde){
52 
53  if(hist_kde.size()>MAX_BANDWIDTHS)
54  throw "Error: input histograms larger than max allowed bandwidths.";
55  if(hist_kde.size()!=fKDEBandwidths.size())
56  throw "Error: input histograms do not have same size as bandwidths.";
57 
58  fPIDATree = tree;
59 
60  hPIDAvalues = hist_vals;
61  hPIDAvalues->SetNameTitle("hPIDAvalues","PIDA Distribution");
63 
64  for(size_t i_hist=0; i_hist<hist_kde.size(); i_hist++){
65  hPIDAKDE[i_hist] = hist_kde[i_hist];
66 
67  std::stringstream hname,htitle;
68  hname << "hPIDAKDE_" << i_hist;
69  htitle << "PIDA KDE-smoothed Distribution, Bandwidth=" << fKDEBandwidths.at(i_hist);
70 
71  hPIDAKDE[i_hist]->SetNameTitle(hname.str().c_str(),htitle.str().c_str());
73  }
74 
76  fPIDATree->Branch("hpida_vals","TH1F",hPIDAvalues);
77  fPIDATree->Branch("n_bandwidths",&(fPIDAProperties.n_bandwidths),"n_bandwidths/i");
78  fPIDATree->Branch("kde_bandwidth",fPIDAProperties.kde_bandwidth,"kde_bandwidth[n_bandwidths]/F");
79  fPIDATree->Branch("kde_mp",fPIDAProperties.kde_mp,"kde_mp[n_bandwidths]/F");
80  fPIDATree->Branch("kde_fwhm",fPIDAProperties.kde_fwhm,"kde_fwhm[n_bandwidths]/F");
81  for(size_t i_hist=0; i_hist<hist_kde.size(); i_hist++){
82  std::stringstream bname;
83  bname << "hpida_kde_" << i_hist;
84  fPIDATree->Branch(bname.str().c_str(),"TH1F",hPIDAKDE[i_hist]);
85  }
86 
87 }
88 
92 
93  return fpida_mean;
94 }
95 
99 
100  return fpida_sigma;
101 }
102 
103 float pid::PIDAAlg::getPIDAKDEMostProbable(const size_t i_b){
104  if(fpida_kde_mp[i_b]==fPIDA_BOGUS)
105  createKDE(i_b);
106 
107  return fpida_kde_mp[i_b];
108 }
109 
111  if(fpida_kde_fwhm[i_b]==fPIDA_BOGUS)
112  createKDE(i_b);
113 
114  return fpida_kde_fwhm[i_b];
115 }
116 
117 const std::vector<float>& pid::PIDAAlg::getPIDAValues(){
118  return fpida_values;
119 }
120 
121 const std::vector<float>& pid::PIDAAlg::getPIDAErrors(){
122  return fpida_errors;
123 }
124 
126  std::vector<float> const& resRange = calo.ResidualRange();
127  std::vector<float> const& dEdx = calo.dEdx();
128  RunPIDAAlg(resRange,dEdx);
129 }
130 
132  float& mean,
133  float& sigma)
134 {
135  RunPIDAAlg(calo);
136  mean = getPIDAMean();
137  sigma = getPIDASigma();
138 }
139 
140 void pid::PIDAAlg::RunPIDAAlg(std::vector<float> const& resRange,
141  std::vector<float> const& dEdx){
142 
144 
145  fpida_values.reserve( resRange.size() );
146  fpida_errors.reserve( resRange.size() );
147 
148  std::map<double,double> range_dEdx_map;
149 
150  for(size_t i_r=0; i_r<resRange.size(); i_r++){
151  if(resRange[i_r]>fMaxResRange || resRange[i_r]<fMinResRange) continue;
152 
153  range_dEdx_map[ resRange[i_r] ] = dEdx[i_r];
154 
155  float val = dEdx[i_r]*std::pow(resRange[i_r],fExponentConstant);
156  if(val < fMaxPIDAValue){
157  fpida_values.push_back(val);
158  //fpida_errors.push_back(0);
159  }
160 
161  }
162 
163  calculatePIDAIntegral(range_dEdx_map);
164 
165  if(fpida_values.size()==0)
166  fpida_values.push_back(-99);
167 
168 }
169 
170 void pid::PIDAAlg::FillPIDATree(unsigned int run,
171  unsigned int event,
172  unsigned int calo_index,
173  anab::Calorimetry const& calo){
174  RunPIDAAlg(calo);
175  FillPIDAProperties(run,event,calo_index,calo);
176 }
177 
179 
180  if(fpida_values.size()==0)
181  throw "pid::PIDAAlg --- PIDA Values not filled!";
182 
183  fpida_mean = 0;
184  for(auto const& val : fpida_values)
185  fpida_mean += val;
186  fpida_mean /= fpida_values.size();
187 
188 }
189 
191 
194 
195  fpida_sigma = 0;
196  for(auto const& val : fpida_values)
197  fpida_sigma += (fpida_mean-val)*(fpida_mean-val);
198 
199  fpida_sigma = std::sqrt(fpida_sigma)/fpida_values.size();
200 }
201 
202 void pid::PIDAAlg::calculatePIDAIntegral(std::map<double,double> const& range_dEdx_map){
203 
204  if(range_dEdx_map.size()<2) return;
205 
207 
208  for( std::map<double,double>::const_iterator map_iter = range_dEdx_map.begin();
209  map_iter != std::prev(range_dEdx_map.end());
210  map_iter++)
211  {
212  double range_width = std::next(map_iter)->first - map_iter->first;
213  fpida_integral_dedx += range_width*( std::next(map_iter)->second + 0.5*(map_iter->second-std::next(map_iter)->second));
214  }
215 
217  std::pow( (std::prev(range_dEdx_map.end())->first - range_dEdx_map.begin()->first),(fExponentConstant-1) );
218 }
219 
220 void pid::PIDAAlg::createKDE(const size_t i_b){
221 
222  if(fpida_values.size()==0)
223  throw "pid::PIDAAlg --- PIDA Values not filled!";
224 
225  //if( fkde_distribution[i_b].size()!=0 ) return;
226 
227  if(fKDEBandwidths[i_b]<=0) {
229  float bandwidth = fpida_sigma*1.06*std::pow((float)(fpida_values.size()),-0.2);
230  fpida_errors = std::vector<float>(fpida_values.size(),bandwidth);
231  fpida_kde_b[i_b] = bandwidth;
232  }
233  else{
234  fpida_errors = std::vector<float>(fpida_values.size(),fKDEBandwidths[i_b]);
235  fpida_kde_b[i_b] = fKDEBandwidths[i_b];
236  }
237 
238  const auto min_pida_iterator = std::min_element(fpida_values.begin(),fpida_values.end());
239  const size_t min_pida_location = std::distance(fpida_values.begin(),min_pida_iterator);
240  fkde_dist_min[i_b] = fpida_values[min_pida_location] - fKDEEvalMaxSigma*fpida_errors[min_pida_location];
241 
242  const auto max_pida_iterator = std::max_element(fpida_values.begin(),fpida_values.end());
243  const size_t max_pida_location = std::distance(fpida_values.begin(),max_pida_iterator);
244  fkde_dist_max[i_b] = fpida_values[max_pida_location] + fKDEEvalMaxSigma*fpida_errors[max_pida_location];
245 
246  //make the kde distribution, and get the max value
247  const size_t kde_dist_size = (size_t)( (fkde_dist_max[i_b] - fkde_dist_min[i_b])/fKDEEvalStepSize ) + 1;
248  fkde_distribution[i_b].resize(kde_dist_size);
249  float kde_max=0;
250  size_t step_max=0;
251  for(size_t i_step=0; i_step<kde_dist_size; i_step++){
252  float pida_val = fkde_dist_min[i_b] + i_step*fKDEEvalStepSize;
253  fkde_distribution[i_b][i_step]=0;
254 
255  for(size_t i_pida=0; i_pida<fpida_values.size(); i_pida++)
256  fkde_distribution[i_b][i_step] += fnormalDist.getValue((fpida_values[i_pida]-pida_val)/fpida_errors[i_pida])/fpida_errors[i_pida];
257 
258  if(fkde_distribution[i_b][i_step]>kde_max){
259  kde_max = fkde_distribution[i_b][i_step];
260  step_max = i_step;
261  fpida_kde_mp[i_b] = pida_val;
262  }
263  }
264 
265  //now get fwhm
266  float half_max = 0.5*fpida_kde_mp[i_b];
267  float low_width=0;
268  for(size_t i_step=step_max; i_step>0; i_step--){
269  if(fkde_distribution[i_b][i_step] < half_max) break;
270  low_width += fKDEEvalStepSize;
271  }
272  float high_width=0;
273  for(size_t i_step=step_max; i_step<kde_dist_size; i_step++){
274  if(fkde_distribution[i_b][i_step] < half_max) break;
275  high_width += fKDEEvalStepSize;
276  }
277  fpida_kde_fwhm[i_b] = low_width+high_width;
278 
279 }
280 
282  for(size_t i_b=0; i_b < fKDEBandwidths.size(); i_b++)
283  createKDE(i_b);
284 }
285 
286 void pid::PIDAAlg::FillPIDAProperties(unsigned int run,
287  unsigned int event,
288  unsigned int calo_index,
289  anab::Calorimetry const& calo){
290 
291  fPIDAProperties.run = run;
292  fPIDAProperties.event = event;
293  fPIDAProperties.calo_index = calo_index;
297 
299  for(size_t i_b=0; i_b<fPIDAProperties.n_bandwidths; i_b++)
300 
305 
308 
309  createKDEs();
310  for(size_t i_b=0; i_b<fPIDAProperties.n_bandwidths; i_b++){
311  fPIDAProperties.kde_mp[i_b] = fpida_kde_mp[i_b];
314  }
315 
316  hPIDAvalues->Reset();
317  for(auto const& val: fpida_values)
318  hPIDAvalues->Fill(val);
319 
320  for(size_t i_b=0; i_b<fPIDAProperties.n_bandwidths; i_b++){
321  hPIDAKDE[i_b]->Reset();
322  for(size_t i_step=0; i_step<fkde_distribution[i_b].size(); i_step++)
323  hPIDAKDE[i_b]->AddBinContent(hPIDAKDE[i_b]->FindBin(i_step*fKDEEvalStepSize+fkde_dist_min[i_b]),
324  fkde_distribution[i_b][i_step]);
325  }
326 
327  fPIDATree->Fill();
328 }
329 
331  for(size_t i_pida=0; i_pida<fpida_values.size(); i_pida++)
332  std::cout << "\tPIDA --- " << i_pida << "\t" << fpida_values[i_pida] << std::endl;
333 }
334 
335 util::NormalDistribution::NormalDistribution(float max_sigma, float step_size){
336 
337  if(step_size==0)
338  throw "util::NormalDistribution --- Cannot have zero step size!";
339 
340  const size_t vector_size = (size_t)(max_sigma / step_size);
341  fValues.resize(vector_size);
342 
343  const float AMPLITUDE = 1. / std::sqrt(2*M_PI);
344 
345  float integral=0;
346  for(size_t i_step=0; i_step<vector_size; i_step++){
347  float diff = i_step*step_size;
348  fValues[i_step] = AMPLITUDE * std::exp(-0.5*diff*diff);
349  integral+= fValues[i_step];
350  }
351 
352  for(size_t i_step=0; i_step<vector_size; i_step++)
353  fValues[i_step] /= (integral*2);
354 
355  fStepSize = step_size;
356  fMaxSigma = fStepSize * vector_size;
357 
358 }
359 
361 
362  x = std::abs(x);
363  if(x > fMaxSigma) return 0;
364 
365  size_t bin_low = x / fStepSize;
366  float remainder = (x - (bin_low*fStepSize)) / fStepSize;
367 
368  return fValues[bin_low]*(1-remainder) + remainder*fValues[bin_low+1];
369 
370 }
Float_t x
Definition: compare.C:6
std::vector< float > fpida_errors
Definition: PIDAAlg.h:86
TH1F * hPIDAKDE[MAX_BANDWIDTHS]
Definition: PIDAAlg.h:115
void SetPIDATree(TTree *, TH1F *, std::vector< TH1F * >)
Definition: PIDAAlg.cxx:51
TH1F * hPIDAvalues
Definition: PIDAAlg.h:114
std::string leaf_structure
Definition: PIDAAlg.h:138
float fPIDAHistMin
Definition: PIDAAlg.h:117
TTree * fPIDATree
Definition: PIDAAlg.h:113
float kde_bandwidth[MAX_BANDWIDTHS]
Definition: PIDAAlg.h:134
float fpida_integral_pida
Definition: PIDAAlg.h:90
void FillPIDATree(unsigned int, unsigned int, unsigned int, anab::Calorimetry const &)
Definition: PIDAAlg.cxx:170
util::NormalDistribution fnormalDist
Definition: PIDAAlg.h:111
float fMaxPIDAValue
Definition: PIDAAlg.h:80
std::vector< float > fkde_dist_min
Definition: PIDAAlg.h:108
const geo::PlaneID & PlaneID() const
Definition: Calorimetry.h:119
float getPIDAMean()
Definition: PIDAAlg.cxx:89
std::vector< float > fpida_kde_mp
Definition: PIDAAlg.h:102
float fMinResRange
Definition: PIDAAlg.h:78
void FillPIDAProperties(unsigned int, unsigned int, unsigned int, anab::Calorimetry const &)
Definition: PIDAAlg.cxx:286
void createKDE(const size_t)
Definition: PIDAAlg.cxx:220
float fMaxResRange
Definition: PIDAAlg.h:79
void createKDEs()
Definition: PIDAAlg.cxx:281
unsigned int n_bandwidths
Definition: PIDAAlg.h:133
const std::vector< float > & ResidualRange() const
Definition: Calorimetry.h:106
float fKDEEvalMaxSigma
Definition: PIDAAlg.h:81
const std::vector< float > & getPIDAValues()
Definition: PIDAAlg.cxx:117
float fKDEEvalStepSize
Definition: PIDAAlg.h:82
std::vector< float > fkde_dist_max
Definition: PIDAAlg.h:109
std::vector< float > fpida_kde_fwhm
Definition: PIDAAlg.h:103
unsigned int fPIDAHistNbins
Definition: PIDAAlg.h:116
float kde_fwhm[MAX_BANDWIDTHS]
Definition: PIDAAlg.h:136
float fPIDAHistMax
Definition: PIDAAlg.h:118
intermediate_table::const_iterator const_iterator
void reconfigure(fhicl::ParameterSet const &p)
Definition: PIDAAlg.cxx:16
float getPIDAKDEMostProbable(const size_t)
Definition: PIDAAlg.cxx:103
T get(std::string const &key) const
Definition: ParameterSet.h:231
const float fPIDA_BOGUS
Definition: PIDAAlg.h:75
float getPIDASigma()
Definition: PIDAAlg.cxx:96
float getValue(float)
Definition: PIDAAlg.cxx:360
void RunPIDAAlg(std::vector< float > const &, std::vector< float > const &)
Definition: PIDAAlg.cxx:140
const unsigned int MAX_BANDWIDTHS
Definition: PIDAAlg.h:41
std::vector< float > fpida_values
Definition: PIDAAlg.h:85
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
void ClearInternalData()
Definition: PIDAAlg.cxx:35
void PrintPIDAValues()
Definition: PIDAAlg.cxx:330
const std::vector< float > & dEdx() const
Definition: Calorimetry.h:104
float getPIDAKDEFullWidthHalfMax(const size_t)
Definition: PIDAAlg.cxx:110
float fpida_sigma
Definition: PIDAAlg.h:88
std::vector< float > fpida_kde_b
Definition: PIDAAlg.h:104
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
Definition: UtilFunc.cxx:15
void calculatePIDAIntegral(std::map< double, double > const &)
Definition: PIDAAlg.cxx:202
float fExponentConstant
Definition: PIDAAlg.h:77
void calculatePIDAMean()
Definition: PIDAAlg.cxx:178
const float & KineticEnergy() const
Definition: Calorimetry.h:108
std::vector< float > fKDEBandwidths
Definition: PIDAAlg.h:83
std::vector< std::vector< float > > fkde_distribution
Definition: PIDAAlg.h:107
const std::vector< float > & getPIDAErrors()
Definition: PIDAAlg.cxx:121
float fpida_mean
Definition: PIDAAlg.h:87
PIDAProperties_t fPIDAProperties
Definition: PIDAAlg.h:143
float fpida_integral_dedx
Definition: PIDAAlg.h:89
const float & Range() const
Definition: Calorimetry.h:109
void calculatePIDASigma()
Definition: PIDAAlg.cxx:190
float kde_mp[MAX_BANDWIDTHS]
Definition: PIDAAlg.h:135
calorimetry
Event finding and building.