Class Reference for E1039 Core & Analysis Software
CalibParam.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <sstream>
3 #include <fstream>
4 #include <TFile.h>
5 #include <TGraph.h>
6 #include <TSpline.h>
7 #include <geom_svc/GeomSvc.h>
8 #include "RTCurve.h"
9 #include "CalibParam.h"
10 using namespace std;
11 
12 const double CalibParam::DT_RT = 2.5;
13 
15  : fix_time_window(false)
16 {
17  memset(m_ana_pl, 0, sizeof(m_ana_pl));
18  memset(m_rtc, 0, sizeof(m_rtc));
19 }
20 
22 {
23  ;
24 }
25 
26 void CalibParam::SetAnaPlanes(const bool d0, const bool d1, const bool d2, const bool d3p, const bool d3m)
27 {
28  for (int ip = 0; ip < 6; ip++) m_ana_pl[ip] = d0;
29  for (int ip = 6; ip < 12; ip++) m_ana_pl[ip] = d1;
30  for (int ip = 12; ip < 18; ip++) m_ana_pl[ip] = d2;
31  for (int ip = 18; ip < 24; ip++) m_ana_pl[ip] = d3p;
32  for (int ip = 24; ip < 30; ip++) m_ana_pl[ip] = d3m;
33 }
34 
35 void CalibParam::Init(const int n_rt_pt)
36 {
37  T_MIN[ 0] = 1200; T_MAX[ 0] = 1330;
38  T_MIN[ 1] = 1200; T_MAX[ 1] = 1330;
39  T_MIN[ 2] = 1200; T_MAX[ 2] = 1330;
40  T_MIN[ 3] = 1200; T_MAX[ 3] = 1330;
41  T_MIN[ 4] = 1200; T_MAX[ 4] = 1330;
42  T_MIN[ 5] = 1200; T_MAX[ 5] = 1330;
43 
44  T_MIN[ 6] = 1700; T_MAX[ 6] = 1800;
45  T_MIN[ 7] = 1700; T_MAX[ 7] = 1800;
46  T_MIN[ 8] = 1700; T_MAX[ 8] = 1800;
47  T_MIN[ 9] = 1700; T_MAX[ 9] = 1800;
48  T_MIN[10] = 1700; T_MAX[10] = 1800;
49  T_MIN[11] = 1700; T_MAX[11] = 1800;
50 
51  T_MIN[12] = 800; T_MAX[12] = 1300;
52  T_MIN[13] = 800; T_MAX[13] = 1300;
53  T_MIN[14] = 800; T_MAX[14] = 1300;
54  T_MIN[15] = 800; T_MAX[15] = 1300;
55  T_MIN[16] = 800; T_MAX[16] = 1300;
56  T_MIN[17] = 800; T_MAX[17] = 1300;
57 
58  T_MIN[18] = 700; T_MAX[18] = 1200;
59  T_MIN[19] = 700; T_MAX[19] = 1200;
60  T_MIN[20] = 700; T_MAX[20] = 1200;
61  T_MIN[21] = 700; T_MAX[21] = 1200;
62  T_MIN[22] = 700; T_MAX[22] = 1200;
63  T_MIN[23] = 700; T_MAX[23] = 1200;
64 
65  T_MIN[24] = 800; T_MAX[24] = 1300;
66  T_MIN[25] = 800; T_MAX[25] = 1300;
67  T_MIN[26] = 800; T_MAX[26] = 1300;
68  T_MIN[27] = 800; T_MAX[27] = 1300;
69  T_MIN[28] = 800; T_MAX[28] = 1300;
70  T_MIN[29] = 800; T_MAX[29] = 1300;
71 
72  GeomSvc* geom = GeomSvc::instance();
73  for (int ip = 0; ip < N_PL; ip++) {
74  R_MAX[ip] = geom->getCellWidth(ip+1) / 2;
75  m_rtc[ip] = new RTCurve(n_rt_pt);
76  m_gr_t2r_in[ip] = new TGraph();
77  }
78 }
79 
82 void CalibParam::ReadRTParam(const string fname)
83 {
84  cout << "ReadRTParam(): input = " << fname << "." << endl;
85  ifstream ifs(fname.c_str());
86  if (! ifs) {
87  cerr << "ERROR: Cannot open the input file. Abort." << endl;
88  exit(1);
89  }
90 
91  char buf[300];
92  while (ifs.getline(buf, 300)) {
93  if (buf[0] != 'D') continue; // comment or PT lines
94  istringstream iss(buf);
95  string det;
96  double t, x, dt, dx;
97  iss >> det >> t >> x >> dt >> dx;
98  int det_id = GeomSvc::instance()->getDetectorID(det);
99  TGraph* gr = m_gr_t2r_in[det_id - 1];
100  gr->SetPoint(gr->GetN(), t, x);
101  }
102  ifs.close();
103 
104  // Invert the points if T is in the descending order.
105  for (int ip = 0; ip < N_PL; ip++) {
106  TGraph* gr = m_gr_t2r_in[ip];
107  int n_pt = gr->GetN();
108  if (n_pt > 0 && gr->GetX()[0] > gr->GetX()[n_pt-1]) {
109  for (int i_pt = 0; i_pt < n_pt / 2; i_pt++) {
110  double t1, x1, t2, x2;
111  gr->GetPoint( i_pt , t1, x1);
112  gr->GetPoint(n_pt-i_pt-1, t2, x2);
113  gr->SetPoint( i_pt , t2, x2);
114  gr->SetPoint(n_pt-i_pt-1, t1, x1);
115  }
116  }
117  }
118 }
119 
120 void CalibParam::WriteRTParam(const string dir_name, const string fname)
121 {
122  cout << "WriteRTParam(): output = " << dir_name << "." << endl;
123  ostringstream oss;
124  oss << dir_name << "/" << fname;
125  ofstream ofs1(oss.str().c_str());
126  ofs1.setf(ios_base::fixed, ios_base::floatfield);
127 
128  ofs1 << "#det\tt\tx\tdt\tdx\n";
129 
130  GeomSvc* geom = GeomSvc::instance();
131  for (int ip = 0; ip < N_PL; ip++) {
132  if (! m_ana_pl[ip]) continue;
133 
134  int n_pt;
135  double t_min, t_max;
136  if (fix_time_window) {
137  t_min = T1[ip];
138  t_max = T0[ip];
139  n_pt = (int)((t_max - t_min) / DT_RT) + 1;
140  } else {
141  n_pt = m_gr_t2r_in[ip]->GetN();
142  t_min = m_gr_t2r_in[ip]->GetX()[0];
143  t_max = m_gr_t2r_in[ip]->GetX()[n_pt-1];
144  }
145 
146  string det = geom->getDetectorName(ip+1);
147  double dx = m_rtc[ip]->GetRWidth();
148  double dt = dx * (t_max - t_min) / R_MAX[ip];
149 
150  for (int i_pt = 0; i_pt < n_pt; i_pt++) {
151  double t = t_min + DT_RT * i_pt;
152  double x = m_rtc[ip]->EvalR(t);
153  ofs1 << det << "\t"
154  << setprecision(1) << t << "\t"
155  << setprecision(4) << x << "\t"
156  << setprecision(3) << dt << "\t"
157  << setprecision(3) << dx << "\n";
158  }
159  }
160 
161  //ifstream ifs( "calibration_prop.txt" );
162  //if(!ifs){
163  // cerr <<"!!ERROR!! Cannot open the input file 'calibration_prop.txt'. Abort." << endl;
164  // exit(1);
165  //}
166  //string buffer;
167  //while (getline(ifs, buffer)) ofs1 << buffer << endl;
168  //ifs .close();
169  ofs1.close();
170 }
171 
172 void CalibParam::WriteRTGraph(const string dir_name, const string fname)
173 {
174  ostringstream oss;
175  oss << dir_name << "/" << fname;
176  TFile* file = new TFile(oss.str().c_str(), "RECREATE");
177 
178  GeomSvc* geom = GeomSvc::instance();
179  for(int ip = 0; ip < N_PL; ip++){
180  if (! m_ana_pl[ip]) continue;
181 
182  oss.str("");
183  oss << "gr_t2r_in_" << geom->getDetectorName(ip+1);
184  m_gr_t2r_in[ip]->SetName(oss.str().c_str());
185  m_gr_t2r_in[ip]->Write();
186 
187  TGraph* gr_t2dr = m_rtc[ip]->GetT2DRGraph();
188  oss.str("");
189  oss << "gr_t2dr_" << geom->getDetectorName(ip+1);
190  gr_t2dr->SetName(oss.str().c_str());
191  gr_t2dr->Write();
192 
193  TSpline3* spl_t2r = m_rtc[ip]->GetT2RSpline();
194  oss.str("");
195  oss << "spl_t2r_" << geom->getDetectorName(ip+1);
196  spl_t2r->SetName(oss.str().c_str());
197  spl_t2r->Write();
198  }
199 
200  file->Write();
201  file->Close();
202 }
203 
204 void CalibParam::ReadTimeWindow(const std::string fname)
205 {
206  cout << "ReadTimeWindow(): input = " << fname << "." << endl;
207  ifstream ifs(fname.c_str());
208  if (! ifs) {
209  cerr << "ERROR: Cannot open the input file. Abort." << endl;
210  exit(1);
211  }
212  int det_id;
213  double t1, t0;
214  string det_name;
215  while (ifs >> det_id >> t1 >> t0 >> det_name) {
216  cout << " " << det_id << "\t" << t1 << "\t" << t0 << "\t" << det_name << "\n";
217  T1[det_id-1] = t1;
218  T0[det_id-1] = t0;
219  }
220  ifs.close();
221  fix_time_window = true;
222 }
223 
224 double CalibParam::ZOfStationID(const int st_id)
225 {
226  if (st_id == 1) return 620; // D0
227  else if (st_id == 3) return 1345; // D2
228  else if (st_id == 4) return 1900; // D3p
229  else if (st_id == 5) return 1900; // D3m
230 
231  cout << "CalibParam::ZOfStationID(): Unsupported station ID (" << st_id << "). Abort." << endl;
232  exit(1);
233 }
void ReadTimeWindow(const std::string fname)
Definition: CalibParam.cc:204
void Init(const int n_rt_pt)
Definition: CalibParam.cc:35
void WriteRTParam(const std::string dir_name, const std::string fname)
Definition: CalibParam.cc:120
void WriteRTGraph(const std::string dir_name, const std::string fname)
Definition: CalibParam.cc:172
void ReadRTParam(const std::string fname)
Definition: CalibParam.cc:82
static double ZOfStationID(const int st_id)
Definition: CalibParam.cc:224
void SetAnaPlanes(const bool d0, const bool d1, const bool d2, const bool d3p, const bool d3m)
Definition: CalibParam.cc:26
virtual ~CalibParam()
Definition: CalibParam.cc:21
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
int getDetectorID(const std::string &detectorName) const
Get the plane position.
Definition: GeomSvc.h:219
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:223
double getCellWidth(int detectorID) const
Definition: GeomSvc.h:238
Class to represent R-T curve.
Definition: RTCurve.h:7
double GetRWidth() const
Definition: RTCurve.cc:48
double EvalR(const double t)
Definition: RTCurve.cc:53
TSpline3 * GetT2RSpline()
Definition: RTCurve.cc:63
TGraph * GetT2DRGraph()
Definition: RTCurve.cc:73