Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFMacro_FineTimeOffsets.C
Go to the documentation of this file.
1 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_FDC_TOF
2 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_FDC_TAGH
3 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_FDC_PSC
4 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_PSC_TAGH
5 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_PSC_TOF
6 // hnamepath: /RF/AverageDeltaT_RF_OtherRFs/RFDeltaT_TAGH_TOF
7 
8 double gBeamSignalPeriod = 2000.0/499.0;
9 int gRebinF1sAmount = 25;
11 
12 Double_t Periodic_Gaussian_Func(Double_t* locXArray, Double_t* locParamArray)
13 {
14  Double_t locValue = locParamArray[0]*TMath::Gaus(locXArray[0], locParamArray[1], locParamArray[2]);
15  if(locParamArray[1] < 0.0)
16  locValue += locParamArray[0]*TMath::Gaus(locXArray[0], locParamArray[1] + gBeamSignalPeriod, locParamArray[2]);
17  else
18  locValue += locParamArray[0]*TMath::Gaus(locXArray[0], locParamArray[1] - gBeamSignalPeriod, locParamArray[2]);
19 
20  return locValue;
21 }
22 
24 {
25  Int_t locMaxBin = locHist->GetMaximumBin();
26  double locMean = locHist->GetBinCenter(locMaxBin);
27 
28  string locFuncName = string(locHist->GetName()) + string("_Func");
29  TF1 *locFunc = new TF1(locFuncName.c_str(), Periodic_Gaussian_Func, -0.5*gBeamSignalPeriod, 0.5*gBeamSignalPeriod, 3);
30  locFunc->SetParameters(locHist->GetBinContent(locMaxBin), locMean, 0.1);
31  locFunc->SetParNames("Height", "Fit #mu", "Fit #sigma");
32 
33  return locFunc;
34 }
35 
36 int RFMacro_FineTimeOffsets(int locRunNumber, string locVariation = "default")
37 {
38  //INPUT RUN NUMBER MUST BE A RUN NUMBER IN THE RUN RANGE YOU ARE TRYING TO COMMIT CONSTANTS TO
39  gDirectory->cd("/"); //return to file base directory
40 
41  //Goto Beam Path
42  TDirectory *locDirectory = (TDirectory*)gDirectory->FindObjectAny("RF");
43  if(!locDirectory)
44  return 0;
45  locDirectory->cd();
46 
47  //Get RF DeltaT Histograms
48  gDirectory->cd("AverageDeltaT_RF_OtherRFs");
49  TH1I* locHist_RFDeltaT_FDC_TOF = (TH1I*)gDirectory->Get("RFDeltaT_FDC_TOF");
50  TH1I* locHist_RFDeltaT_FDC_TAGH = (TH1I*)gDirectory->Get("RFDeltaT_FDC_TAGH");
51  TH1I* locHist_RFDeltaT_FDC_PSC = (TH1I*)gDirectory->Get("RFDeltaT_FDC_PSC");
52  TH1I* locHist_RFDeltaT_PSC_TAGH = (TH1I*)gDirectory->Get("RFDeltaT_PSC_TAGH");
53  TH1I* locHist_RFDeltaT_PSC_TOF = (TH1I*)gDirectory->Get("RFDeltaT_PSC_TOF");
54  TH1I* locHist_RFDeltaT_TAGH_TOF = (TH1I*)gDirectory->Get("RFDeltaT_TAGH_TOF");
55 
56  //Get/Make Canvas
57  TCanvas *locCanvas = NULL;
58  if(TVirtualPad::Pad() == NULL)
59  locCanvas = new TCanvas("RF_FineTimeOffsets", "RF_FineTimeOffsets", 1200, 800); //for testing
60  else
61  locCanvas = gPad->GetCanvas();
62  locCanvas->Divide(3, 2);
63 
64  //Draw
65  locCanvas->cd(1);
66  gPad->SetTicks();
67  gPad->SetGrid();
68  TF1* locTOFFitFunc_TAGH = NULL;
69  if(locHist_RFDeltaT_TAGH_TOF != NULL)
70  {
72  locHist->Rebin(gRebinTOFAmount);
73  locTOFFitFunc_TAGH = Create_FitFunc(locHist);
74  locHist->Fit(locTOFFitFunc_TAGH, "QR");
75 
76  locHist->GetXaxis()->SetTitleSize(0.05);
77  locHist->GetYaxis()->SetTitleSize(0.05);
78  locHist->GetXaxis()->SetLabelSize(0.05);
79  locHist->GetYaxis()->SetLabelSize(0.05);
80  locHist->Draw();
81  }
82 
83  locCanvas->cd(2);
84  gPad->SetTicks();
85  gPad->SetGrid();
86  TF1* locTOFFitFunc_PSC = NULL;
87  if(locHist_RFDeltaT_PSC_TOF != NULL)
88  {
90  locHist->Rebin(gRebinTOFAmount);
91  locTOFFitFunc_PSC = Create_FitFunc(locHist);
92  locHist->Fit(locTOFFitFunc_PSC, "QR");
93 
94  locHist->GetXaxis()->SetTitleSize(0.05);
95  locHist->GetYaxis()->SetTitleSize(0.05);
96  locHist->GetXaxis()->SetLabelSize(0.05);
97  locHist->GetYaxis()->SetLabelSize(0.05);
98  locHist->Draw();
99  }
100 
101  locCanvas->cd(3);
102  gPad->SetTicks();
103  gPad->SetGrid();
104  TF1* locTOFFitFunc_FDC = NULL;
105  if(locHist_RFDeltaT_FDC_TOF != NULL)
106  {
108  locTOFFitFunc_FDC = Create_FitFunc(locHist);
109  locHist->Rebin(gRebinTOFAmount);
110  locHist->Fit(locTOFFitFunc_FDC, "QR");
111 
112  locHist->GetXaxis()->SetTitleSize(0.05);
113  locHist->GetYaxis()->SetTitleSize(0.05);
114  locHist->GetXaxis()->SetLabelSize(0.05);
115  locHist->GetYaxis()->SetLabelSize(0.05);
116  locHist->Draw();
117  }
118 
119  locCanvas->cd(4);
120  gPad->SetTicks();
121  gPad->SetGrid();
122  if(locHist_RFDeltaT_FDC_TAGH != NULL)
123  {
125  locHist->Rebin(gRebinF1sAmount);
126  TF1* locFunc = Create_FitFunc(locHist);
127  locHist->Fit(locFunc, "QR");
128 
129  locHist->GetXaxis()->SetTitleSize(0.05);
130  locHist->GetYaxis()->SetTitleSize(0.05);
131  locHist->GetXaxis()->SetLabelSize(0.05);
132  locHist->GetYaxis()->SetLabelSize(0.05);
133  locHist->Draw();
134  }
135 
136  locCanvas->cd(5);
137  gPad->SetTicks();
138  gPad->SetGrid();
139  if(locHist_RFDeltaT_FDC_PSC != NULL)
140  {
142  locHist->Rebin(gRebinF1sAmount);
143  TF1* locFunc = Create_FitFunc(locHist);
144  locHist->Fit(locFunc, "QR");
145 
146  locHist->GetXaxis()->SetTitleSize(0.05);
147  locHist->GetYaxis()->SetTitleSize(0.05);
148  locHist->GetXaxis()->SetLabelSize(0.05);
149  locHist->GetYaxis()->SetLabelSize(0.05);
150  locHist->Draw();
151  }
152 
153  locCanvas->cd(6);
154  gPad->SetTicks();
155  gPad->SetGrid();
156  if(locHist_RFDeltaT_PSC_TAGH != NULL)
157  {
159  locHist->Rebin(gRebinF1sAmount);
160  TF1* locFunc = Create_FitFunc(locHist);
161  locHist->Fit(locFunc, "QR");
162 
163  locHist->GetXaxis()->SetTitleSize(0.05);
164  locHist->GetYaxis()->SetTitleSize(0.05);
165  locHist->GetXaxis()->SetLabelSize(0.05);
166  locHist->GetYaxis()->SetLabelSize(0.05);
167  locHist->Draw();
168  }
169 
170  //Fine time offsets (with respect to TOF)
171  double locTimeOffset_TOF = 0.0;
172 
173  double locTimeOffset_FDC = locTOFFitFunc_FDC->GetParameter(1);
174  if(!(locHist_RFDeltaT_FDC_TOF->GetEntries() > 0.0))
175  locTimeOffset_FDC = 0.0;
176 
177  double locTimeOffset_PSC = locTOFFitFunc_PSC->GetParameter(1);
178  if(!(locHist_RFDeltaT_PSC_TOF->GetEntries() > 0.0))
179  locTimeOffset_PSC = 0.0;
180 
181  double locTimeOffset_TAGH = locTOFFitFunc_TAGH->GetParameter(1);
182  if(!(locHist_RFDeltaT_TAGH_TOF->GetEntries() > 0.0))
183  locTimeOffset_TAGH = 0.0;
184 
185  //Fine time offset variances
186  double locTimeOffsetVariance_TOF = 0.0;
187  double locTimeOffsetVariance_FDC = locTOFFitFunc_FDC->GetParError(1)*locTOFFitFunc_FDC->GetParError(1);
188  double locTimeOffsetVariance_PSC = locTOFFitFunc_PSC->GetParError(1)*locTOFFitFunc_PSC->GetParError(1);
189  double locTimeOffsetVariance_TAGH = locTOFFitFunc_TAGH->GetParError(1)*locTOFFitFunc_TAGH->GetParError(1);
190 
191  //Print Fine offsets to screen:
192  cout << "Fine offsets for TOF/TAGH/PSC/FDC are: 0.0, " << locTimeOffset_TAGH << ", " << locTimeOffset_PSC << ", " << locTimeOffset_FDC << endl;
193 
194  //Pipe the current constants into this macro
195  //NOTE: This dumps the "LATEST" values. If you need something else, modify this script.
196  ostringstream locCommandStream;
197  locCommandStream << "ccdb -v " << locVariation << " dump PHOTON_BEAM/RF/time_offset -r " << locRunNumber;
198  FILE* locInputFile = gSystem->OpenPipe(locCommandStream.str().c_str(), "r");
199  if(locInputFile == NULL)
200  return 0;
201 
202  //get the first (comment) line
203  char buff[1024]; // I HATE char buffers
204  if(fgets(buff, sizeof(buff), locInputFile) == NULL)
205  return 0;
206 
207  //get the second (data) line
208  if(fgets(buff, sizeof(buff), locInputFile) == NULL)
209  return 0;
210  istringstream locConstantsStream(buff);
211 
212  //Close the pipe
213  gSystem->ClosePipe(locInputFile);
214 
215  //extract the numbers
216  double locCoarseTimeOffset_TOF = 0.0, locCoarseTimeOffset_TAGH = 0.0, locCoarseTimeOffset_PSC = 0.0, locCoarseTimeOffset_FDC = 0.0;
217  if(!(locConstantsStream >> locCoarseTimeOffset_TOF))
218  return 0;
219  if(!(locConstantsStream >> locCoarseTimeOffset_TAGH))
220  return 0;
221  if(!(locConstantsStream >> locCoarseTimeOffset_PSC))
222  return 0;
223  if(!(locConstantsStream >> locCoarseTimeOffset_FDC))
224  return 0;
225 
226  //Print old Coarse offsets to screen:
227  cout << "Old coarse offsets for TOF/TAGH/PSC/FDC are: " << locCoarseTimeOffset_TOF << ", " << locCoarseTimeOffset_TAGH << ", " << locCoarseTimeOffset_PSC << ", " << locCoarseTimeOffset_FDC << endl;
228 
229  //Compute new time offsets
230  locTimeOffset_TOF += locCoarseTimeOffset_TOF;
231  locTimeOffset_TAGH += locCoarseTimeOffset_TAGH;
232  locTimeOffset_PSC += locCoarseTimeOffset_PSC;
233  locTimeOffset_FDC += locCoarseTimeOffset_FDC;
234 
235  //Print Final offsets to screen:
236  cout << "Final offsets for TOF/TAGH/PSC/FDC are: 0.0, " << locTimeOffset_TAGH << ", " << locTimeOffset_PSC << ", " << locTimeOffset_FDC << endl;
237 
238  //Print Final variances to screen:
239  cout << "Final offset variances for TOF/TAGH/PSC/FDC are: 0.0, " << locTimeOffsetVariance_TAGH << ", " << locTimeOffsetVariance_PSC << ", " << locTimeOffsetVariance_FDC << endl;
240 
241  //Print Fine offsets to file:
242  ofstream locOutputFileStream;
243  locOutputFileStream.open("rf_fine_time_offsets.txt");
244  locOutputFileStream << "0.0 " << std::setprecision(8) << locTimeOffset_TAGH << " " << locTimeOffset_PSC << " " << locTimeOffset_FDC << endl;
245  locOutputFileStream.close();
246 
247  //Print Fine offset variances to file:
248  ofstream locOutputFileStream2;
249  locOutputFileStream2.open("rf_time_offset_vars.txt");
250  locOutputFileStream2 << "0.0 " << std::setprecision(8) << locTimeOffsetVariance_TAGH << " " << locTimeOffsetVariance_PSC << " " << locTimeOffsetVariance_FDC << endl;
251  locOutputFileStream2.close();
252 
253  return 1;
254 }
TH1I * locHist_RFDeltaT_PSC_TOF
TH1I * locHist_RFDeltaT_TAGH_TOF
Double_t Periodic_Gaussian_Func(Double_t *locXArray, Double_t *locParamArray)
char string[256]
int RFMacro_FineTimeOffsets(int locRunNumber, string locVariation="default")
TH1I * locHist_RFDeltaT_PSC_TAGH
TF1 * Create_FitFunc(TH1I *locHist)
TH1I * locHist_RFDeltaT_FDC_TAGH
TDirectory * locDirectory
TH1I * locHist_RFDeltaT_FDC_PSC
TH2D * locHist
double gBeamSignalPeriod
int gRebinTOFAmount
int gRebinF1sAmount
TCanvas * locCanvas
TH1I * locHist_RFDeltaT_FDC_TOF