Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:43

0001 #include <string>
0002 #include <vector>
0003 #include "TF1.h"
0004 #include "TH1D.h"
0005 #include "TH1F.h"
0006 #include "TH2.h"
0007 #include "TCanvas.h"
0008 #include "TCutG.h"
0009 #include "TFile.h"
0010 #include "TLegend.h"
0011 #include "TMath.h"
0012 #include "TNtuple.h"
0013 #include "TPad.h"
0014 #include "TPaveText.h"
0015 #include "TROOT.h"
0016 #include "TString.h"
0017 #include "TSystem.h"
0018 #include <MuonAnalysis/MomentumScaleCalibration/test/Macros/RooFit/TkAlStyle.cc>
0019 
0020 
0021 using namespace ROOT::Math;
0022 
0023 struct Plot_1D{
0024   TString directory;
0025   TString histoname;
0026   TString plotname;
0027   TString xtitle;
0028   TString ytitle;
0029   float xmin;
0030   float xmax;
0031   float ymin;
0032   float ymax;
0033   float absMin;
0034   float absMax;
0035 
0036   vector<TH1D*> histogram;
0037   vector<TString> histogram_label;
0038   vector<pair<float, float>> mean_sigma_Y;
0039   vector<pair<int, float>> chi2Y;
0040 
0041   Plot_1D(TString directory_, TString histoname_, TString plotname_, TString xtitle_, TString ytitle_, float xmin_, float xmax_, float ymin_, float ymax_) :
0042     directory(directory_),
0043     histoname(histoname_),
0044     plotname(plotname_),
0045     xtitle(xtitle_),
0046     ytitle(ytitle_),
0047     xmin(xmin_),
0048     xmax(xmax_),
0049     ymin(ymin_),
0050     ymax(ymax_),
0051     absMin(9e9),
0052     absMax(-9e9)
0053   {}
0054   ~Plot_1D(){}
0055 
0056   void addMeanSigmaChi2Y(TH1D* htmp){
0057     double mean=0;
0058     double sigma=0;
0059     for (int bin=1; bin<=htmp->GetNbinsX(); bin++){
0060       float bincontent = htmp->GetBinContent(bin);
0061       float binerror = htmp->GetBinError(bin);
0062       if (binerror==0. && bincontent==0.) continue;
0063       absMin = min(absMin, bincontent - binerror);
0064       absMax = max(absMax, bincontent + binerror);
0065       mean += bincontent/pow(binerror, 2);
0066       sigma += 1./pow(binerror, 2);
0067     }
0068     mean /= sigma;
0069     sigma = sqrt(1./sigma);
0070     mean_sigma_Y.push_back(pair<float, float>(mean, sigma));
0071 
0072     float chi2=0; int ndof=-1;
0073     for (int bin=1; bin<=htmp->GetNbinsX(); bin++){
0074       float bincontent = htmp->GetBinContent(bin);
0075       float binerror = htmp->GetBinError(bin);
0076       if (binerror==0. && bincontent==0.) continue;
0077       ndof++;
0078       chi2 += pow((bincontent-mean)/binerror, 2);
0079     }
0080     chi2Y.push_back(pair<int, float>(ndof, chi2));
0081   }
0082 
0083   void addHistogramFromFile(TFile* finput, TString label, int color, int linestyle, int markerstyle){
0084     if (finput!=0 && !finput->IsZombie()){
0085       TH1D* htmp = (TH1D*)finput->Get(histoname);
0086       if (htmp!=0){
0087         htmp->SetTitle("");
0088         htmp->GetXaxis()->SetTitle(xtitle);
0089         htmp->GetYaxis()->SetTitle(ytitle);
0090 
0091         htmp->SetLineWidth(1);
0092         htmp->SetLineColor(color);
0093         htmp->SetLineStyle(linestyle);
0094         htmp->SetMarkerSize(1.2);
0095         htmp->SetMarkerColor(color);
0096         htmp->SetMarkerStyle(markerstyle);
0097 
0098         histogram.push_back(htmp);
0099         histogram_label.push_back(label);
0100         addMeanSigmaChi2Y(htmp);
0101       }
0102     }
0103   }
0104 
0105   void plot(TString fitFormula="", bool AutoSetRange=false){
0106     // Determine the plot ranges
0107     unsigned int nValidHistos = histogram.size();
0108     const unsigned int nValidHistosLimit = 4;
0109     float rangeFactor[2]={ 1.01, 1.05 };
0110     float dampingFactor = 0;
0111     float deviationThreshold = 1.04;
0112     if (nValidHistos>nValidHistosLimit) dampingFactor = 0.07/nValidHistosLimit*(nValidHistos-nValidHistosLimit);
0113     double rangeMaxReduction = 0.02;
0114     if (nValidHistos>nValidHistosLimit) rangeMaxReduction = rangeMaxReduction*nValidHistosLimit/nValidHistos;
0115     double dampingFactorEff = dampingFactor;
0116 
0117     float overallAvg=0;
0118     float overallSigmaSqInv=0;
0119     for (unsigned int f=0; f<nValidHistos; f++){
0120       float avgY = mean_sigma_Y.at(f).first;
0121       float devY = mean_sigma_Y.at(f).second;
0122       devY = 1./pow(devY, 2);
0123       avgY *= devY;
0124       overallAvg += avgY;
0125       overallSigmaSqInv += devY;
0126     }
0127     overallAvg /= overallSigmaSqInv;
0128 
0129     // Update August 2018, implemented autorange for histograms
0130     float MaxY=-1000.;
0131     float MinY=1000.;
0132 
0133 
0134     for (unsigned int f=0; f<nValidHistos; f++){
0135       for (int bin=1; bin<=histogram.at(f)->GetNbinsX(); bin++){
0136         float bincontent = histogram.at(f)->GetBinContent(bin);
0137         float binerror = histogram.at(f)->GetBinError(bin);
0138         if (binerror==0 && bincontent==0) continue;
0139         if ((bincontent + binerror)>deviationThreshold*overallAvg) rangeMaxReduction = 0;
0140     if(bincontent>=MaxY)MaxY=bincontent;
0141     if(bincontent<=MinY)MinY=bincontent;
0142 
0143       }
0144     }
0145     if (nValidHistos>nValidHistosLimit && rangeMaxReduction!=0) dampingFactorEff = dampingFactorEff*0.7;
0146     
0147     float Yrange = MaxY-MinY; 
0148     if(AutoSetRange){ // Overwtite fixed range for Y axis with the one computed for the current histogram
0149       ymin=MinY-Yrange/10.;
0150       ymax=Yrange*1.6+ymin;
0151     }
0152 
0153    
0154     if (ymin>=ymax){ // If range is not already set
0155       ymin = absMin/rangeFactor[0];
0156       ymax = absMax*(rangeFactor[1]+dampingFactorEff-rangeMaxReduction);
0157     }
0158     for (unsigned int f=0; f<nValidHistos; f++){
0159       histogram.at(f)->GetXaxis()->SetRangeUser(xmin, xmax);
0160       histogram.at(f)->GetYaxis()->SetRangeUser(ymin, ymax);
0161     }
0162 
0163     // Begin plotting
0164     vector<TF1*> hfit;
0165     TCanvas* canvas;
0166     TLegend* legend;
0167     TString cname = Form("c%s", plotname.Data());
0168     canvas = new TCanvas(cname, cname, 8, 30, 800, 800);
0169     canvas->cd();
0170     legend = TkAlStyle::legend("topleft", nValidHistos);
0171     for (unsigned int f=0; f<nValidHistos; f++){
0172       TString legendlabel = histogram_label.at(f);
0173       if (f==0) histogram.at(f)->Draw();
0174       else histogram.at(f)->Draw("same");
0175 
0176       if (fitFormula==""){
0177         if (plotname.Contains("mean")){
0178           float normchi2=0.;
0179           if (chi2Y.at(f).first>0) normchi2=chi2Y.at(f).first/chi2Y.at(f).second;
0180           legendlabel = Form("%s (#chi^{2}/dof=%.2f)", legendlabel.Data(), normchi2);
0181         }
0182       }
0183       else{
0184         TF1* ftmp = new TF1(Form("fit_%i", f), fitFormula, xmin, xmax);
0185         ftmp->SetParameter(0, 91.0);
0186         ftmp->SetParameter(1, 0);
0187         if (fitFormula.Contains("[2]")) ftmp->SetParameter(2, 0);
0188         ftmp->SetLineColor(1);
0189 
0190         histogram.at(f)->Fit(ftmp, "R", "same", xmin, xmax);
0191         ftmp->Draw("same");
0192         hfit.push_back(ftmp);
0193       }
0194 
0195       legend->AddEntry(histogram.at(f), legendlabel, "lp");
0196     }
0197     legend->Draw("same");
0198     TkAlStyle::drawStandardTitle();
0199     canvas->RedrawAxis();
0200     canvas->Modified();
0201     canvas->Update();
0202     canvas->SaveAs(Form("%s/%s%s", directory.Data(), plotname.Data(), ".png"));
0203     canvas->SaveAs(Form("%s/%s%s", directory.Data(), plotname.Data(), ".pdf"));
0204     canvas->SaveAs(Form("%s/%s%s", directory.Data(), plotname.Data(), ".root"));
0205 
0206     for (unsigned int hf=0; hf<hfit.size(); hf++) delete hfit.at(hf);
0207     delete legend;
0208     canvas->Close();
0209     delete canvas;
0210   }
0211 
0212 };
0213 
0214 
0215 void splitOption(string rawoption, string& wish, string& value, char delimiter){
0216   size_t posEq = rawoption.find(delimiter);
0217   if (posEq!=string::npos){
0218     wish=rawoption;
0219     value=rawoption.substr(posEq+1);
0220     wish.erase(wish.begin()+posEq, wish.end());
0221   }
0222   else{
0223     wish="";
0224     value=rawoption;
0225   }
0226 }
0227 void splitOptionRecursive(string rawoption, vector<string>& splitoptions, char delimiter){
0228   string suboption=rawoption, result=rawoption;
0229   string remnant;
0230   while (result!=""){
0231     splitOption(suboption, result, remnant, delimiter);
0232     if (result!="") splitoptions.push_back(result);
0233     suboption = remnant;
0234   }
0235   if (remnant!="") splitoptions.push_back(remnant);
0236 }
0237 
0238 void getCustomRanges(TString type, TString resonance, int iP, float minmax_plot[2]);
0239 void MultiHistoOverlapAll_Base_one(string files, string labels, string colors, string linestyles, string markerstyles, TString directory, TString resonance, TString type, bool switchONfit, bool AutoSetRange=false, float CustomMinY=90.85, float CustomMaxY=91.4);
0240 
0241 void MultiHistoOverlapAll_Base(string files, string labels, string colors, string linestyles, string markerstyles, TString directory, TString resonance, bool switchONfit=false, bool AutoSetRange=false, float CustomMinY=90.85, float CustomMaxY=91.4){
0242   gSystem->mkdir(directory, true);
0243   MultiHistoOverlapAll_Base_one(files, labels, colors, linestyles, markerstyles, directory, resonance, "mean", switchONfit, AutoSetRange, CustomMinY, CustomMaxY);
0244   MultiHistoOverlapAll_Base_one(files, labels, colors, linestyles, markerstyles, directory, resonance, "sigma", switchONfit, AutoSetRange, CustomMinY, CustomMaxY);
0245 }
0246 
0247 void MultiHistoOverlapAll_Base_one(string files, string labels, string colors, string linestyles, string markerstyles, TString directory, TString resonance, TString type, bool switchONfit, bool AutoSetRange, float CustomMinY, float CustomMaxY){
0248   gROOT->Reset();
0249   if (TkAlStyle::status() == NO_STATUS) TkAlStyle::set(INTERNAL);
0250   gROOT->ForceStyle();
0251 
0252   vector<string> strValidation_file;
0253   vector<string> strValidation_label;
0254   vector<string> strValidation_color;
0255   vector<string> strValidation_linestyle;
0256   vector<string> strValidation_markerstyle;
0257   splitOptionRecursive(files, strValidation_file, ',');
0258   splitOptionRecursive(labels, strValidation_label, ',');
0259   splitOptionRecursive(colors, strValidation_color, ',');
0260   splitOptionRecursive(linestyles, strValidation_linestyle, ',');
0261   splitOptionRecursive(markerstyles, strValidation_markerstyle, ',');
0262   int nfiles = strValidation_file.size();
0263   int nlabels = strValidation_label.size();
0264   int ncolors = strValidation_color.size();
0265   int nlinestyles = strValidation_linestyle.size();
0266   int nmarkerstyles = strValidation_markerstyle.size();
0267   if (nlabels!=nfiles){
0268     cout << "nlabels!=nfiles" << endl;
0269     return;
0270   }
0271   if (ncolors!=0 && ncolors!=nfiles){
0272     cout << "ncolors!=nfiles" << endl;
0273     return;
0274   }
0275   if (nlinestyles!=0 && nlinestyles!=nfiles){
0276     cout << "nlinestyles!=nfiles" << endl;
0277     return;
0278   }
0279   if (nmarkerstyles!=0 && nmarkerstyles!=nfiles){
0280     cout << "nmarkerstyles!=nfiles" << endl;
0281     return;
0282   }
0283 
0284   const int nhistos=8;
0285   const int n2Dhistos = 2;
0286   TFile** file = new TFile*[nfiles];
0287   for (int f=0; f<nfiles; f++) file[f] = TFile::Open((strValidation_file[f]).c_str(), "read");
0288 
0289   TString histoname[nhistos] ={
0290     TString("MassVsPt/allHistos/") + type + TString("Histo"),
0291     TString("MassVsPhiPlus/allHistos/") + type + TString("Histo"),
0292     TString("MassVsPhiMinus/allHistos/") + type + TString("Histo"),
0293     TString("MassVsEtaPlus/allHistos/") + type + TString("Histo"),
0294     TString("MassVsEtaMinus/allHistos/") + type + TString("Histo"),
0295     TString("MassVsEtaPlusMinusDiff/allHistos/") + type + TString("Histo"),
0296     TString("MassVsCosThetaCS/allHistos/") + type + TString("Histo"),
0297     TString("MassVsPhiCS/allHistos/") + type + TString("Histo")
0298   };
0299   TString histo2Dname[n2Dhistos] ={
0300     TString("MassVsEtaPhiPlus/allHistos/") + type + TString("Histo"),
0301     TString("MassVsEtaPhiMinus/allHistos/") + type + TString("Histo"),
0302   };
0303   TString xtitle[nhistos] ={
0304     "p^{T}_{#mu} (GeV)",
0305     "#phi_{#mu+}",
0306     "#phi_{#mu-}",
0307     "#eta_{#mu+}",
0308     "#eta_{#mu-}",
0309     "#eta_{#mu+} - #eta_{#mu-}",
0310     "cos #theta_{CS}",
0311     "#phi_{CS}"
0312   };
0313   TString x2Dtitle[n2Dhistos] ={
0314     "#phi_{#mu+}",
0315     "#phi_{#mu-}"
0316   };
0317   TString y2Dtitle[n2Dhistos] ={
0318     "#eta_{#mu+}",
0319     "#eta_{#mu-}"
0320   };
0321   TString ytitle;
0322   if (type=="mean") ytitle = "M_{#mu#mu} (GeV)";
0323   else if (type=="sigma") ytitle = "#sigma_{#mu#mu} (GeV)";
0324   TString plotname[nhistos] ={
0325     type + "MassVsPt_ALL",
0326     type + "MassVsPhiPlus_ALL",
0327     type + "MassVsPhiMinus_ALL",
0328     type + "MassVsEtaPlus_ALL",
0329     type + "MassVsEtaMinus_ALL",
0330     type + "MassVsDeltaEta_ALL",
0331     type + "MassVsCosThetaCS_ALL",
0332     type + "MassVsPhiCS_ALL"
0333   };
0334   TString plotname2D[n2Dhistos] ={
0335     type + "MassVsEtaPhiPlus",
0336     type + "MassVsEtaPhiMinus"
0337   };
0338   TString fitFormula[nhistos]={
0339     "[0]+[1]*cos(x+[2])",
0340     "[0]+[1]*cos(x+[2])",
0341     "[0]+[1]*cos(x+[2])",
0342     "[0]+[1]*x",
0343     "[0]+[1]*x",
0344     "[0]+[1]*x",
0345     "[0]+[1]*cos(x+[2])",
0346     "[0]+[1]*cos(x+[2])"
0347   };
0348   float plot_xmax[nhistos]={
0349     static_cast<float>((resonance=="Z" ? 100. : (resonance=="Y1S" ? 20. : (resonance=="Y2S" ? 20. : (resonance=="Y3S" ? 20. : (resonance=="JPsi" ? 20. : 20.)))))),
0350     static_cast<float>(TMath::Pi()),
0351     static_cast<float>(TMath::Pi()),
0352     2.4,
0353     2.4,
0354     4.8,
0355     1,
0356     static_cast<float>(TMath::Pi())
0357   };
0358   float plot_xmin[nhistos]={
0359     0.,
0360     static_cast<float>(-TMath::Pi()),
0361     static_cast<float>(-TMath::Pi()),
0362     -2.4,
0363     -2.4,
0364     -4.8,
0365     -1,
0366     static_cast<float>(-TMath::Pi())
0367   };
0368 
0369 
0370   // Plot 2D hystograms
0371   // TO BE FIXED: something strange appears in the 2D plots for the sigma, this should be further investigated
0372   Double_t zMass(91.1876);
0373   Double_t Y1SMass(9.46);
0374   Double_t maxDist(0.);
0375   Double_t Mass=zMass;
0376   Double_t fixrange=1.5;
0377   if (resonance=="Y1S"){
0378     Mass=Y1SMass;
0379     fixrange=0.15;
0380   }
0381   for (int f=0; f<nfiles; f++){
0382  
0383     if (file[f]!=0 && !file[f]->IsZombie()){
0384 
0385       
0386       TH2D* histo2D[2];
0387       histo2D[0]= (TH2D*)file[f]->Get(histo2Dname[0]); //histoMassVsEtaPhiPlus
0388       histo2D[1]= (TH2D*)file[f]->Get(histo2Dname[1]); //histoMassVsEtaPhiMinus
0389       for(int s=0;s<2;s++){
0390     TCanvas dummycanvas;
0391     dummycanvas.SetFillColor(0);  
0392     dummycanvas.cd()->SetTopMargin(0.07);
0393     dummycanvas.cd()->SetRightMargin(0.15);
0394     dummycanvas.cd()->SetLeftMargin(0.12);
0395     histo2D[s]->SetTitle("");
0396     histo2D[s]->GetZaxis()->SetLabelFont(42);
0397     histo2D[s]->GetXaxis()->SetTitle(x2Dtitle[s]);
0398     histo2D[s]->GetYaxis()->SetTitle(y2Dtitle[s]);
0399     histo2D[s]->GetYaxis()->SetTitleOffset(0.9);
0400     if(type=="mean"){//This is made to ensure that the range options are not used for the sigma
0401       Double_t zMin=Mass;
0402       Double_t zMax=Mass;
0403  
0404       for(int nbinX=1;nbinX<=histo2D[s]->GetNbinsX();nbinX++){
0405         for(int nbinY=1;nbinY<=histo2D[s]->GetNbinsY();nbinY++){
0406           Double_t value = histo2D[s]->GetBinContent(nbinX,nbinY);
0407           if(value<zMin)zMin=value;
0408           if(value>zMax)zMax=value;
0409         }
0410       }
0411       maxDist=fabs(zMax-Mass);
0412       if(fabs(Mass-zMin) >= maxDist) maxDist=fabs(Mass-zMin);
0413       histo2D[s]->GetZaxis()->SetRangeUser(Mass-fixrange,Mass+fixrange);//Default range Zmass +- 1.5 GeV (Ymass +- 0.15 GeV), NOTE: this will create empty bins when the bin content is either lower or higher than the fixed range
0414       if(AutoSetRange)histo2D[s]->GetZaxis()->SetRangeUser(zMin,zMax); //Set range automatically
0415       
0416     }
0417     histo2D[s]->Draw("COLZ");
0418     TString alignment_label = strValidation_label.at(f);
0419     alignment_label.ReplaceAll(" ","_");
0420     dummycanvas.SaveAs(Form("%s/%sc%s%s%s", directory.Data(),alignment_label.Data(),"_", plotname2D[s].Data(), ".png"));
0421     dummycanvas.SaveAs(Form("%s/%sc%s%s%s", directory.Data(),alignment_label.Data(),"_", plotname2D[s].Data(), ".pdf"));
0422     dummycanvas.SaveAs(Form("%s/%sc%s%s%s", directory.Data(),alignment_label.Data(),"_", plotname2D[s].Data(), ".root"));
0423       }
0424       
0425     }
0426 
0427   };
0428   
0429 
0430   for (int iP=0; iP<nhistos; iP++){
0431     TString theFitFormula="";
0432     if (switchONfit && type=="mean") theFitFormula = fitFormula[iP];
0433 
0434     float minmax_plot[2]={ -1, -1 };
0435     //getCustomRanges(type, resonance, iP, minmax_plot);
0436     if (type=="mean") {
0437       minmax_plot[0]=CustomMinY;
0438       minmax_plot[1]=CustomMaxY;
0439     }
0440     Plot_1D thePlot(directory, histoname[iP], plotname[iP], xtitle[iP], ytitle, plot_xmin[iP], plot_xmax[iP], minmax_plot[0], minmax_plot[1]);
0441     for (int f=0; f<nfiles; f++){
0442       int color=0, linestyle=0, markerstyle=0;
0443       if (ncolors>f) color = stoi(strValidation_color.at(f));
0444       else if (f==0) color = (int)kBlack;
0445       else if (f==(nfiles-1)) color = (int)kViolet;
0446       else if (f==1) color = (int)kBlue;
0447       else if (f==2) color = (int)kRed;
0448       else if (f==3) color = (int)(kGreen+2);
0449       else if (f==4) color = (int)(kOrange+3);
0450       else if (f==5) color = (int)kGreen;
0451       else if (f==6) color = (int)kYellow;
0452       else if (f==7) color = (int)(kPink+9);
0453       else if (f==8) color = (int)kCyan;
0454       else if (f==9) color = (int)(kGreen+3);
0455       else color = (int)kBlack;
0456       if (nlinestyles>f){
0457         linestyle = stoi(strValidation_linestyle.at(f));
0458         if (linestyle > 100) markerstyle = (linestyle / 100);
0459         linestyle = linestyle % 100;
0460       }
0461       else linestyle=1;
0462       if (nmarkerstyles>f) markerstyle = stoi(strValidation_markerstyle.at(f));
0463       else if (markerstyle==0){
0464         if (strValidation_label.at(f).find("reference")!=string::npos || strValidation_label.at(f).find("Reference")!=string::npos) markerstyle=1;
0465         else markerstyle=20;
0466       }
0467       thePlot.addHistogramFromFile(file[f], strValidation_label.at(f), color, linestyle, markerstyle);
0468     }
0469     if(type=="mean")thePlot.plot(theFitFormula,AutoSetRange);
0470     else thePlot.plot(theFitFormula,true);
0471     
0472   }
0473   
0474   for (int f=nfiles-1; f>=0; f--){ if (file[f]!=0 && file[f]->IsOpen()) file[f]->Close(); }
0475   delete[] file;
0476 }
0477 
0478 void getCustomRanges(TString type, TString resonance, int iP, float minmax_plot[2]){
0479   if (type=="mean"){ // Custom ranges
0480 
0481     if (resonance=="Z"){
0482       if (iP==7){ // PhiCS
0483         minmax_plot[0] = 90.9;
0484         minmax_plot[1] = 91.3;
0485       }
0486       else if (iP==6){ // cosThetaCS
0487         minmax_plot[0] = 90.95;
0488         minmax_plot[1] = 91.4;
0489       }
0490       else if (iP==5){ // deta
0491         minmax_plot[0] = 90.9;
0492         minmax_plot[1] = 91.6;
0493       }
0494       else if (iP==4){ // eta-
0495         minmax_plot[0] = 90.8;
0496         minmax_plot[1] = 91.5;
0497       }
0498       else if (iP==3){ // eta+
0499         minmax_plot[0] = 90.8;
0500         minmax_plot[1] = 91.5;
0501       }
0502       else if (iP==2){ // phi-
0503         minmax_plot[0] = 90.8;
0504         minmax_plot[1] = 91.7;
0505       }
0506       else if (iP==1){ // phi+
0507         minmax_plot[0] = 90.8;
0508         minmax_plot[1] = 91.7;
0509       }
0510       else{ // pt
0511         minmax_plot[0] = 90.;
0512         minmax_plot[1] = 94.;
0513       }
0514     }
0515     else if (resonance=="Y1S"){
0516       if (iP==7){ // PhiCS
0517         minmax_plot[0] = 9.42;
0518         minmax_plot[1] = 9.50;
0519       }
0520       else if (iP==6){ // cosThetaCS
0521         minmax_plot[0] = 9.43;
0522         minmax_plot[1] = 9.49;
0523       }
0524       else if (iP==5){ // deta
0525         minmax_plot[0] = 9.43;
0526         minmax_plot[1] = 9.49;
0527       }
0528       else if (iP==4){ // eta-
0529         minmax_plot[0] = 9.40;
0530         minmax_plot[1] = 9.52;
0531       }
0532       else if (iP==3){ // eta+
0533         minmax_plot[0] = 9.41;
0534         minmax_plot[1] = 9.52;
0535       }
0536       else if (iP==2){ // phi-
0537         minmax_plot[0] = 9.435;
0538         minmax_plot[1] = 9.475;
0539       }
0540       else if (iP==1){ // phi+
0541         minmax_plot[0] = 9.435;
0542         minmax_plot[1] = 9.475;
0543       }
0544       else{ // pt
0545         minmax_plot[0] = 9.4;
0546         minmax_plot[1] = 9.65;
0547       }
0548     }
0549 
0550   }
0551   else if (type=="sigma"){
0552 
0553     if (resonance=="Z"){
0554       if (iP==7){ // PhiCS
0555         minmax_plot[0] = 1.1;
0556         minmax_plot[1] = 1.7;
0557       }
0558       else if (iP==6){ // cosThetaCS
0559         minmax_plot[0] = 1.0;
0560         minmax_plot[1] = 1.65;
0561       }
0562       else if (iP==5){ // deta
0563         minmax_plot[0] = 0.8;
0564         minmax_plot[1] = 2.3;
0565       }
0566       else if (iP==4){ // eta-
0567         minmax_plot[0] = 0.8;
0568         minmax_plot[1] = 3.2;
0569       }
0570       else if (iP==3){ // eta+
0571         minmax_plot[0] = 0.8;
0572         minmax_plot[1] = 3.2;
0573       }
0574       else if (iP==2){ // phi-
0575         minmax_plot[0] = 1.0;
0576         minmax_plot[1] = 1.7;
0577       }
0578       else if (iP==1){ // phi+
0579         minmax_plot[0] = 1.0;
0580         minmax_plot[1] = 1.7;
0581       }
0582       else{ // pt
0583         minmax_plot[0] = 0;
0584         minmax_plot[1] = 3.0;
0585       }
0586     }
0587     else if (resonance=="Y1S"){
0588       if (iP==7){ // PhiCS
0589         minmax_plot[0] = 0.05;
0590         minmax_plot[1] = 0.16;
0591       }
0592       else if (iP==6){ // cosThetaCS
0593         minmax_plot[0] = 0.06;
0594         minmax_plot[1] = 0.12;
0595       }
0596       else if (iP==5){ // deta
0597         minmax_plot[0] = 0.04;
0598         minmax_plot[1] = 0.16;
0599       }
0600       else if (iP==4){ // eta-
0601         minmax_plot[0] = 0.05;
0602         minmax_plot[1] = 0.22;
0603       }
0604       else if (iP==3){ // eta+
0605         minmax_plot[0] = 0.05;
0606         minmax_plot[1] = 0.22;
0607       }
0608       else if (iP==2){ // phi-
0609         minmax_plot[0] = 0.06;
0610         minmax_plot[1] = 0.11;
0611       }
0612       else if (iP==1){ // phi+
0613         minmax_plot[0] = 0.06;
0614         minmax_plot[1] = 0.11;
0615       }
0616       else{ // pt
0617         minmax_plot[0] = 0.0;
0618         minmax_plot[1] = 0.26;
0619       }
0620     }
0621 
0622   } // End custom ranges
0623 
0624 }