Script to plot changes between two calib bin files.

This commit is contained in:
redford_s
2019-02-26 12:21:07 +01:00
parent 9bf8924dc8
commit 2a9a583d17
3 changed files with 135 additions and 1 deletions

1
.gitignore vendored
View File

@ -11,6 +11,7 @@ JFMC_CalibWriter_wBP
JFMC_BadPixels
LP_analysis
LC_analysis
JFMC_CalibComp
data/
plots/

130
JFMC_CalibComp.cpp Normal file
View File

@ -0,0 +1,130 @@
#include "../sls_detector_calibration/jungfrauCommonHeader.h"
#include "../sls_detector_calibration/jungfrauCommonFunctions.h"
#include <fstream>
#include "TH1F.h"
#include "TCanvas.h"
#include "TPaveText.h"
struct GainMaps
{
double g0VALs[NCH];
double g1VALs[NCH];
double g2VALs[NCH];
double hg0VALs[NCH];
};
int main(int argc, char* argv[]) {
jungfrauStyle();
if (argc != 3) {
cout << "Correct usage:" << endl;
cout << "arg 1: specify 1st bin file to compare" << endl;
cout << "arg 2: specify 2nd bin file to compare" << endl;
cout << " " << endl;
exit(1);
}
string file1_str = argv[1];
string file2_str = argv[2];
char savename[128];
GainMaps *mapsObject1 = new GainMaps;
fstream infile1;
sprintf(savename,"%s",file1_str.c_str());
infile1.open((char*)savename, ios::in | ios::binary);
if (not infile1.is_open()) {
cout << "file 1 doesn't exist" << endl;
exit(1);
}
infile1.read((char*)mapsObject1->g0VALs, sizeof(mapsObject1->g0VALs));
infile1.read((char*)mapsObject1->g1VALs, sizeof(mapsObject1->g1VALs));
infile1.read((char*)mapsObject1->g2VALs, sizeof(mapsObject1->g2VALs));
infile1.read((char*)mapsObject1->hg0VALs, sizeof(mapsObject1->hg0VALs));
infile1.close();
GainMaps *mapsObject2 = new GainMaps;
fstream infile2;
sprintf(savename,"%s",file2_str.c_str());
infile2.open((char*)savename, ios::in | ios::binary);
if (not infile2.is_open()) {
cout << "file 2 doesn't exist" << endl;
exit(1);
}
infile2.read((char*)mapsObject2->g0VALs, sizeof(mapsObject2->g0VALs));
infile2.read((char*)mapsObject2->g1VALs, sizeof(mapsObject2->g1VALs));
infile2.read((char*)mapsObject2->g2VALs, sizeof(mapsObject2->g2VALs));
infile2.read((char*)mapsObject2->hg0VALs, sizeof(mapsObject2->hg0VALs));
infile2.close();
TH1F* g0_diffpc_hist = new TH1F("g0_diffpc_hist","",100,-5-0.05,5-0.05);
TH1F* g1_diffpc_hist = new TH1F("g1_diffpc_hist","",100,-5-0.05,5-0.05);
TH1F* g2_diffpc_hist = new TH1F("g2_diffpc_hist","",100,-5-0.05,5-0.05);
TH1F* hg0_diffpc_hist = new TH1F("hg0_diffpc_hist","",100,-5-0.05,5-0.05);
for (int i = 0; i < NCH; i++) {
g0_diffpc_hist->Fill(100.*(mapsObject2->g0VALs[i] - mapsObject1->g0VALs[i])/mapsObject2->g0VALs[i]);
g1_diffpc_hist->Fill(100.*(mapsObject2->g1VALs[i] - mapsObject1->g1VALs[i])/mapsObject2->g1VALs[i]);
g2_diffpc_hist->Fill(100.*(mapsObject2->g2VALs[i] - mapsObject1->g2VALs[i])/mapsObject2->g2VALs[i]);
hg0_diffpc_hist->Fill(100.*(mapsObject2->hg0VALs[i] - mapsObject1->hg0VALs[i])/mapsObject2->hg0VALs[i]);
}
TCanvas* c1 = new TCanvas("c1","");
TPaveText *pave = new TPaveText(0.86,0.97,0.91,0.99,"blNDC");
pave->SetBorderSize(0);
pave->SetFillStyle(0);
pave->SetTextSize(0.06);
pave->SetTextAlign(32);
g0_diffpc_hist->GetXaxis()->SetTitle("% change G0");
g0_diffpc_hist->Draw();
pave->Clear();
sprintf(savename,"Outliers %0.0f",g0_diffpc_hist->GetBinContent(0)+g0_diffpc_hist->GetBinContent(101));
pave->AddText(savename);
pave->Draw();
c1->SaveAs("plots/CalibComp_g0.png");
c1->SetLogy();
c1->SaveAs("plots/CalibComp_g0_log.png");
c1->SetLogy(0);
g1_diffpc_hist->GetXaxis()->SetTitle("% change G1");
g1_diffpc_hist->Draw();
pave->Clear();
sprintf(savename,"Outliers %0.0f",g1_diffpc_hist->GetBinContent(0)+g1_diffpc_hist->GetBinContent(101));
pave->AddText(savename);
pave->Draw();
c1->SaveAs("plots/CalibComp_g1.png");
c1->SetLogy();
c1->SaveAs("plots/CalibComp_g1_log.png");
c1->SetLogy(0);
g2_diffpc_hist->GetXaxis()->SetTitle("% change G2");
g2_diffpc_hist->Draw();
pave->Clear();
sprintf(savename,"Outliers %0.0f",g2_diffpc_hist->GetBinContent(0)+g2_diffpc_hist->GetBinContent(101));
pave->AddText(savename);
pave->Draw();
c1->SaveAs("plots/CalibComp_g2.png");
c1->SetLogy();
c1->SaveAs("plots/CalibComp_g2_log.png");
c1->SetLogy(0);
hg0_diffpc_hist->GetXaxis()->SetTitle("% change HG0");
hg0_diffpc_hist->Draw();
pave->Clear();
sprintf(savename,"Outliers %0.0f",hg0_diffpc_hist->GetBinContent(0)+hg0_diffpc_hist->GetBinContent(101));
pave->AddText(savename);
pave->Draw();
c1->SaveAs("plots/CalibComp_hg0.png");
c1->SetLogy();
c1->SaveAs("plots/CalibComp_hg0_log.png");
c1->SetLogy(0);
}

View File

@ -33,4 +33,7 @@ LP_analysis: LP_analysis.cpp
g++ -Wall -O3 -m64 -I$(ROOTSYS)/include -L$(ROOTSYS)/lib -lGui -lCore -lCint -lRIO -lNet -lHist -lGraf -lGraf3d -lGpad -lTree -lRint -lPostscript -lMatrix -lPhysics -lMathCore -lThread -pthread -lm -ldl -rdynamic LP_analysis.cpp -o LP_analysis
LC_analysis: LC_analysis.cpp
g++ -Wall -O3 -m64 -I$(ROOTSYS)/include -L$(ROOTSYS)/lib -lGui -lCore -lCint -lRIO -lNet -lHist -lGraf -lGraf3d -lGpad -lTree -lRint -lPostscript -lMatrix -lPhysics -lMathCore -lThread -pthread -lm -ldl -rdynamic LC_analysis.cpp -o LC_analysis
g++ -Wall -O3 -m64 -I$(ROOTSYS)/include -L$(ROOTSYS)/lib -lGui -lCore -lCint -lRIO -lNet -lHist -lGraf -lGraf3d -lGpad -lTree -lRint -lPostscript -lMatrix -lPhysics -lMathCore -lThread -pthread -lm -ldl -rdynamic LC_analysis.cpp -o LC_analysis
JFMC_CalibComp: JFMC_CalibComp.cpp
g++ -Wall -O3 -m64 -I$(ROOTSYS)/include -L$(ROOTSYS)/lib -lGui -lCore -lCint -lRIO -lNet -lHist -lGraf -lGraf3d -lGpad -lTree -lRint -lPostscript -lMatrix -lPhysics -lMathCore -lThread -pthread -lm -ldl -rdynamic JFMC_CalibComp.cpp -o JFMC_CalibComp