ATLAS Offline Software
Loading...
Searching...
No Matches
dqutils Namespace Reference

Classes

class  CoolMdt
class  CoolRpc
class  CoolTgc
class  HanOutputFile
class  HistogramDataCOOL
class  MonitoringFile
class  StatusFlagCommentCOOL
class  StatusFlagCOOL
class  StatusFlagCOOLBase

Typedefs

typedef std::map< std::string, std::vector< int > > keycyclemap
typedef std::map< std::string, std::set< std::string > > fileLBMap_t

Enumerations

enum  debugLevel_t { none = 0 , DEBUG , VERBOSE }

Functions

void populateKeyMapping (TDirectory *, keycyclemap &)
std::string getInputDirectory (const std::string &outputDirName, TFile *input, bool has_multiple_runs, std::map< TFile *, std::string > *prefixes)
void getImageBuffer ATLAS_NOT_THREAD_SAFE (TImage **img, TCanvas *myC, char **x, int *y)
void plotResolution (const std::string &coordinate, const std::string &versus)
void plotEfficiency ()
double error_func (float x, const Double_t *par)
double scaleFactorFitFcn (double *x, double *par)
std::vector< float > stableGaussianFit (TH1 *histo)
int updateHists (const std::string &inFileName, const std::string &inStem, const std::string &outFileName="", const std::string &outStem="")
bool makeDirectories (const std::string &dirName)
bool makeDir (const std::string &dirName)
void Copy (TFile *source, TFile *target, const std::string &inDir, const std::string &outDir, const std::string &inHist="", const std::string &outHist="")
void CopyHist (TFile *source, TFile *target, const std::string &inDir, const std::string &outDir, const std::string &inHist, const std::string &outHist)

Variables

std::vector< int > root_color_choices
static const bool fdbg = true
static const bool fpdbg = false
std::atomic< int > padding = 0
static const bool rno_debug = false
static const float TGCChamberLowOccupancyCut = 1.0e-6
static const float TGCChamberHighOccupancyCut = 0.005
static const float TGCChannelOccupancyCut = 0.01
static const float TGCChamberEfficiencyCut = 0.7
static const float TGCChamberTimingCut = 0.95
static const bool tgc_debug = false

Typedef Documentation

◆ fileLBMap_t

typedef std::map<std::string, std::set<std::string> > dqutils::fileLBMap_t

Definition at line 52 of file MonitoringFile.h.

◆ keycyclemap

typedef std::map<std::string, std::vector<int> > dqutils::keycyclemap

Definition at line 51 of file MonitoringFile.h.

Enumeration Type Documentation

◆ debugLevel_t

Enumerator
none 
DEBUG 
VERBOSE 

Definition at line 55 of file MonitoringFile.h.

Function Documentation

◆ ATLAS_NOT_THREAD_SAFE()

void getImageBuffer dqutils::ATLAS_NOT_THREAD_SAFE ( TImage ** img,
TCanvas * myC,
char ** x,
int * y )

Definition at line 1128 of file HanOutputFile.cxx.

1128 {
1129 gVirtualPS->Open(myC->GetName(), 114);
1130 myC->Paint();
1131 auto pImgDump = dynamic_cast<TImageDump*>(gVirtualPS);
1132 if (pImgDump) {
1133 (*img) = pImgDump->GetImage();
1134 if (*img) {
1135 (*img)->GetImageBuffer(x, y, TImage::kPng);
1136 }
1137 }
1138 }
#define y
#define x

◆ Copy()

void dqutils::Copy ( TFile * source,
TFile * target,
const std::string & inDir,
const std::string & outDir,
const std::string & inHist = "",
const std::string & outHist = "" )

Definition at line 247 of file MonitoringFile_MoveVertexMonitoring.cxx.

248 {
249 padding += 3;
250
251 if (!inHist.empty()) {
252 CopyHist(source, target, inDir, outDir, inHist, outHist);
253 } else {
254 TDirectory* sourceDir = source->GetDirectory(inDir.c_str());
255 TDirectory* targetDir = target->GetDirectory(outDir.c_str());
256
257 TIter nextKey(sourceDir->GetListOfKeys());
258
259 TKey* key;
260 while ((key = (TKey*) nextKey())) {
261 std::string keyName = key->GetName();
262 std::string className(key->GetClassName());
263
264 if (className == "TDirectoryFile") {
265 std::string newInDir = inDir;
266 newInDir += '/';
267 newInDir += keyName;
268
269 std::string newOutDir = outDir;
270 newOutDir += '/';
271 newOutDir += keyName;
272
273 if (!targetDir->FindKey(keyName.c_str())) {
274 /*
275 std::cout << std::setw(padding) << " ";
276 std::cout << "creating Dir " << newOutDir << std::endl;
277 */
278 targetDir->mkdir(keyName.c_str());
279 }
280 /*
281 std::cout << std::setw(padding) << " ";
282 std::cout << "moving to " << newInDir << ", a " << className << std::endl;
283 */
284 Copy(source, target, newInDir, newOutDir);
285 } else {
286 CopyHist(source, target, inDir, outDir, keyName, keyName);
287 }
288 }
289 }
290
291 padding -= 3;
292 }
void Copy(TFile *source, TFile *target, const std::string &inDir, const std::string &outDir, const std::string &inHist="", const std::string &outHist="")
void CopyHist(TFile *source, TFile *target, const std::string &inDir, const std::string &outDir, const std::string &inHist, const std::string &outHist)

◆ CopyHist()

void dqutils::CopyHist ( TFile * source,
TFile * target,
const std::string & inDir,
const std::string & outDir,
const std::string & inHist,
const std::string & outHist )

Definition at line 294 of file MonitoringFile_MoveVertexMonitoring.cxx.

295 {
296 TDirectory* sourceDir = source->GetDirectory(inDir.c_str());
297 TDirectory* targetDir = target->GetDirectory(outDir.c_str());
298
299 targetDir->cd();
300 TKey* key = sourceDir->FindKey(inHist.c_str());
301 TObject* object = key->ReadObj();
302 bool permission = true;
303
304 if (targetDir->FindKey(outHist.c_str())) permission = false;
305 if (permission) object->Write(outHist.c_str(), TObject::kOverwrite);
306 }

◆ error_func()

double dqutils::error_func ( float x,
const Double_t * par )

Definition at line 447 of file MonitoringFile_IDEnhancedPrimaryVertex.cxx.

447 {
448//calculating the square of the propagated error on the fit values
449 return(TMath::Power(par[0], 2) + x * TMath::Power(par[1], 2) + TMath::Power(x * par[2], 2));
450 }

◆ getInputDirectory()

std::string dqutils::getInputDirectory ( const std::string & outputDirName,
TFile * input,
bool has_multiple_runs,
std::map< TFile *, std::string > * prefixes )

◆ makeDir()

bool dqutils::makeDir ( const std::string & dirName)

Definition at line 233 of file MonitoringFile_MoveVertexMonitoring.cxx.

233 {
234 padding += 3;
235 std::cout << std::setw(padding) << " ";
236
237 if (!gDirectory->FindKey(dirName.c_str())) {
238 gDirectory->mkdir(dirName.c_str());
239// std::cout << "makeDir=" << dirName << std::endl;
240 } else
241// std::cout << "object " << dirName << " already exists in directory " << gDirectory->GetName() <<
242// std::endl;
243 padding -= 3;
244 return gDirectory->cd(dirName.c_str());
245 }

◆ makeDirectories()

bool dqutils::makeDirectories ( const std::string & dirName)

Definition at line 215 of file MonitoringFile_MoveVertexMonitoring.cxx.

215 {
216 bool success = true;
217
218 if (!dirName.empty()) {
219 std::string::size_type firstSlash = dirName.find('/');
220 if (firstSlash == std::string::npos) {
221 success &= makeDir(dirName);
222 } else {
223 std::string subdir(dirName, 0, firstSlash);
224 if (!subdir.empty()) success &= makeDir(subdir);
225
226 std::string newSubdir(dirName, firstSlash + 1, dirName.size() - firstSlash);
227 success &= makeDirectories(newSubdir);
228 }
229 }
230 return success;
231 }
bool makeDirectories(const std::string &dirName)
bool makeDir(const std::string &dirName)

◆ plotEfficiency()

void dqutils::plotEfficiency ( )

Definition at line 371 of file MonitoringFile_IDEnhancedPrimaryVertex.cxx.

371 {
372// cout << "Creating and writing histos for efficiency" << endl;
373
374 // get histos but return if any histo is not there in input
375 TH1F* h_Vrt_split_tag_ntrk = (TH1F*) gDirectory->Get("Vrt_split_tag_ntrk");
376
377 if (h_Vrt_split_tag_ntrk == 0) return;
378
379 TH1F* h_Vrt_split_probe_ntrk = (TH1F*) gDirectory->Get("Vrt_split_probe_ntrk");
380 if (h_Vrt_split_probe_ntrk == 0) return;
381
382 TH1F* h_Vrt_split_matched_tag_ntrk = (TH1F*) gDirectory->Get("Vrt_split_matched_tag_ntrk");
383 if (h_Vrt_split_matched_tag_ntrk == 0) return;
384
385 TH1F* h_Vrt_split_matched_probe_ntrk = (TH1F*) gDirectory->Get("Vrt_split_matched_probe_ntrk");
386 if (h_Vrt_split_matched_probe_ntrk == 0) return;
387
388 TH1F* h_Vrt_split_dist_tag = (TH1F*) gDirectory->Get("Vrt_split_dist_tag");
389 if (h_Vrt_split_dist_tag == 0) return;
390
391 TH1F* h_Vrt_split_dist_probe = (TH1F*) gDirectory->Get("Vrt_split_dist_probe");
392 if (h_Vrt_split_dist_probe == 0) return;
393
394 // Use BayesDivide routing of TGraphAsymmErrors
395 TGraphAsymmErrors* g_Vrt_rec_eff_m1_split_vs_ntrk = new TGraphAsymmErrors();
396
397 g_Vrt_rec_eff_m1_split_vs_ntrk->BayesDivide(h_Vrt_split_probe_ntrk, h_Vrt_split_tag_ntrk);
398 g_Vrt_rec_eff_m1_split_vs_ntrk->SetName("g_RecEff_M1");
399
400 TGraphAsymmErrors* g_Vrt_sel_eff_m1_split_vs_ntrk = new TGraphAsymmErrors();
401 g_Vrt_sel_eff_m1_split_vs_ntrk->BayesDivide(h_Vrt_split_matched_probe_ntrk, h_Vrt_split_matched_tag_ntrk);
402 g_Vrt_sel_eff_m1_split_vs_ntrk->SetName("g_SelEff_M1");
403
404 // formatting and writing out
405 g_Vrt_rec_eff_m1_split_vs_ntrk->GetHistogram()->GetXaxis()->SetTitle("Number of tracks");
406 g_Vrt_rec_eff_m1_split_vs_ntrk->GetHistogram()->GetYaxis()->SetTitle("Reconstruction efficiency");
407 g_Vrt_rec_eff_m1_split_vs_ntrk->SetMarkerStyle(20);
408 g_Vrt_rec_eff_m1_split_vs_ntrk->Write("", TObject::kOverwrite);
409 delete g_Vrt_rec_eff_m1_split_vs_ntrk;
410
411 g_Vrt_sel_eff_m1_split_vs_ntrk->GetHistogram()->GetXaxis()->SetTitle("Number of tracks");
412 g_Vrt_sel_eff_m1_split_vs_ntrk->GetHistogram()->GetYaxis()->SetTitle("Selection Efficiency");
413 g_Vrt_sel_eff_m1_split_vs_ntrk->SetMarkerStyle(20);
414 g_Vrt_sel_eff_m1_split_vs_ntrk->Write("", TObject::kOverwrite);
415 delete g_Vrt_sel_eff_m1_split_vs_ntrk;
416
417 return;
418 }

◆ plotResolution()

void dqutils::plotResolution ( const std::string & coordinate = "Z",
const std::string & versus = "Ntrk" )

Definition at line 99 of file MonitoringFile_IDEnhancedPrimaryVertex.cxx.

99 {
100// cout << "Creating and writing histos for resolution " << coordinate << " versus " << versus << endl;
101
102 TH2F* h_Vrt_pullVsSomething_split(0);
103 TH2F* h_Vrt_err_vs_Something(0);
104// TH2F* h_Vrt_Tag_err_vs_Something(0);
105 std::string xAxisLabel("");
106
107 if (versus == "Ntrk") {
108 h_Vrt_pullVsSomething_split = (TH2F*) gDirectory->Get(("Vrt_" + coordinate + "pullVsNtrkAverage_split").c_str());
109 h_Vrt_err_vs_Something = (TH2F*) gDirectory->Get(("Vrt_" + coordinate + "err_vs_ntrk").c_str());
110// h_Vrt_Tag_err_vs_Something = (TH2F*)gDirectory->Get("Vrt_Tag_"+coordinate+"err_vs_ntrk");
111 xAxisLabel = "Number of fitted tracks";
112 } else if (versus == "SumPt2") {
113 h_Vrt_pullVsSomething_split = (TH2F*) gDirectory->Get(("Vrt_" + coordinate + "pullVsPt2Average_split").c_str());
114 h_Vrt_err_vs_Something = (TH2F*) gDirectory->Get(("Vrt_" + coordinate + "err_vs_pt2").c_str());
115// h_Vrt_Tag_err_vs_Something = (TH2F*)gDirectory->Get("Vrt_Tag_"+coordinate+"err_vs_pt2");
116 xAxisLabel = "#sqrt{#sum p_{T}^{2}} [GeV]";
117 } else return;
118
119 //if (h_Vrt_pullVsSomething_split == 0) std::cout << "h_Vrt_pullVsSomething_split has zero pointer!" << std::endl;
120 //if (h_Vrt_err_vs_Something == 0) std::cout << "h_Vrt_err_vs_Something has zero pointer!" << std::endl;
121 if (h_Vrt_pullVsSomething_split == 0 or h_Vrt_err_vs_Something == 0) return;
122
123 int n_bins = h_Vrt_pullVsSomething_split->GetNbinsX();
124 std::vector<float> rms_z;
125 std::vector<float> rms_z_er;
126 std::vector<float> sigma_z;
127 std::vector<float> sigma_z_er;
128 std::vector<float> bins_z_nt;
129 std::vector<float> bins_z_nt_er;
130
131// root starts counting the bins at 1, i.e. bin 1 holds NTrk = 0. or sqrt(sumpt2) = 0. - 0.25. GeV
132// std::cout << "n bins: " << n_bins << "\tTotal entries: " << h_Vrt_pullVsSomething_split->GetEntries() << std::endl;
133// TH1D * nTrksPerVertex = h_Vrt_pullVsSomething_split->ProjectionX("projectionNTrks_"+coordinate+"_"+versus);
134
135// TH1D * profileZFull = h_Vrt_pullVsSomething_split->ProjectionY("projectionPullsFull_"+coordinate+"_"+versus);
136
137 Int_t startBin = 0;
138 TH1D* profileZ = 0;
139 const Int_t minEntriesForKFactorBin = 1000;
140 for (int bin_count = 1; bin_count < n_bins + 1; bin_count++) {
141 //Fixed binning
142// TH1D *profileZ = h_Vrt_pullVsSomething_split->ProjectionY("projectionPulls", bin_count, bin_count,"e");
143// Double_t binCenter = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinCenter(bin_count);
144// Double_t binWidth = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinWidth(bin_count)/2.;
145
146 //Variable binning
147 TH1D* profileZTmp = h_Vrt_pullVsSomething_split->ProjectionY("projectionPulls", bin_count, bin_count, "e");
148 //cout << "Bin: " << bin_count << ", Entries: " << profileZTmp->GetEntries() << endl;
149 if (profileZ == 0) {
150 startBin = bin_count;
151 profileZ = (TH1D*) profileZTmp->Clone("projectionPulls_Integrated");
152 //cout << "StartBin = " << startBin << endl;
153 } else {
154 profileZ->Add(profileZTmp);
155 }
156 delete profileZTmp;
157 profileZTmp = 0;
158 if ((profileZ->GetEntries() < minEntriesForKFactorBin) && (bin_count < n_bins)) //not enough entries, not last bin
159 continue;
160
161 Double_t lowEdge = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinLowEdge(startBin);
162 Double_t highEdge = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinLowEdge(bin_count) +
163 h_Vrt_pullVsSomething_split->GetXaxis()->GetBinWidth(bin_count);
164 Double_t binCenter = (lowEdge + highEdge) / 2;
165 Double_t binWidth = (highEdge - lowEdge) / 2; //half of the bin width
166 //cout << "Bin done: " << binCenter << " +/- " << binWidth << ", Entries: " << profileZ->GetEntries() << endl;
167 // variable binning end
168
169 bins_z_nt.push_back(binCenter);
170 bins_z_nt_er.push_back(binWidth); // dummy error of binwidth for now
171
172 rms_z.push_back(profileZ->GetRMS());
173 rms_z_er.push_back(profileZ->GetRMSError());
174
175 //making a gaussian fit if there is anough entries
176 if (profileZ->GetEntries() > 100.) {
177 std::vector<float> fit_res = stableGaussianFit(profileZ);
178 sigma_z.push_back(fit_res[0]);
179 sigma_z_er.push_back(fit_res[1]);
180 } else {
181 sigma_z.push_back(0.);
182 sigma_z_er.push_back(0.);
183 }//end of good number of bins selection
184
185 delete profileZ; // must keep this to delete the projection from memory (next one has same name!)
186 profileZ = 0;
187 }//end of loop over all the ntrk bins
188
189 TGraphErrors* krms_z_vs_ntrk = new TGraphErrors(
190 bins_z_nt.size(), &(bins_z_nt[0]), &(rms_z[0]), &(bins_z_nt_er[0]), &(rms_z_er[0]));
191 krms_z_vs_ntrk->GetYaxis()->SetTitle((coordinate + " scale factor from RMS").c_str());
192 krms_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel.c_str());
193 krms_z_vs_ntrk->SetTitle(("scaleFactor" + coordinate + "_RMS").c_str());
194 krms_z_vs_ntrk->SetName(("scaleFactor" + coordinate + "_" + versus + "_RMS").c_str());
195
196 TGraphErrors* kgs_z_vs_ntrk = new TGraphErrors(
197 bins_z_nt.size(), &(bins_z_nt[0]), &(sigma_z[0]), &(bins_z_nt_er[0]), &(sigma_z_er[0]));
198 kgs_z_vs_ntrk->GetYaxis()->SetTitle((coordinate + " scale factor from gauss fit").c_str());
199 kgs_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel.c_str());
200 kgs_z_vs_ntrk->SetTitle(("scaleFactor" + coordinate + "_Fit").c_str());
201 kgs_z_vs_ntrk->SetName(("scaleFactor_" + coordinate + "_" + versus + "_Fit").c_str());
202
203// approximating the graph with 2nd order polynomial.
204 float maxFitRange(100.);
205 float minFitRange(2.);
206 if (versus == "SumPt2") {
207 minFitRange = 0.5;
208 maxFitRange = 20.;
209 }
210 TF1* kgs_z_ntrk_fit;
211 const Double_t* kgs_z_ntrk_fit_er;
212 int fitResKFactorMethod = 2; // set by hand for now
213 if (fitResKFactorMethod == 1) {
214 //Fit with a pol2
215 //coverity[DEADCODE]
216 kgs_z_vs_ntrk->Fit("pol2", "Q", "", minFitRange, maxFitRange);
217 kgs_z_vs_ntrk->GetFunction("pol2")->SetLineColor(kRed);
218 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction("pol2");
219 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
220 } else if (fitResKFactorMethod == 2) {
221 //Fit with a pol1
222 kgs_z_vs_ntrk->Fit("pol1", "Q", "", minFitRange, maxFitRange);
223 kgs_z_vs_ntrk->GetFunction("pol1")->SetLineColor(kRed);
224 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction("pol1");
225 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
226 //coverity[DEADCODE]
227 } else if (fitResKFactorMethod == 3) {
228 TF1* kgsFitFcn = new TF1("kgsFitFcn", scaleFactorFitFcn, minFitRange, maxFitRange, 3);
229 kgsFitFcn->SetParameter(0, minFitRange);
230 kgsFitFcn->SetParameter(1, 1.0);
231 kgsFitFcn->SetParameter(2, 1.0);
232 for (int ifit = 0; ifit < 1; ifit++) //initial estimation of best parameters
233 kgs_z_vs_ntrk->Fit(kgsFitFcn, "Q");
234 kgs_z_vs_ntrk->Fit(kgsFitFcn, "Q"); //perform actual fit
235 kgs_z_ntrk_fit = kgsFitFcn;
236 kgs_z_ntrk_fit_er = kgsFitFcn->GetParErrors();
237/* cout << "ScaleFactor fit for " << coordinate << " vs " << versus << " (method " << fitResKFactorMethod << ")= "
238 << kgsFitFcn->GetParameter(0) << " +/- " << kgsFitFcn->GetParError(0) << " "
239 << kgsFitFcn->GetParameter(1) << " +/- " << kgsFitFcn->GetParError(1) << " "
240 << kgsFitFcn->GetParameter(2) << " +/- " << kgsFitFcn->GetParError(2) << endl; */
241 } else if (fitResKFactorMethod == 4) {
242 //constant fit
243 kgs_z_vs_ntrk->Fit("pol0", "Q", "", minFitRange, maxFitRange);
244 kgs_z_vs_ntrk->GetFunction("pol0")->SetLineColor(kRed);
245 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction("pol0");
246 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
247/* cout << "ScaleFactor fit for " << coordinate << " vs " << versus << " (method " << fitResKFactorMethod << ")= "
248 << kgs_z_ntrk_fit->GetParameter(0) << " +/- " << kgs_z_ntrk_fit->GetParError(0) << endl; */
249 }
250
251// plotting the fit error of the unconstrained primary vertex and correcting them
252 int nbins_z_err_ntrk = h_Vrt_err_vs_Something->GetNbinsX();
253
254 std::vector<float> av_err_z;
255 std::vector<float> av_err_z_er;
256// std::vector<float> av_err_tag_z;
257// std::vector<float> av_err_tag_z_er;
258 std::vector<float> err_bins_z_nt;
259 std::vector<float> err_bins_z_nt_er;
260 std::vector<float> res_z;
261 std::vector<float> res_z_er;
262// std::vector<float> res_tag_z;
263// std::vector<float> res_tag_z_er;
264
265 for (int bin_count = 1; bin_count <= nbins_z_err_ntrk; ++bin_count) {
266 err_bins_z_nt.push_back(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count));
267 err_bins_z_nt_er.push_back(h_Vrt_err_vs_Something->GetXaxis()->GetBinWidth(bin_count) / 2.);
268
269 TH1D* profileY = h_Vrt_err_vs_Something->ProjectionY("projectionErrors", bin_count, bin_count, "e");
270// TH1D * profileYTag(0);
271// if (h_Vrt_Tag_err_vs_Something)
272// profileYTag = h_Vrt_Tag_err_vs_Something->ProjectionY("projectionErrorsTag",bin_count,bin_count,"e");
273
274 float mean = profileY->GetMean();
275 float mean_error = profileY->GetMeanError();
276// float meanTag(0);
277// float mean_errorTag(0);
278// if (profileYTag) {
279// meanTag = profileYTag->GetMean();
280// mean_errorTag = profileYTag->GetMeanError();
281// }
282 delete profileY;
283// delete profileYTag;
284
285 av_err_z.push_back(mean);
286 av_err_z_er.push_back(mean_error);
287// av_err_tag_z.push_back(meanTag);
288// av_err_tag_z_er.push_back(mean_errorTag);
289
290 //estimating the approximate k-factor and the error value
291 double pr_er = 0.0;
292 float val(0.);
293 if (fitResKFactorMethod == 1) {
294 //coverity[DEADCODE]
295 pr_er = error_func(bin_count, kgs_z_ntrk_fit_er);
296 } else if (fitResKFactorMethod == 2) {
297 val = h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count);
298 pr_er = TMath::Power(kgs_z_ntrk_fit_er[1] * val, 2) + TMath::Power(kgs_z_ntrk_fit_er[0], 2);
299 pr_er = TMath::Sqrt(pr_er);
300// cout << "val = " << val << ", pr_er = " << pr_er << ", p0er = " << kgs_z_ntrk_fit_er[0] << ", p1er = "<<
301// kgs_z_ntrk_fit_er[1] << endl;
302 //coverity[DEADCODE]
303 } else if (fitResKFactorMethod == 3) {
304 val = h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count);
305 //approximately the error on the plateau
306 pr_er = kgs_z_ntrk_fit_er[2];
307 } else if (fitResKFactorMethod == 4) {
308 pr_er = kgs_z_ntrk_fit_er[0];
309 }
310
311 res_z.push_back(mean * kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count)));
312 res_z_er.push_back(TMath::Sqrt(TMath::Power(mean_error *
313 kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(
314 bin_count)),
315 2) + TMath::Power(pr_er * mean, 2)));
316// res_tag_z.push_back(meanTag * kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count)));
317// res_tag_z_er.push_back(TMath::Sqrt(TMath::Power(mean_errorTag *
318// kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count)),2) + TMath::Power( pr_er * mean
319// ,2)));
320 }
321 TGraphErrors* res_z_vs_ntrk =
322 new TGraphErrors(err_bins_z_nt.size(), &(err_bins_z_nt[0]), &(res_z[0]), &(err_bins_z_nt_er[0]), &(res_z_er[0]));
323 res_z_vs_ntrk->GetYaxis()->SetTitle((coordinate + " Vertex Resolution [mm]").c_str());
324 res_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel.c_str());
325 res_z_vs_ntrk->SetTitle((coordinate + " Vertex Resolution").c_str());
326 res_z_vs_ntrk->SetName(("resolution_" + coordinate + "_" + versus).c_str());
327
328// TGraphErrors * res_tag_z_vs_ntrk = new TGraphErrors(err_bins_z_nt.size(),
329// &(err_bins_z_nt[0]),&(res_tag_z[0]),&(err_bins_z_nt_er[0]), &(res_tag_z_er[0]) );
330// res_tag_z_vs_ntrk->GetYaxis()->SetTitle(coordinate+" Tagged Vertex Resolution [mm]");
331// res_tag_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel);
332// res_tag_z_vs_ntrk->SetTitle(coordinate+" Tagged Vertex Resolution");
333// res_tag_z_vs_ntrk->SetName("resolution_tag_"+coordinate+"_"+versus);
334
335// TGraphErrors * mean_err_z_vs_ntrk = new TGraphErrors(err_bins_z_nt.size(),
336// &(err_bins_z_nt[0]),&(av_err_z[0]),&(err_bins_z_nt_er[0]), &(av_err_z_er[0]) );
337// mean_err_z_vs_ntrk->GetYaxis()->SetTitle(coordinate+" mean vertex error [mm]");
338// mean_err_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel);
339// mean_err_z_vs_ntrk->SetTitle(coordinate+" Mean Vertex Error");
340// mean_err_z_vs_ntrk->SetName("resolution_"+coordinate+"_"+versus+"_Uncorrected"); //not corrected with k-factor
341
343 // Writing out
345 if (versus == "Ntrk") res_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
346 else res_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
347 res_z_vs_ntrk->GetYaxis()->SetRangeUser(0.0025, 1.);
348 res_z_vs_ntrk->Write("", TObject::kOverwrite);
349 delete res_z_vs_ntrk;
350
351 if (versus == "Ntrk") krms_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
352 else krms_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
353 krms_z_vs_ntrk->GetYaxis()->SetRangeUser(0.5, 1.3);
354 krms_z_vs_ntrk->Write("", TObject::kOverwrite);
355 delete krms_z_vs_ntrk;
356
357 if (versus == "Ntrk") kgs_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
358 else kgs_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
359 kgs_z_vs_ntrk->GetYaxis()->SetRangeUser(0.5, 1.3);
360 kgs_z_vs_ntrk->Write("", TObject::kOverwrite);
361 delete kgs_z_vs_ntrk;
362
363// nTrksPerVertex->GetYaxis()->SetTitle("Entries");
364// nTrksPerVertex->Draw();
365// nTrksPerVertex->Write("", TObject::kOverwrite);
366
367 return;
368 }
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
void binWidth(TH1 *h)
Definition listroot.cxx:80
std::vector< float > stableGaussianFit(TH1 *histo)
double error_func(float x, const Double_t *par)
double scaleFactorFitFcn(double *x, double *par)

◆ populateKeyMapping()

void dqutils::populateKeyMapping ( TDirectory * ,
keycyclemap &  )

◆ scaleFactorFitFcn()

double dqutils::scaleFactorFitFcn ( double * x,
double * par )

Definition at line 452 of file MonitoringFile_IDEnhancedPrimaryVertex.cxx.

452 {
453 // Truncated gaus-smeared step function
454 // par 0: mean of turn-on (and truncation point)
455 // par 1: slope of turn-on
456 // par 2: plateau
457 return (x[0] >= par[0]) * TMath::Erf(par[1] * x[0] - par[0]) * par[2];
458 }

◆ stableGaussianFit()

std::vector< float > dqutils::stableGaussianFit ( TH1 * histo)

Definition at line 420 of file MonitoringFile_IDEnhancedPrimaryVertex.cxx.

420 {
421 std::vector<float> results;
422 if (histo) {
423 double sigmas = 2.;
424
425 histo->Fit("gaus", "Q0", "", -2, 2);
426 TF1* func = histo->GetFunction("gaus");
427 double actualSigma = func->GetParameter(2);
428 double actualSigmaErr = func->GetParError(2);
429
430 for (int u = 0; u < 10; u++) {
431 histo->Fit("gaus", "Q0", "", -sigmas * actualSigma, sigmas * actualSigma);
432 func = histo->GetFunction("gaus");
433 actualSigma = func->GetParameter(2);
434 actualSigmaErr = func->GetParError(2);
435 }//end of the fitting loop
436
437 results.push_back(actualSigma);
438 results.push_back(actualSigmaErr);
439 } else {
440 results.push_back(0.);
441 results.push_back(0.);
442 }//end of protection against an empty histogram
443
444 return results;
445 }//end of stableGaussian Fit function

◆ updateHists()

int dqutils::updateHists ( const std::string & inFileName,
const std::string & inStem,
const std::string & outFileName = "",
const std::string & outStem = "" )

Definition at line 143 of file MonitoringFile_MoveVertexMonitoring.cxx.

144 {
145 std::string outFileName = fileName;
146 std::string outStem = stem;
147 //open original file
148 TFile* source = TFile::Open(inFileName.c_str());
149 if (!source) {
150 //std::cout << "Couldn't open input file, " << inFileName << std::endl;
151 return 1;
152 }//else std::cout << "opening input file, " << inFileName << std::endl;
153
154 //find out whether inStem is a histogram or a directory
155 bool isDir = true;
156 std::string path = inStem;
157 std::string hist;
158
159 if (inStem[inStem.size() - 1] != '/') {
160 std::string::size_type lastSlash = inStem.find_last_of('/');
161 hist = inStem.substr(lastSlash + 1, inStem.size() - lastSlash - 1);
162 path = inStem.substr(0, lastSlash + 1);
163 isDir = (source->FindObjectAny(hist.c_str()))->InheritsFrom("TDirectory");
164 if (isDir) {
165 path = inStem;
166 hist.clear();
167 }
168 }
169
170 if (path[path.size() - 1] == '/') path.resize(path.size() - 1);
171 if (!source->GetDirectory(path.c_str())) {
172 //std::cout << "path " << path << " directory doesn't exist in input file" << std::endl;
173 return 1;
174 }
175
176 //open target file
177 std::cout.fill(' ');
178 if (outFileName.empty()) {
179 outFileName = std::string(inFileName, 0, inFileName.find(".root"));
180 outFileName += "_trimmed.root";
181 }
182
183 TFile* target = TFile::Open(outFileName.c_str(), "update");
184 if (!target) {
185 //std::cout << "coundn't open output file " << outFileName << std::endl;
186 return 1;
187 }//else std::cout << "opening output file " << outFileName << std::endl;
188
189 if (outStem.empty()) outStem = inStem;
190 std::string targetPath = outStem;
191 std::string targetHist;
192 if (!isDir) {
193 std::string::size_type lastSlash = outStem.find_last_of('/');
194 targetPath = outStem.substr(0, lastSlash + 1);
195 targetHist = outStem.substr(lastSlash + 1, outStem.size() - lastSlash - 1);
196 }
197
198 if (targetPath[targetPath.size() - 1] == '/') targetPath.resize(targetPath.size() - 1);
199
200 target->cd();
201 if (!makeDirectories(targetPath)) {
202 //std::cout << "couldn't create outStem directories in output" << std::endl;
203 return 1;
204 }
205 //copy relevant objects
206 if (!target->IsWritable()) return 1;
207
208 Copy(source, target, path, targetPath, hist, targetHist);
209
210 delete target;
211 delete source;
212 return 0;
213 }
path
python interpreter configuration --------------------------------------—
Definition athena.py:126

Variable Documentation

◆ fdbg

const bool dqutils::fdbg = true
static

Definition at line 40 of file MonitoringFile_HLTMuonHistogramDivision.cxx.

◆ fpdbg

const bool dqutils::fpdbg = false
static

Definition at line 38 of file MonitoringFile_HLTMuonPostProcess.cxx.

◆ padding

std::atomic<int> dqutils::padding = 0

Definition at line 16 of file MonitoringFile_MoveVertexMonitoring.cxx.

◆ rno_debug

const bool dqutils::rno_debug = false
static

Definition at line 38 of file MonitoringFile_PixelPostProcess.cxx.

◆ root_color_choices

std::vector<int> dqutils::root_color_choices
Initial value:
= {
kBlue, kRed, kGray, kOrange, kViolet, kGreen + 1
}

Definition at line 84 of file HanOutputFile.cxx.

84 {
85 kBlue, kRed, kGray, kOrange, kViolet, kGreen + 1
86 };

◆ tgc_debug

const bool dqutils::tgc_debug = false
static

Definition at line 43 of file MonitoringFile_TGCPostProcess.cxx.

◆ TGCChamberEfficiencyCut

const float dqutils::TGCChamberEfficiencyCut = 0.7
static

Definition at line 40 of file MonitoringFile_TGCPostProcess.cxx.

◆ TGCChamberHighOccupancyCut

const float dqutils::TGCChamberHighOccupancyCut = 0.005
static

Definition at line 38 of file MonitoringFile_TGCPostProcess.cxx.

◆ TGCChamberLowOccupancyCut

const float dqutils::TGCChamberLowOccupancyCut = 1.0e-6
static

Definition at line 37 of file MonitoringFile_TGCPostProcess.cxx.

◆ TGCChamberTimingCut

const float dqutils::TGCChamberTimingCut = 0.95
static

Definition at line 41 of file MonitoringFile_TGCPostProcess.cxx.

◆ TGCChannelOccupancyCut

const float dqutils::TGCChannelOccupancyCut = 0.01
static

Definition at line 39 of file MonitoringFile_TGCPostProcess.cxx.