LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
CreateHybridLibrary_module.cc
Go to the documentation of this file.
1 
3 // Chris Backhouse, UCL, Nov 2017
5 
14 
15 #include <iostream>
16 
18 
19 #include "TF1.h"
20 #include "TFile.h"
21 #include "TGraph.h"
22 #include "TGraphErrors.h"
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TPad.h"
26 #include "TTree.h"
27 #include "TVectorD.h"
28 #include "TCanvas.h"
29 #include "TStyle.h"
30 #include "TMultiGraph.h"
31 
32 #define PI 3.14159265
33 
34 namespace phot
35 {
37  {
38  public:
39  explicit CreateHybridLibrary(const fhicl::ParameterSet& p);
40 
41  // Plugins should not be copied or assigned.
46 
47  void analyze(const art::Event& e) override;
48  };
49 
50  //--------------------------------------------------------------------
52  : EDAnalyzer(p)
53  {
54 
56 
58  sim::PhotonVoxelDef voxdef = pvs->GetVoxelDef();
59 
60  TFile* fout_full = new TFile("full.root", "RECREATE");
61  TFile* fout_fit = new TFile("fit.root", "RECREATE");
62 
63  std::cout << voxdef.GetNVoxels() << " voxels for each of " << geom->NOpDets() << " OpDets" << std::endl;
64  std::cout << std::endl;
65 
66  //EP = Exception points: the parameterization is not a good description of the visibility and the value of the Photon Library is kept.
67  long totExceptions = 0;
68  long totPts = 0;
69 
70  for(unsigned int opdetIdx =0; opdetIdx < geom->NOpDets(); ++opdetIdx){
71  std::cout << opdetIdx << " / " << geom->NOpDets() << std::endl;
72 
73  TDirectory* d_full = fout_full->mkdir(TString::Format("opdet_%d", opdetIdx).Data());
74  TDirectory* d_fit = fout_fit->mkdir(TString::Format("opdet_%d", opdetIdx).Data());
75 
76  d_full->cd();
77  TTree* tr_full = new TTree("tr", "tr");
78  int vox, taken;
79  float dist, vis, psi, theta, xpos;
80  tr_full->Branch("vox", &vox);
81  tr_full->Branch("dist", &dist);
82  tr_full->Branch("vis", &vis);
83  tr_full->Branch("taken", &taken);
84  tr_full->Branch("psi", &psi); //not needed to parameterize the visibilities, useful for tests in SP
85  tr_full->Branch("theta", &theta); //not needed to parameterize the visibilities, useful for tests in SP
86  tr_full->Branch("xpos", &xpos); //not needed to parameterize the visibilities, useful for tests in SP
87 
88  const geo::OpDetGeo& opdet = geom->OpDetGeoFromOpDet(opdetIdx);
89  double xyzopdet[3];
90  opdet.GetCenter(xyzopdet);
91  //std::cout << "OpDet position " << opdetIdx << " x = " <<xyzopdet[0] << " y = " <<xyzopdet[1] << " z = " <<xyzopdet[2]<< std::endl;
92  const TVector3 opdetvec(xyzopdet);
93 
94  struct Visibility{
95  Visibility(int vx, int t, float d, float v, float p, float th, float xp) : vox(vx), taken(t), dist(d), vis(v), psi(p), theta(th), xpos(xp) {}
96  int vox;
97  int taken;
98  float dist;
99  float vis;
100  float psi;
101  float theta;
102  float xpos;
103  };
104  TCanvas *c1=new TCanvas("c1", "c1");
105  //c1->SetCanvasSize(1500, 1500);
106  c1->SetWindowSize(600, 600);
107  //c1->Divide(1,2);
108  TGraph g;
109  TGraph g2;
110 
111  std::vector<Visibility> viss;
112  viss.reserve(voxdef.GetNVoxels());
113 
114  for(int voxIdx = 0; voxIdx < voxdef.GetNVoxels(); ++voxIdx){
115 
116  const TVector3 voxvec = voxdef.GetPhotonVoxel(voxIdx).GetCenter();
117  const double xyzvox[] = {voxvec.X(), voxvec.Y(), voxvec.Z()};
118  const double fc_y = 600; //624cm is below the center of the first voxel outside the Field Cage
119  const double fc_z = 1394;
120  const double fc_x = 350;
121  taken = 0;
122  //DP does not need variable "taken" because all voxels are inside the Field Cage for the Photon Library created in LightSim.
123  //DP taken = 1;
124  dist = opdet.DistanceToPoint(xyzvox);
125  vis = pvs->GetVisibility(xyzvox, opdetIdx);
126  // all voxels outside the Field Cage would be assigned these values of psi and theta
127  psi = 100;
128  theta = 200;
129  xpos = xyzvox[0];
130 
131  if((xyzvox[0] - xyzopdet[0])<0){
132  psi = atan((xyzvox[1] - xyzopdet[1])/(-xyzvox[0] +xyzopdet[0]));
133  }
134 
135  if(xyzvox[0]<fc_x && xyzvox[0]>-fc_x && xyzvox[1]<fc_y && xyzvox[1]>-fc_y && xyzvox[2]> -9 && xyzvox[2]<fc_z){
136  g.SetPoint(g.GetN(), dist, vis*dist*dist);
137  taken = 1;
138  psi = atan((xyzvox[1] - xyzopdet[1])/(xyzvox[0] - xyzopdet[0]))* 180.0/PI ; // psi takes values within (-PI/2, PI/2)
139  theta = acos((xyzvox[2] - xyzopdet[2])/dist) * 180.0/PI; // theta takes values within (0 (beam direction, z), PI (-beam direction, -z))
140 
141  if((xyzvox[0] - xyzopdet[0])<0){
142  psi = atan((xyzvox[1] - xyzopdet[1])/(-xyzvox[0] +xyzopdet[0]));
143  }
144 
145  }
146  tr_full->Fill();
147  viss.emplace_back(voxIdx, taken, dist, vis, psi, theta, xpos);
148  } // end for voxIdx
149 
150  d_full->cd();
151  tr_full->Write();
152  delete tr_full;
153 
154  g.SetMarkerStyle(7);
155  c1->cd();
156  g.Draw("ap");
157  g.GetXaxis()->SetTitle("Distance (cm)");
158  g.GetYaxis()->SetTitle("Visibility #times r^{2}");
159  g.Fit("expo");
160  TF1* fit = g.GetFunction("expo");
161 
162  d_fit->cd();
163  TVectorD fitres(2);
164  fitres[0] = fit->GetParameter(0);
165  fitres[1] = fit->GetParameter(1);
166  fitres.Write("fit");
167 
168  gPad->SetLogy();
169 
170  TH1F h("", "", 200, 0, 20);
171 
172  d_fit->cd();
173  TTree* tr_fit = new TTree("tr", "tr");
174  tr_fit->Branch("vox", &vox);
175  tr_fit->Branch("dist", &dist);
176  tr_fit->Branch("vis", &vis);
177  tr_fit->Branch("taken", &taken);
178  tr_fit->Branch("psi", &psi);
179  tr_fit->Branch("theta", &theta);
180  tr_fit->Branch("xpos", &xpos);
181 
182 
183  for(const Visibility& v: viss){
184  // 2e-5 is the magic scaling factor to get back to integer photon
185  // counts. TODO this will differ for new libraries, should work out a
186  // way to communicate it or derive it.
187  const double obs = v.vis / 2e-5; //taken from the Photon Library
188  const double pred = fit->Eval(v.dist) / (v.dist*v.dist) / 2e-5; //calculated with parameterization
189 
190  //DP const double obs = v.vis / 1e-7; //magic scaling factor for DP library created in LightSim
191  //Minimal amount of detected photons is 50, bc of Landau dustribution
192  //Those voxels with detected photons < 50 were set to 0
193  //DP const double pred = fit->Eval(v.dist) / (v.dist*v.dist) / 1e-7; //calculated with parametrisation
194  //std::cout << "observed = "<<obs<<" predicted (by parametrization) = "<<pred <<std::endl;
195 
196 
197  // Log-likelihood ratio for poisson statistics
198  double chisq = 2*(pred-obs);
199  if(obs) chisq += 2*obs*log(obs/pred);
200 
201  vox = v.vox;
202  dist = v.dist;
203  vis = pred *2e-5;
204  psi = v.psi;
205  theta = v.theta;
206  xpos = v.xpos;
207  //DP vis = pred *1e-7;
208 
209 
210  if (v.taken==1){
211  h.Fill(chisq);
212  }
213 
214 
215  if(chisq > 9){ //equivalent to more than 9 chisquare = 3 sigma //maybe play around with this cutoff
216  g2.SetPoint(g2.GetN(), v.dist, v.vis*v.dist*v.dist);
217  vis = obs *2e-5;
218  //DP vis = obs *1e-7;
219  tr_fit->Fill();
220  }
221 
222  }
223 
224  d_fit->cd();
225  tr_fit->Write();
226  delete tr_fit;
227  g2.SetMarkerSize(1);
228  g2.SetLineWidth(3);
229  g2.SetMarkerStyle(7);
230  g2.SetMarkerColor(kRed);
231  gStyle->SetLabelSize(.04, "XY");
232  gStyle->SetTitleSize(.04,"XY");
233  c1->cd();
234  g2.Draw("p same");
235  c1->SaveAs(TString::Format("plots/Chris_vis_vs_dist_%d.png", opdetIdx).Data());
236 
237 
238  h.GetXaxis()->SetTitle("#chi^{2}");
239  h.GetYaxis()->SetTitle("Counts");
240  gStyle->SetLabelSize(.04, "XY");
241  gStyle->SetTitleSize(.04,"XY");
242  h.Draw("hist");
243  std::cout <<"Integral(90,-1)/Integral(all) = "<< h.Integral(90, -1) / h.Integral(0, -1) << std::endl;
244  std::cout <<"*****************************" << std::endl;
245  totExceptions += h.Integral(90, -1);
246  totPts += h.Integral(0, -1);
247  gPad->Print(TString::Format("plots/chisq_opdet_%d.eps", opdetIdx).Data());
248 
249  delete c1;
250  } // end for opdetIdx
251 
252  std::cout << totExceptions << " exceptions from " << totPts << " points = " << (100.*totExceptions)/totPts << "%" << std::endl;
253 
254  delete fout_full;
255  delete fout_fit;
256  exit(0); // We're done :)
257  }
258 
259  //--------------------------------------------------------------------
261  {
262  }
263 
265  } // namespace
CreateHybridLibrary(const fhicl::ParameterSet &p)
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Number of OpDets in the whole detector.
Encapsulate the construction of a single cyostat.
Float_t vox
Definition: plot.C:40
void analyze(const art::Event &e) override
void GetCenter(double *xyz, double localz=0.0) const
Definition: OpDetGeo.cxx:47
TVector3 GetCenter() const
int GetNVoxels() const
TCanvas * c1
Definition: plotHisto.C:7
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
Float_t d
Definition: plot.C:237
EDAnalyzer(Table< Config > const &config)
Definition: EDAnalyzer.h:100
unsigned int NOpDets() const
Number of OpDets in the whole detector.
Encapsulate the geometry of an optical detector.
const sim::PhotonVoxelDef & GetVoxelDef() const
double DistanceToPoint(geo::Point_t const &point) const
Returns the distance of the specified point from detector center [cm].
Definition: OpDetGeo.cxx:97
General LArSoft Utilities.
CreateHybridLibrary & operator=(const CreateHybridLibrary &)=delete
Float_t e
Definition: plot.C:34
art framework interface to geometry description
float GetVisibility(double const *xyz, unsigned int OpChannel, bool wantReflected=false) const
PhotonVoxel GetPhotonVoxel(int ID) const