BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
ElectronInterface.cc
Go to the documentation of this file.
1#include "ElectronInterface.h"
2
4
5 // if no files are specified, assume sample is MC (1) and not data (0)
6 m_mcFlag = 1;
7
8 // if no files are specified, assume sample is BESIII (0) and not Belle II (1)
9 m_type = 0;
10
11 // initialize various parameters
12 m_filename = "";
13 m_constfilename = "";
14 m_nevents = 0;
15
16 // default truncation
17 m_removeLowest = 0.05;
18 m_removeHighest = 0.25;
19
20 // specify whether to take truncated means (0) or perform fits (1)
21 m_fits = 0;
22}
23
24void
26
27 cout << "ElectronInterface: setting up from " << configfile << endl;
28
29 std::ifstream input(configfile.c_str());
30
31 std::string varname, vartype, line;
32 while( std::getline(input, line) ){
33
34 // skip comment lines
35 line.erase( std::find( line.begin(), line.end(), '#' ), line.end() );
36 if( line.empty() ) continue;
37
38 // parse the lines in the configuration file
39 std::istringstream iss(line);
40 iss >> vartype;
41
42 // check for a mcFlag
43 if( vartype == "mcFlag" ){
44 iss >> m_mcFlag;
45 continue;
46 }
47 else if( vartype == "type" ){
48 iss >> m_type;
49 continue;
50 }
51 else if( vartype == "fit" ){
52 iss >> m_fits;
53 continue;
54 }
55 else if( vartype == "constfile" ){
56 iss >> m_constfilename;
57 continue;
58 }
59 else if( vartype == "lowest" ){
60 iss >> m_removeLowest;
61 continue;
62 }
63 else if( vartype == "highest" ){
64 iss >> m_removeHighest;
65 continue;
66 }
67 else if( vartype != "electron" ){ // the first word in the line should be the particle type
68 std::cout << "ElectronInterface: ERROR: unknown particle type " << vartype << " - this line will not be read" << endl;
69 continue;
70 }
71
72 // the next word should be the variable type
73 iss >> varname;
74 if( varname == "filename")
75 iss >> m_filename;
76 else if( varname == "nevents")
77 iss >> m_nevents;
78 else if( varname == "docabins")
79 iss >> m_docabins;
80 else if( varname == "upperdoca")
81 iss >> m_upperdoca;
82 else if( varname == "lowerdoca")
83 iss >> m_lowerdoca;
84 else if( varname == "bgbins")
85 iss >> m_bgbins;
86 else if( varname == "upperp"){
87 iss >> m_upperbg;
88 m_upperbg = m_upperbg/Widget::melectron;
89 }
90 else if( varname == "lowerp"){
91 iss >> m_lowerbg;
92 m_lowerbg = m_lowerbg/Widget::melectron;
93 }
94 else if( varname == "entabins")
95 iss >> m_entabins;
96 else if( varname == "upperenta")
97 iss >> m_upperenta;
98 else if( varname == "lowerenta")
99 iss >> m_lowerenta;
100 else if( varname == "costhbins")
101 iss >> m_costhbins;
102 else{
103 std::cout << "ElectronInterface: ERROR: unknown variable type " << vartype << " - this line will not be read" << endl;
104 continue;
105 }
106 }
107
108 std::cout << "ElectronInterface: Done parsing configuration file..." << endl;
109}
110
111
112void
113ElectronInterface::SetFileName( TString newfilename ){
114 m_filename = newfilename;
115}
116
117
118void
120
121 std::cout << "ElectronInterface: generating electron sample..." << std::endl;
122
123 TFile* outfile = new TFile(outfilename,"RECREATE");
124
125 // bin the data in run numbers and fit for the mean in each bin
126 ElectronGenerator egen(100000,2.0,0.2);
127 egen.generateEvents(outfile);
128
129 // if we are generating a sample, use it for calibration
130 m_filename = outfilename;
131}
132
133
134
135// -------------------- ELECTRON CORRECTIONS --------------------
136
137void
138ElectronInterface::RunGains( TString outfilename ) {
139
140 std::cout << "ElectronInterface: getting run gains..." << std::endl;
141
142 TFile* outfile = new TFile(outfilename,"RECREATE");
143 // make sure the samples exist
144 TFile* testfile = new TFile(m_filename);
145 if( testfile->GetListOfKeys()->IsEmpty() ){
146 cout << "\t No electron sample" << endl;
147 exit;
148 }
149 testfile->Close();
150
151 // bin the data in run numbers and fit for the mean in each bin
152 ElectronCalibration ecalib(m_filename, m_constfilename, m_mcFlag, m_type, m_fits, m_docabins, m_upperdoca, m_lowerdoca, m_entabins, m_upperenta, m_lowerenta, m_costhbins);
153 ecalib.fitRunGains(outfile);
154
155 // plot the run gains
156 ecalib.plotRunGains(outfilename);
157}
158
159void
160ElectronInterface::TwoDCorrection( TString outfilename ) {
161
162 std::cout << "WidgetInterface: getting two dimensional correction..." << std::endl;
163
164 TFile* outfile = new TFile(outfilename,"RECREATE");
165 // make sure the samples exist
166 TFile* testfile = new TFile(m_filename);
167 if( testfile->GetListOfKeys()->IsEmpty() ){
168 cout << "\t No electron sample" << endl;
169 exit;
170 }
171 testfile->Close();
172
173 // bin the data in doca and entrance angle and fit for the mean in each bin
174 ElectronCalibration ecalib(m_filename, m_constfilename, m_mcFlag, m_type, m_fits, m_docabins, m_upperdoca, m_lowerdoca, m_entabins, m_upperenta, m_lowerenta, m_costhbins);
175 ecalib.TwoDCorrection(outfile);
176
177 outfile->Close();
178}
179
180void
182
183 std::cout << "WidgetInterface: doing electron saturation correction..." << std::endl;
184
185 TFile* outfile = new TFile(outfilename,"RECREATE");
186 // make sure the samples exist
187 TFile* testfile = new TFile(m_filename);
188 if( testfile->GetListOfKeys()->IsEmpty() ){
189 cout << "\t No electron sample" << endl;
190 exit;
191 }
192 testfile->Close();
193
194 // bin the data in doca and entrance angle and fit for the mean in each bin
195 ElectronCalibration ecalib(m_filename, m_constfilename, m_mcFlag, m_type, m_fits, m_docabins, m_upperdoca, m_lowerdoca, m_entabins, m_upperenta, m_lowerenta, m_costhbins);
196 ecalib.SaturationCorrection(outfile);
197
198 outfile->Close();
199}
200
201void
203
204 std::cout << "WidgetInterface: applying electron corrections..." << std::endl;
205
206 TFile* outfile = new TFile(outfilename,"RECREATE");
207 // make sure the samples exist
208 TFile* testfile = new TFile(m_filename);
209 if( testfile->GetListOfKeys()->IsEmpty() ){
210 cout << "\t No electron sample" << endl;
211 exit;
212 }
213 testfile->Close();
214
215 // bin the data in doca and entrance angle and fit for the mean in each bin
216 ElectronCorrection ecor(m_constfilename,m_filename,m_removeLowest,m_removeHighest);
217 ecor.process(outfile);
218
219 outfile->Close();
220}
void generateEvents(TFile *outfile)
void SetupFromConfigFile(std::string configfile)
void RunGains(TString outfilename)
void ApplyCorrections(TString outfile)
void SaturationCorrection(TString outfilename)
void SetFileName(TString newfilename)
void TwoDCorrection(TString outfilename)
void GenerateSample(TString outfile)