Class Reference for E1039 Core & Analysis Software
MakeRTCurve.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <sstream>
3 #include <fstream>
4 #include <TFile.h>
5 #include <TTree.h>
6 #include <TH2D.h>
7 #include <TCanvas.h>
8 #include <TClonesArray.h>
9 #include <TSystem.h>
10 #include <TStyle.h>
11 #include <TGraph.h>
12 #include <TLegend.h>
13 #include <geom_svc/GeomSvc.h>
14 //#include <SRecEvent.h>
15 #include <ktracker/FastTracklet.h>
16 #include "RTCurve.h"
17 #include "CalibParam.h"
18 #include "CalibData.h"
19 #include "FitRTDist.h"
20 #include "MakeRTCurve.h"
21 using namespace std;
22 
23 MakeRTCurve::MakeRTCurve(const int iter)
24  : m_iter(iter)
25  , m_dir_name_out("")
26  , m_n_evt_all(0)
27  , m_n_evt_ana(0)
28  , m_n_trk_all(0)
29  , m_n_trk_ana(0)
30  , m_verb (0)
31 {
32  if (iter <= 0) {
33  cout << "'iter' must be a positive integer. Abort." << endl;
34  exit(1);
35  }
36  //gErrorIgnoreLevel = 1111;
37  //GeomSvc* geom = GeomSvc::instance();
38  //geom->init();
39 }
40 
42 {
43  ;
44 }
45 
47 {
48  ostringstream oss;
49  oss << "calib/" << m_iter;
50  m_dir_name_out = oss.str();
51  gSystem->mkdir(m_dir_name_out.c_str(), kTRUE);
52 
53  cal_par.SetAnaPlanes(false, false, false, false, true); // (d0, d1, d2, d3p, d3m)
54  cal_par.Init(N_RT_PT);
55 
56  oss.str("");
57  oss << "calib/" << (m_iter-1) << "/param.tsv";
58  cal_par.ReadRTParam(oss.str());
59 
60  //const char* env_str = gSystem->Getenv("FIX_TIME_WINDOW");
61  //if (! env_str) {
62  // cout << "!!ERROR!! FIX_TIME_WINDOW is not defined. Abort." << endl;
63  // exit(1);
64  //}
65  //if (strcmp(env_str, "true" ) == 0) cal_par.ReadTimeWindow("calib/time_window.txt");
66 
67  cal_dat.Init(&cal_par);
68 }
69 
70 void MakeRTCurve::AnalyzeListOfFiles(const char* fn_list)
71 {
72  ifstream ifs(fn_list);
73  string fn_file;
74  while (ifs >> fn_file) AnalyzeFile(fn_file.c_str());
75  ifs.close();
76 }
77 
78 void MakeRTCurve::AnalyzeFile(const char* fname)
79 {
80  if (gSystem->AccessPathName(fname)) return;
81  cout << "Input: " << fname << endl;
82  TFile* dataFile = new TFile(fname, "READ");
83  TTree* dataTree = (TTree*)dataFile->Get("eval");
84 
85  // SRawEvent is not available in the evaluation file. Need to use DST.
86  //static SRawEvent* raw = new SRawEvent();
87  //static SRecEvent* rec = 0; // new SRecEvent();
88  static TClonesArray* tracklets = new TClonesArray("Tracklet");
89  tracklets->Clear();
90  //dataTree->SetBranchAddress("rawEvent" , &raw );
91  //dataTree->SetBranchAddress("recEvent" , &rec );
92  dataTree->SetBranchAddress("tracklets", &tracklets);
93 
94  int n_evt_ana = 0;
95  int n_trk_all = 0;
96  int n_trk_ana = 0;
97  int n_evt = dataTree->GetEntries();
98  for (int i_evt = 0; i_evt < n_evt; i_evt++) {
99  dataTree->GetEntry(i_evt);
100  //bool nim1 = raw->isTriggeredBy(SRawEvent::NIM1); // H1234
101  //bool nim2 = raw->isTriggeredBy(SRawEvent::NIM2); // H12
102  //bool nim4 = raw->isTriggeredBy(SRawEvent::NIM4); // H24
103  //if (! nim4) continue;
104  n_evt_ana++;
105 
106  int n_trk = tracklets->GetEntries();
107  n_trk_all += n_trk;
108 
109  map<int, int> list_n_trk; // <station ID, N of tracks>
110  for (int i_trk = 0; i_trk < n_trk; i_trk++) {
111  Tracklet* trk = (Tracklet*)tracklets->At(i_trk);
112  //cal_dat.FillTracklet(trk); // Use this to fill all tracklets
113  list_n_trk[trk->stationID]++;
114  }
115  if (list_n_trk[ST_ID_D3M] > 64) continue; // Exclude noisy events
116 
117  int rec_st = 0; // rec->getRecStatus()
118  cal_dat.FillEventInfo(rec_st, list_n_trk);
119 
120  int i_trk_best = -1;
121  double rchi2_best = 0;
122  for (int i_trk = 0; i_trk < n_trk; i_trk++) {
123  Tracklet* trk = (Tracklet*)tracklets->At(i_trk);
124  int st_id = trk->stationID;
125  if (st_id != ST_ID_D3M) continue;
126 
127  int n_hit = trk->getNHits();
128  int ndf = n_hit - 4; // Correct only when KMag is off
129  double rchi2 = trk->chisq / ndf;
130  if (n_hit != 6) continue;
131  if (rchi2 > 3.0) continue;
132 
133  double x_d2, y_d2;
134  double x_d3, y_d3;
135  EvalD2XY(trk, x_d2, y_d2);
136  EvalD3XY(trk, x_d3, y_d3);
137 
138  if (fabs(trk->tx) > 0.3 || fabs(trk->ty) > 0.4 ||
139  fabs(x_d3) > 131 || fabs(y_d3+85) > 85 ) continue; // Track in the D3m acceptance, page 15 of doc9856.
140 
141  // Cut L1 = st_id & tx & ty & x_d3 & y_d3
142  // Cut L2 = L1 & n_hit & n_trk <= 64
143  // Cut L3 = L2 & rchi2 < 2
144 
145  if (i_trk_best < 0 || rchi2 < rchi2_best) {
146  i_trk_best = i_trk;
147  rchi2_best = rchi2;
148  }
149  }
150 
151  if (i_trk_best >= 0) {
152  Tracklet* trk = (Tracklet*)tracklets->At(i_trk_best);
153  cal_dat.FillTracklet(trk); // Use this to fill only good tracklets
154  n_trk_ana++;
155 
156  for(std::list<SignedHit>::iterator it = trk->hits.begin(); it != trk->hits.end(); ++it) {
157  if(it->hit.index < 0) continue; // Probably not necessary.
158  //int sign = it->hit.sign;
159  int det_id = it->hit.detectorID;
160  int ele_id = it->hit.elementID;
161  double drift_dist = it->hit.driftDistance;
162  double tdc_time = it->hit.tdcTime;
163  double track_pos = trk->getExpPositionW(det_id);
164  double wire_pos = it->hit.pos; // GeomSvc::instance()->getMeasurement(det_id, ele_id);
165  cal_dat.FillHit(det_id, drift_dist, tdc_time, track_pos, wire_pos);
166  //cout << "D " << det_id << " " << it->hit.pos << " " << tdc_time << " " << drift_dist << " " << track_dist << endl;
167  }
168  }
169 
170  tracklets->Clear();
171  }
172  dataFile->Close();
173  cout << " N of all events = " << n_evt << ", analyzed events = " << n_evt_ana << ", tracks = " << n_trk_ana << endl;
174  m_n_evt_all += n_evt;
175  m_n_evt_ana += n_evt_ana;
176  m_n_trk_all += n_trk_all;
177  m_n_trk_ana += n_trk_ana;
178 }
179 
181 {
182  cal_dat.DrawHistEvent(m_dir_name_out);
183  string fname = m_dir_name_out + "/info_event.txt";
184  ofstream ofs(fname.c_str());
185  ofs << m_n_evt_all << "\n"
186  << m_n_evt_ana << "\n"
187  << m_n_trk_all << "\n"
188  << m_n_trk_ana << "\n";
189  ofs.close();
190 }
191 
193 {
194  cal_dat.DrawHistHit(m_dir_name_out);
195 }
196 
198 {
199  cout << "ExtractRT()" << endl;
200  FitRTDist* fit = new FitRTDist();
201  fit->Verbosity(m_verb);
202  for (int ip = 0; ip < cal_par.GetNumPlanes(); ip++) {
203  if (! cal_par.GetAnaPlane(ip)) continue;
204 
205  cout << " Plane " << ip+1 << endl;
206  if (cal_par.TimeWindowIsFixed()) {
207  double t1 = cal_par.GetT1(ip);
208  double t0 = cal_par.GetT0(ip);
209  fit->FixT1T0(t1, t0);
210  }
211  TH2* h2_rt = cal_dat.GetHistRT (ip);
212  double r_max = cal_par.GetRMax (ip);
213  TGraph* gr_t2r = cal_par.GetGraphT2R(ip);
214  fit->DoFit(N_RT_PT, h2_rt, r_max, gr_t2r, cal_par.GetRTCurve(ip));
215  }
216  delete fit;
217 }
218 
220 {
221  cal_par.WriteRTParam(m_dir_name_out, "param.tsv");
222  cal_par.WriteRTGraph(m_dir_name_out, "rt_graph.root");
223 }
224 
226 {
227  GeomSvc* geom = GeomSvc::instance();
228  cout << "DrawCalibResult()" << endl;
229  TCanvas* c1 = new TCanvas("c1", "", 800, 600);
230  c1->SetGrid();
231  gStyle->SetOptStat(0);
232  TLegend* leg = new TLegend(0.7, 0.7, 0.9, 0.9);
233  leg->SetTextFont(22);
234  leg->SetBorderSize(1);
235  leg->SetFillColor(0);
236 
237  ostringstream oss;
238  oss << setfill('0');
239  for(int ip = 0; ip < cal_par.GetNumPlanes(); ip++) {
240  if (! cal_par.GetAnaPlane(ip)) continue;
241 
242  leg->Clear();
243 
244  cal_dat.GetHistRT(ip)->Draw("colz");
245 
246  TGraph* gr_t2r_in = cal_par.GetGraphT2R(ip);
247  TSpline3* spl_t2r = cal_par.GetRTCurve(ip)->GetT2RSpline();
248  TGraph* gr_t2dr = cal_par.GetRTCurve(ip)->GetT2DRGraph();
249 
250  spl_t2r->SetMarkerStyle(10);
251  spl_t2r->SetLineWidth(2);
252  spl_t2r->Draw("sameCP"); // LPE1
253 
254  gr_t2dr->SetMarkerColor(kGreen);
255  gr_t2dr->SetMarkerStyle(10);
256  gr_t2dr->SetLineWidth(2);
257  gr_t2dr->Draw("sameP");
258 
259  gr_t2r_in->SetLineStyle(2);
260  gr_t2r_in->SetLineWidth(2);
261  gr_t2r_in->Draw("sameC"); // LP
262 
263  leg->AddEntry(spl_t2r , "Output", "l");
264  leg->AddEntry(gr_t2r_in, "Input" , "l");
265  leg->Draw();
266 
267  oss.str("");
268  oss << m_dir_name_out << "/rt_result_" << setw(2) << (ip+1) << "_" << geom->getDetectorName(ip+1) << ".png";
269 
270  c1->SaveAs(oss.str().c_str());
271  }
272  delete leg;
273 
274  int n_pl = cal_par.GetNumPlanes();
275  TH1* h1_res = new TH1D("h1_res", ";;Resolution (cm)", n_pl, 0, n_pl);
276  for (int ip = 0; ip < n_pl; ip++) {
277  if (! cal_par.GetAnaPlane(ip)) continue;
278  string det_name = geom->getDetectorName(ip+1);
279  h1_res->SetBinContent(ip+1, cal_par.GetRTCurve(ip)->GetRWidth());
280  h1_res->GetXaxis()->SetBinLabel(ip+1, det_name.c_str());
281  }
282  h1_res->GetYaxis()->SetRangeUser(0, 0.08);
283  h1_res->SetMarkerColor(kRed);
284  h1_res->SetMarkerStyle(10);
285  h1_res->GetXaxis()->LabelsOption("v");
286  h1_res->Draw("P");
287  oss.str("");
288  oss << m_dir_name_out << "/h1_res.png";
289  c1->SaveAs(oss.str().c_str());
290  delete c1;
291 
292  //const char* NAME_GROUP[4] = {"D1", "D2", "D3p", "D3m"};
293  //const double RES_MAX[4] = {0.1, 0.1, 0.1, 0.1};
294  //for (int ig = 0; ig < 4; ig++) {
295  // TH1* fr = c1->DrawFrame(TMIN[ig*6], 0, TMAX[ig*6], RES_MAX[ig], "Resolution vs TDC ;TDC time (ns);Resolution (cm)");
296  // for (int ip = 0; ip < 6; ip++) {
297  // gr_t2dr[ ig*6 +ip ]->SetLineColor (ip+1);
298  // gr_t2dr[ ig*6 +ip ]->SetMarkerColor(ip+1);
299  // gr_t2dr[ ig*6 +ip ]->Draw("sameCP"); // LPE1
300  //
301  // oss.str(""); oss << geom->getDetectorName( ig*6 + ip + 1);
302  // }
303  // c1->SetGrid();
304  // oss.str(""); oss << m_dir_name_out << "/res_" << NAME_GROUP[ig] <<".png";
305  // c1->SaveAs(oss.str().c_str());
306  //}
307 }
308 
309 void MakeRTCurve::EvalD2XY(const Tracklet* trk, double& x, double& y)
310 {
311  x = trk->x0 + 1345 * trk->tx;
312  y = trk->y0 + 1345 * trk->ty;
313 }
314 
315 void MakeRTCurve::EvalD3XY(const Tracklet* trk, double& x, double& y)
316 {
317  x = trk->x0 + 1900 * trk->tx;
318  y = trk->y0 + 1900 * trk->ty;
319 }
320 
void FillTracklet(const Tracklet *trk)
Definition: CalibData.cc:123
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
TH2 * GetHistRT(const int ipl) const
Definition: CalibData.h:49
void Init(CalibParam *ptr)
Definition: CalibData.cc:28
bool GetAnaPlane(const int i_pl)
Definition: CalibParam.h:31
void Init(const int n_rt_pt)
Definition: CalibParam.cc:35
double GetRMax(const int i_pl) const
Definition: CalibParam.h:38
void WriteRTParam(const std::string dir_name, const std::string fname)
Definition: CalibParam.cc:120
double GetT0(const int i_pl) const
Definition: CalibParam.h:39
void WriteRTGraph(const std::string dir_name, const std::string fname)
Definition: CalibParam.cc:172
TGraph * GetGraphT2R(const int i_pl) const
Definition: CalibParam.h:42
double GetT1(const int i_pl) const
Definition: CalibParam.h:40
int GetNumPlanes() const
Definition: CalibParam.h:28
void ReadRTParam(const std::string fname)
Definition: CalibParam.cc:82
RTCurve * GetRTCurve(const int i_pl) const
Definition: CalibParam.h:34
void SetAnaPlanes(const bool d0, const bool d1, const bool d2, const bool d3p, const bool d3m)
Definition: CalibParam.cc:26
bool TimeWindowIsFixed() const
Definition: CalibParam.h:33
Class for fitting R-T histogram.
Definition: FitRTDist.h:9
void Verbosity(const int verb)
Definition: FitRTDist.h:34
int DoFit(const int n_pt, TH2 *h2, double r_max, TGraph *gr_init, RTCurve *rtc)
Definition: FitRTDist.cc:88
void FixT1T0(const double t1, const double t0)
Definition: FitRTDist.cc:81
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:223
void DrawHistHit()
Definition: MakeRTCurve.cc:192
void Init()
Definition: MakeRTCurve.cc:46
MakeRTCurve(const int iter)
Definition: MakeRTCurve.cc:23
void AnalyzeFile(const char *fname)
Definition: MakeRTCurve.cc:78
void ExtractRT()
Definition: MakeRTCurve.cc:197
void DrawHistEvent()
Definition: MakeRTCurve.cc:180
void DrawCalibResult()
Definition: MakeRTCurve.cc:225
void AnalyzeListOfFiles(const char *fn_list)
Definition: MakeRTCurve.cc:70
virtual ~MakeRTCurve()
Definition: MakeRTCurve.cc:41
void WriteRT()
Definition: MakeRTCurve.cc:219
double GetRWidth() const
Definition: RTCurve.cc:48
TSpline3 * GetT2RSpline()
Definition: RTCurve.cc:63
TGraph * GetT2DRGraph()
Definition: RTCurve.cc:73
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