Class Reference for E1039 Core & Analysis Software
CalibData.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <sstream>
3 #include <TH2D.h>
4 #include <TCanvas.h>
5 #include <TSystem.h>
6 #include <TStyle.h>
7 #include <THStack.h>
8 #include <TGraph.h>
9 #include <TLegend.h>
10 #include <geom_svc/GeomSvc.h>
11 #include <UtilAna/UtilHist.h>
12 #include <ktracker/FastTracklet.h>
13 #include "CalibParam.h"
14 #include "CalibData.h"
15 using namespace std;
16 
18  : cal_par(0)
19 {
20  ;
21 }
22 
24 {
25  ;
26 }
27 
29 {
30  cal_par = ptr;
31 
32  h1_stat = new TH1D("h1_stat" , ";Reco. status;N of tracks", 25, -24.5, 0.5);
33 
34  int id_n = ST_ID_MAX - ST_ID_MIN + 1;
35  double id_l = ST_ID_MIN - 0.5;
36  double id_h = ST_ID_MAX + 0.5;
37  h1_st_id = new TH1D("h1_st_id", ";Station ID;N of tracks", id_n, id_l, id_h);
38  h2_ntrk = new TH2D("h2_ntrk" , ";N of tracks/event;Station ID;N of tracks", 200, -0.5, 199.5, id_n, id_l, id_h);
39  h2_nhit = new TH2D("h2_nhit" , ";N of hits/track;Station ID;N of tracks", 18, 0.5, 18.5, id_n, id_l, id_h);
40  h2_mom = new TH2D("h2_mom" , ";Momentum (GeV);Station ID;N of tracks", 22, 0.0, 120.0, id_n, id_l, id_h);
41  h2_rchi2 = new TH2D("h2_rchi2", ";#chi^{2}/NDF;Station ID;N of tracks", 50, 0.0, 10.0, id_n, id_l, id_h);
42  h2_x0 = new TH2D("h2_x0" , ";x_{0} (cm);Station ID;N of tracks", 100, -500, 500, id_n, id_l, id_h);
43  h2_y0 = new TH2D("h2_y0" , ";y_{0} (cm);Station ID;N of tracks", 100, -500, 500, id_n, id_l, id_h);
44  h2_tx = new TH2D("h2_tx" , ";t_{x};Station ID;N of tracks", 100, -1.0, 1.0, id_n, id_l, id_h);
45  h2_ty = new TH2D("h2_ty" , ";t_{y};Station ID;N of tracks", 100, -1.0, 1.0, id_n, id_l, id_h);
46  h2_x_st = new TH2D("h2_x_st" , ";x_{St} (cm);Station ID;N of tracks", 100, -250, 250, id_n, id_l, id_h);
47  h2_y_st = new TH2D("h2_y_st" , ";y_{St} (cm);Station ID;N of tracks", 100, -250, 250, id_n, id_l, id_h);
48 
49  GeomSvc* geom = GeomSvc::instance();
50  ostringstream oss;
51  ostringstream oss2;
52  oss2 << setfill('0');
53  for (int ip = 0; ip < cal_par->GetNumPlanes(); ip++) {
54  string det_name = geom->getDetectorName(ip+1);
55  double det_dy = geom->getPlaneScaleX (ip+1);
56  double t_min = cal_par->GetTMin(ip);
57  double t_max = cal_par->GetTMax(ip);
58  double r_max = cal_par->GetRMax(ip);
59 
60  const double t_width = 4.0/9.0 * 3; // 4/9 ns * integer
61  oss.str("");
62  oss << "h1_time_wide_" << det_name;
63  oss2.str("");
64  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";TDC time (ns);Yields";
65  h1_time_wide[ip] = new TH1D(oss.str().c_str(), oss2.str().c_str(), 2000, 0, 2000*t_width);
66 
67  int bin_t_lo = (int)(t_min / t_width);
68  int bin_t_hi = (int)(t_max / t_width) + 1;
69  oss.str("");
70  oss << "h1_time_" << det_name;
71  oss2.str("");
72  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";TDC time (ns);Yields";
73  h1_time[ip] = new TH1D(oss.str().c_str(), oss2.str().c_str(), bin_t_hi-bin_t_lo, t_width*(bin_t_lo)-0.01, t_width*(bin_t_hi)-0.01);
74 
75  oss.str("");
76  oss << "h2_tp_wp_" << det_name;
77  oss2.str("");
78  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";Track position (cm);Wire position (cm);";
79  h2_tp_wp[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), 120, -0.6*det_dy, 0.6*det_dy, 120, -0.6*det_dy, 0.6*det_dy);
80 
81  oss.str("");
82  oss << "h2_rt_" << det_name;
83  oss2.str("");
84  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";TDC time (ns);R (cm);";
85  h2_RT[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), T_BIN, t_min, t_max, R_BIN, 0, 1.2*r_max );
86 
87  oss.str("");
88  oss << "h2_resi_" << det_name;
89  oss2.str("");
90  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";TDC time (ns);Residual (cm);";
91  h2_Resi[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), T_BIN, t_min, t_max, R_BIN, -0.5, 0.5);
92 
93  oss.str("");
94  oss << "h2_tr_" << det_name;
95  oss2.str("");
96  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";TDC time (ns);R (cm);";
97  h2_TR[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), R_BIN, 0, 1.2*r_max, T_BIN, t_min, t_max );
98 
99  oss.str("");
100  oss << "h2_tr_ex_" << det_name;
101  oss2.str("");
102  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";Track distance (cm);TDC time (ns);";
103  h2_TR_ex[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), R_BIN*2, -1.2*r_max, 1.2*r_max, T_BIN, t_min, t_max );
104 
105  oss.str("");
106  oss << "h2_td_dd_" << det_name;
107  oss2.str("");
108  oss2 << "#" << setw(2) << (ip+1) << ": " << det_name << ";Track distance (cm);Drift distance (cm);Yields";
109  h2_td_dd[ip] = new TH2D(oss.str().c_str(), oss2.str().c_str(), R_BIN*2, -1.2*r_max, 1.2*r_max, R_BIN, 0, 1.2*r_max);
110  }
111 }
112 
113 void CalibData::FillEventInfo(const int rec_stat, const std::map<int, int> list_n_trk)
114 {
115  h1_stat->Fill(rec_stat);
116 
117  for (int st_id = ST_ID_MIN; st_id <= ST_ID_MAX; st_id++) {
118  int n_trk = list_n_trk.find(st_id) != list_n_trk.end() ? list_n_trk.at(st_id) : 0;
119  h2_ntrk->Fill(n_trk, st_id);
120  }
121 }
122 
124 {
125  int st_id = trk->stationID;
126  int n_hit = trk->getNHits();
127  int ndf = n_hit - (st_id == 7 ? 5 : 4); // Incorrect, when KMag is off
128  double rchi2 = ndf > 0 ? trk->chisq / ndf : 0;
129 
130  double z_st = CalibParam::ZOfStationID(st_id);
131  double x_st = trk->x0 + z_st * trk->tx;
132  double y_st = trk->y0 + z_st * trk->ty;
133 
134  h1_st_id->Fill(st_id);
135  h2_nhit ->Fill(n_hit , st_id);
136  h2_mom ->Fill(1/trk->invP, st_id);
137  h2_rchi2->Fill(rchi2 , st_id);
138  h2_x0 ->Fill(trk->x0 , st_id);
139  h2_y0 ->Fill(trk->y0 , st_id);
140  h2_tx ->Fill(trk->tx , st_id);
141  h2_ty ->Fill(trk->ty , st_id);
142  h2_x_st ->Fill(x_st , st_id);
143  h2_y_st ->Fill(y_st , st_id);
144 }
145 
147 {
148  for(auto it = trk->hits.begin(); it != trk->hits.end(); ++it) {
149  if(it->hit.index < 0) continue; // Probably not necessary.
150  //int sign = it->hit.sign;
151  int det_id = it->hit.detectorID;
152  int ele_id = it->hit.elementID;
153  double drift_dist = it->hit.driftDistance;
154  double tdc_time = it->hit.tdcTime;
155  double track_pos = trk->getExpPositionW(det_id);
156  double wire_pos = it->hit.pos; // GeomSvc::instance()->getMeasurement(det_id, ele_id);
157  FillHit(det_id, drift_dist, tdc_time, track_pos, wire_pos);
158  //cout << "D " << det_id << " " << it->hit.pos << " " << tdc_time << " " << drift_dist << " " << track_dist << endl;
159  }
160 }
161 
163 
166 void CalibData::FillHit(const int det_id, const double drift_dist, const double tdc_time, const double track_pos, const double wire_pos)
167 {
168  double track_dist = track_pos - wire_pos;
169  h1_time_wide[det_id-1]->Fill(tdc_time);
170  h1_time [det_id-1]->Fill(tdc_time);
171  h2_tp_wp [det_id-1]->Fill(track_pos, wire_pos);
172  h2_RT [det_id-1]->Fill(tdc_time, fabs(track_dist));
173  h2_Resi [det_id-1]->Fill(tdc_time, fabs(drift_dist) - fabs(track_dist));
174  h2_TR [det_id-1]->Fill(fabs(track_dist), tdc_time);
175  h2_TR_ex [det_id-1]->Fill(track_dist, tdc_time);
176  h2_td_dd [det_id-1]->Fill(track_dist, drift_dist);
177 }
178 
179 void CalibData::DrawHistEvent(const string dir_out)
180 {
181  cout << "DrawHistEvent()" << endl;
182  ostringstream oss;
183  TCanvas* c1 = new TCanvas("c1", "", 800, 600);
184  c1->SetGrid();
185  gStyle->SetOptStat(0);
186 
187  //h1_stat->Draw();
188  //oss.str("");
189  //oss << dir_out << "/h1_rec_stat.png";
190  //c1->SaveAs(oss.str().c_str());
191 
192  h1_st_id->Draw();
193  oss.str("");
194  oss << dir_out << "/h1_st_id.png";
195  c1->SaveAs(oss.str().c_str());
196 
197  const int iy_lo = 1; // 4 = D3p, 5 = D3m, 6 = D23
198  const int iy_hi = 5;
199  DrawIn1D(h2_ntrk , "ntrk" , dir_out, iy_lo, iy_hi);
200  DrawIn1D(h2_nhit , "nhit" , dir_out, iy_lo, iy_hi);
201  //DrawIn1D(h2_mom , "mom" , dir_out, iy_lo, iy_hi);
202  DrawIn1D(h2_rchi2, "rchi2", dir_out, iy_lo, iy_hi);
203  //DrawIn1D(h2_x0 , "x0" , dir_out, iy_lo, iy_hi);
204  //DrawIn1D(h2_y0 , "y0" , dir_out, iy_lo, iy_hi);
205  DrawIn1D(h2_tx , "tx" , dir_out, iy_lo, iy_hi);
206  DrawIn1D(h2_ty , "ty" , dir_out, iy_lo, iy_hi);
207  DrawIn1D(h2_x_st , "x_st" , dir_out, iy_lo, iy_hi);
208  DrawIn1D(h2_y_st , "y_st" , dir_out, iy_lo, iy_hi);
209 
210  delete c1;
211 }
212 
213 void CalibData::DrawHistHit(const string dir_out)
214 {
215  cout << "DrawHistHit()" << endl;
216  ostringstream oss;
217  TCanvas* c1 = new TCanvas("c1", "", 800, 600);
218  c1->SetGrid();
219  gStyle->SetOptStat(0);
220 
221  GeomSvc* geom = GeomSvc::instance();
222 
223  oss << setfill('0');
224  for (int ip = 0; ip < cal_par->GetNumPlanes(); ip++) {
225  if (! cal_par->GetAnaPlane(ip)) continue;
226  string det_name = geom->getDetectorName(ip+1);
227  //cout << " Plane " << ip+1 << endl;
228 
229  oss.str("");
230  oss << setw(2) << (ip+1) << "_" << det_name << ".png";
231  string suffix = oss.str();
232 
233  UtilHist::AutoSetRange(h1_time_wide[ip]);
234  h1_time_wide[ip]->Draw();
235  c1->SaveAs( (dir_out + "/h1_time_wide_" + suffix).c_str() );
236 
237  h1_time[ip]->Draw();
238  c1->SaveAs( (dir_out + "/h1_time_" + suffix).c_str() );
239 
240  //h2_tp_wp[ip]->Draw("colz"); // To check if any outlier exists.
241  //c1->SaveAs( (dir_out + "/h2_tp_wp_" + suffix).c_str() );
242 
243  TH1* h1_tp = h2_tp_wp[ip]->ProjectionX("h1_tp");
244  h1_tp->Draw();
245  c1->SaveAs( (dir_out + "/h1_tp_" + suffix).c_str() );
246  delete h1_tp;
247 
248  //TH1* h1_wp = h2_tp_wp[ip]->ProjectionY("h1_wp");
249  //h1_wp->Draw();
250  //c1->SaveAs( (dir_out + "/h1_wp_" + suffix).c_str() );
251  //delete h1_wp;
252 
253  h2_TR_ex[ip]->Draw("colz");
254  c1->SaveAs( (dir_out + "/h2_tr_ex_" + suffix).c_str() );
255 
256  h2_Resi[ip]->Draw("colz");
257  c1->SaveAs( (dir_out + "/h2_resi_" + suffix).c_str() );
258 
259  h2_td_dd[ip]->Draw("colz");
260  c1->SaveAs( (dir_out + "/h2_td_dd_" + suffix).c_str() );
261 
262  TH1* h1_td = h2_td_dd[ip]->ProjectionX("h1_td");
263  h1_td->Draw();
264  c1->SaveAs( (dir_out + "/h1_td_" + suffix).c_str() );
265  delete h1_td;
266 
267  TH1* h1_dd = h2_td_dd[ip]->ProjectionY("h1_dd");
268  h1_dd->Draw();
269  c1->SaveAs( (dir_out + "/h1_dd_" + suffix).c_str() );
270  delete h1_dd;
271  }
272  delete c1;
273 }
274 
275 void CalibData::DrawIn1D(TH2* h2, const string label, const string dir_out, int iy_lo, int iy_hi)
276 {
277  if (iy_lo <= 0) iy_lo = 1;
278  if (iy_hi <= 0) iy_hi = h2->GetNbinsY();
279 
280  ostringstream oss;
281  oss << h2->GetTitle() << ";" << h2->GetXaxis()->GetTitle() << ";" << h2->GetZaxis()->GetTitle();
282  THStack hs("hs", oss.str().c_str());
283  TLegend leg(0.90, 0.90-0.03*(iy_hi-iy_lo+2), 0.99, 0.99);
284  leg.SetHeader("St. ID");
285  for (int iy = iy_lo; iy <= iy_hi; iy++) {
286  oss.str("");
287  oss << h2->GetName() << "_" << iy;
288  TH1* h1 = h2->ProjectionX(oss.str().c_str(), iy, iy);
289  h1->SetLineColor(iy - iy_lo + 1);
290  hs.Add(h1);
291  oss.str("");
292  //oss << h2->GetYaxis()->GetTitle() << h2->GetYaxis()->GetBinCenter(iy);
293  oss << h2->GetYaxis()->GetBinCenter(iy);
294  leg.AddEntry(h1, oss.str().c_str(), "l");
295  }
296 
297  TCanvas* c1 = new TCanvas("c1_din1d", "", 800, 600);
298  c1->SetGrid();
299  c1->SetLogy();
300  gStyle->SetOptStat(0);
301  hs.Draw("nostack");
302  leg.Draw();
303 
304  oss.str("");
305  oss << dir_out << "/h1_" << label << ".png";
306  c1->SaveAs(oss.str().c_str());
307  delete c1;
308 }
void FillTracklet(const Tracklet *trk)
Definition: CalibData.cc:123
virtual ~CalibData()
Definition: CalibData.cc:23
void FillTrackletHits(const Tracklet *trk)
Definition: CalibData.cc:146
void DrawHistHit(const std::string dir_out)
Definition: CalibData.cc:213
void DrawHistEvent(const std::string dir_out)
Definition: CalibData.cc:179
void FillEventInfo(const int rec_stat, const std::map< int, int > list_n_trk)
Definition: CalibData.cc:113
void FillHit(const int det_id, const double drift_dist, const double tdc_time, const double track_pos, const double wire_pos)
Fill the hit information.
Definition: CalibData.cc:166
void Init(CalibParam *ptr)
Definition: CalibData.cc:28
Class to hold the calibration parameters.
Definition: CalibParam.h:8
bool GetAnaPlane(const int i_pl)
Definition: CalibParam.h:31
double GetTMax(const int i_pl) const
Definition: CalibParam.h:37
double GetRMax(const int i_pl) const
Definition: CalibParam.h:38
int GetNumPlanes() const
Definition: CalibParam.h:28
static double ZOfStationID(const int st_id)
Definition: CalibParam.cc:224
double GetTMin(const int i_pl) const
Definition: CalibParam.h:36
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
double getPlaneScaleX(int detectorID) const
Definition: GeomSvc.h:242
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:223
double invP
Definition: FastTracklet.h:239
std::list< SignedHit > hits
Definition: FastTracklet.h:228
double y0
Definition: FastTracklet.h:238
int stationID
Definition: FastTracklet.h:214
double ty
Definition: FastTracklet.h:236
double tx
Definition: FastTracklet.h:235
double x0
Definition: FastTracklet.h:237
int getNHits() const
Definition: FastTracklet.h:145
double chisq
Definition: FastTracklet.h:222
double getExpPositionW(int detectorID) const
void AutoSetRange(TH1 *h1, const int margin_lo=5, const int margin_hi=5)
Adjust the axis range via "h1->GetXaxis()->SetRange(bin_lo, bin_hi)" to zoom up non-empty bins.
Definition: UtilHist.cc:17