Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFMacro_CoarseTimeOffsets.C
Go to the documentation of this file.
1 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_FDC_TOF
2 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_FDC_TAGH
3 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_FDC_PSC
4 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_PSC_TAGH
5 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_PSC_TOF
6 // hnamepath: /RF/AbsoluteDeltaT_RF_OtherRFs/RFDeltaT_TAGH_TOF
7 
8 double Calc_Mean(TH1I* locHist)
9 {
10  double locTotal = 0.0;
11  Int_t locEntries = 0;
12  for(Int_t loc_i = 1; loc_i <= locHist->GetNbinsX(); ++loc_i)
13  {
14  Int_t locThisBinEntries = locHist->GetBinContent(loc_i);
15  locEntries += locThisBinEntries;
16  locTotal += double(locThisBinEntries)*locHist->GetBinCenter(loc_i);
17  }
18  if(locEntries == 0)
19  return 0.0;
20  return (locTotal / double(locEntries));
21 }
22 
23 int RFMacro_CoarseTimeOffsets(int locRunNumber, string locVariation = "default")
24 {
25  //INPUT RUN NUMBER MUST BE A RUN NUMBER IN THE RUN RANGE YOU ARE TRYING TO COMMIT CONSTANTS TO
26  gDirectory->cd("/"); //return to file base directory
27 
28  //Goto Beam Path
29  TDirectory *locDirectory = (TDirectory*)gDirectory->FindObjectAny("RF");
30  if(!locDirectory)
31  return 0;
32  locDirectory->cd();
33 
34  //Get RF DeltaT Histograms
35  gDirectory->cd("AbsoluteDeltaT_RF_OtherRFs");
36  TH1I* locHist_RFDeltaT_FDC_TOF = (TH1I*)gDirectory->Get("RFDeltaT_FDC_TOF");
37  TH1I* locHist_RFDeltaT_FDC_TAGH = (TH1I*)gDirectory->Get("RFDeltaT_FDC_TAGH");
38  TH1I* locHist_RFDeltaT_FDC_PSC = (TH1I*)gDirectory->Get("RFDeltaT_FDC_PSC");
39  TH1I* locHist_RFDeltaT_PSC_TAGH = (TH1I*)gDirectory->Get("RFDeltaT_PSC_TAGH");
40  TH1I* locHist_RFDeltaT_PSC_TOF = (TH1I*)gDirectory->Get("RFDeltaT_PSC_TOF");
41  TH1I* locHist_RFDeltaT_TAGH_TOF = (TH1I*)gDirectory->Get("RFDeltaT_TAGH_TOF");
42 
43  //New coarse time offsets (with respect to TOF)
44  double locTimeOffset_TOF = 0.0;
45  double locTimeOffset_FDC = Calc_Mean(locHist_RFDeltaT_FDC_TOF);
46  double locTimeOffset_PSC = Calc_Mean(locHist_RFDeltaT_PSC_TOF);
47  double locTimeOffset_TAGH = Calc_Mean(locHist_RFDeltaT_TAGH_TOF);
48 
49  //Print extracted coarse offsets to screen:
50  cout << "Extracted offsets for TOF/TAGH/PSC/FDC are: 0.0, " << locTimeOffset_TAGH << ", " << locTimeOffset_PSC << ", " << locTimeOffset_FDC << endl;
51 
52  //Pipe the current constants into this macro
53  //NOTE: This dumps the "LATEST" values. If you need something else, modify this script.
54  ostringstream locCommandStream;
55  locCommandStream << "ccdb -v " << locVariation << " dump PHOTON_BEAM/RF/time_offset -r " << locRunNumber;
56  FILE* locInputFile = gSystem->OpenPipe(locCommandStream.str().c_str(), "r");
57  if(locInputFile == NULL)
58  return 0;
59 
60  //get the first (comment) line
61  char buff[1024]; // I HATE char buffers
62  if(fgets(buff, sizeof(buff), locInputFile) == NULL)
63  return 0;
64 
65  //get the second (data) line
66  if(fgets(buff, sizeof(buff), locInputFile) == NULL)
67  return 0;
68  istringstream locConstantsStream(buff);
69 
70  //Close the pipe
71  gSystem->ClosePipe(locInputFile);
72 
73  //extract the numbers
74  double locOldTimeOffset_TOF = 0.0, locOldTimeOffset_TAGH = 0.0, locOldTimeOffset_PSC = 0.0, locOldTimeOffset_FDC = 0.0;
75  if(!(locConstantsStream >> locOldTimeOffset_TOF))
76  return 0;
77  if(!(locConstantsStream >> locOldTimeOffset_TAGH))
78  return 0;
79  if(!(locConstantsStream >> locOldTimeOffset_PSC))
80  return 0;
81  if(!(locConstantsStream >> locOldTimeOffset_FDC))
82  return 0;
83 
84  //Print old time offsets to screen:
85  cout << "Old time offsets for TOF/TAGH/PSC/FDC are: " << locOldTimeOffset_TOF << ", " << locOldTimeOffset_TAGH << ", " << locOldTimeOffset_PSC << ", " << locOldTimeOffset_FDC << endl;
86 
87  //Compute new time offsets
88  locTimeOffset_TOF += locOldTimeOffset_TOF;
89  locTimeOffset_TAGH += locOldTimeOffset_TAGH;
90  locTimeOffset_PSC += locOldTimeOffset_PSC;
91  locTimeOffset_FDC += locOldTimeOffset_FDC;
92 
93  //Print new coarse time offsets to screen:
94  cout << "New coarse time offsets for TOF/TAGH/PSC/FDC are: " << locTimeOffset_TOF << ", " << locTimeOffset_TAGH << ", " << locTimeOffset_PSC << ", " << locTimeOffset_FDC << endl;
95 
96  //Print new coarse offsets to file:
97  ofstream locOutputFileStream;
98  locOutputFileStream.open("rf_coarse_time_offsets.txt");
99  locOutputFileStream << "0.0 " << std::setprecision(8) << locTimeOffset_TAGH << " " << locTimeOffset_PSC << " " << locTimeOffset_FDC << endl;
100  locOutputFileStream.close();
101 
102  //Get/Make Canvas
103  TCanvas *locCanvas = NULL;
104  if(TVirtualPad::Pad() == NULL)
105  locCanvas = new TCanvas("RF_CoarseTimeOffsets", "RF_CoarseTimeOffsets", 1200, 800); //for testing
106  else
107  locCanvas = gPad->GetCanvas();
108  locCanvas->Divide(3, 2);
109 
110  //Draw
111  locCanvas->cd(1);
112  gPad->SetTicks();
113  gPad->SetGrid();
114  if(locHist_RFDeltaT_TAGH_TOF != NULL)
115  {
116  locHist_RFDeltaT_TAGH_TOF->GetXaxis()->SetTitleSize(0.05);
117  locHist_RFDeltaT_TAGH_TOF->GetYaxis()->SetTitleSize(0.05);
118  locHist_RFDeltaT_TAGH_TOF->GetYaxis()->SetLabelSize(0.05);
119  locHist_RFDeltaT_TAGH_TOF->Draw();
120  }
121 
122  locCanvas->cd(2);
123  gPad->SetTicks();
124  gPad->SetGrid();
125  if(locHist_RFDeltaT_PSC_TOF != NULL)
126  {
127  locHist_RFDeltaT_PSC_TOF->GetXaxis()->SetTitleSize(0.05);
128  locHist_RFDeltaT_PSC_TOF->GetYaxis()->SetTitleSize(0.05);
129  locHist_RFDeltaT_PSC_TOF->GetYaxis()->SetLabelSize(0.05);
130  locHist_RFDeltaT_PSC_TOF->Draw();
131  }
132 
133  locCanvas->cd(3);
134  gPad->SetTicks();
135  gPad->SetGrid();
136  if(locHist_RFDeltaT_FDC_TOF != NULL)
137  {
138  locHist_RFDeltaT_FDC_TOF->GetXaxis()->SetTitleSize(0.05);
139  locHist_RFDeltaT_FDC_TOF->GetYaxis()->SetTitleSize(0.05);
140  locHist_RFDeltaT_FDC_TOF->GetYaxis()->SetLabelSize(0.05);
141  locHist_RFDeltaT_FDC_TOF->Draw();
142  }
143 
144  locCanvas->cd(4);
145  gPad->SetTicks();
146  gPad->SetGrid();
147  if(locHist_RFDeltaT_FDC_TAGH != NULL)
148  {
149  locHist_RFDeltaT_FDC_TAGH->GetXaxis()->SetTitleSize(0.05);
150  locHist_RFDeltaT_FDC_TAGH->GetYaxis()->SetTitleSize(0.05);
151  locHist_RFDeltaT_FDC_TAGH->GetYaxis()->SetLabelSize(0.05);
152  locHist_RFDeltaT_FDC_TAGH->Draw();
153  }
154 
155  locCanvas->cd(5);
156  gPad->SetTicks();
157  gPad->SetGrid();
158  if(locHist_RFDeltaT_FDC_PSC != NULL)
159  {
160  locHist_RFDeltaT_FDC_PSC->GetXaxis()->SetTitleSize(0.05);
161  locHist_RFDeltaT_FDC_PSC->GetYaxis()->SetTitleSize(0.05);
162  locHist_RFDeltaT_FDC_PSC->GetYaxis()->SetLabelSize(0.05);
163  locHist_RFDeltaT_FDC_PSC->Draw();
164  }
165 
166  locCanvas->cd(6);
167  gPad->SetTicks();
168  gPad->SetGrid();
169  if(locHist_RFDeltaT_PSC_TAGH != NULL)
170  {
171  locHist_RFDeltaT_PSC_TAGH->GetXaxis()->SetTitleSize(0.05);
172  locHist_RFDeltaT_PSC_TAGH->GetYaxis()->SetTitleSize(0.05);
173  locHist_RFDeltaT_PSC_TAGH->GetYaxis()->SetLabelSize(0.05);
174  locHist_RFDeltaT_PSC_TAGH->Draw();
175  }
176 
177  return 1;
178 }
TH1I * locHist_RFDeltaT_PSC_TOF
TH1I * locHist_RFDeltaT_TAGH_TOF
TH1I * locHist_RFDeltaT_PSC_TAGH
int RFMacro_CoarseTimeOffsets(int locRunNumber, string locVariation="default")
TH1I * locHist_RFDeltaT_FDC_TAGH
TDirectory * locDirectory
TH1I * locHist_RFDeltaT_FDC_PSC
TH2D * locHist
double Calc_Mean(TH1I *locHist)
TCanvas * locCanvas
TH1I * locHist_RFDeltaT_FDC_TOF