Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-30 09:14:58

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #include "Mat_map_detector_plot.C"
0010 
0011 /// Draw and save the ratio plots.
0012 
0013 void plot_ratio(std::vector<TH2F*> Map_prop, std::vector<TH2F*> Map_geant, std::vector<int> detectors, const std::string& name){
0014 
0015 
0016   TH2F *Unit_Map_prop = (TH2F*) Map_prop[2]->Clone();
0017   Unit_Map_prop->Divide(Map_prop[2]);
0018   TH2F *Unit_Map_geant = (TH2F*) Map_geant[2]->Clone();
0019   Unit_Map_geant->Divide(Map_geant[2]);
0020 
0021   TH1D *Proj_eta_prop = (TH1D*) Map_prop[0]->ProjectionX()->Clone();
0022   Proj_eta_prop->Divide(Unit_Map_prop->ProjectionX());
0023   TH1D *Proj_eta_geant = (TH1D*) Map_geant[0]->ProjectionX()->Clone();
0024   Proj_eta_geant->Divide(Unit_Map_geant->ProjectionX());
0025 
0026   TH1D *Proj_phi_prop = (TH1D*) Map_prop[0]->ProjectionY()->Clone();
0027   Proj_phi_prop->Divide(Unit_Map_prop->ProjectionY());
0028   TH1D *Proj_phi_geant = (TH1D*) Map_geant[0]->ProjectionY()->Clone();
0029   Proj_phi_geant->Divide(Unit_Map_geant->ProjectionY());
0030 
0031   std::string sVol = "Detector volumes :";
0032   for(auto const& det: detectors) {
0033     sVol += " ";
0034     sVol += std::to_string(det);
0035   }
0036   TText *vol = new TText(.1, .95, sVol.c_str());
0037   vol->SetNDC();
0038 
0039   TCanvas *c1 = new TCanvas("c1","mat_X0",1200,1200);
0040   c1->SetRightMargin(0.14);
0041   c1->SetTopMargin(0.14);
0042   c1->SetLeftMargin(0.14);
0043   c1->SetBottomMargin(0.14);
0044   Map_prop[0]->Divide(Map_geant[0]);
0045   Map_prop[0]->GetZaxis()->SetTitle("X0 Val/Geant");
0046   Map_prop[0]->SetMaximum(2.);
0047   Map_prop[0]->Draw("COLZ");
0048   vol->Draw();
0049   c1->Print( (name+"ratio_X0.pdf").c_str());
0050 
0051   TCanvas *c2 = new TCanvas("c2", "mat_X0/eta", 1200, 1200);
0052   c2->SetRightMargin(0.14);
0053   c2->SetTopMargin(0.14);
0054   c2->SetLeftMargin(0.14);
0055   c2->SetBottomMargin(0.14);
0056 
0057   Proj_eta_prop->Divide(Proj_eta_geant);
0058   Proj_eta_prop->GetYaxis()->SetTitle("X0 Val/Geant");
0059   Proj_eta_prop->SetMarkerStyle(7);
0060   Proj_eta_prop->Draw("HIST PC");
0061   c2->Print((name + "ratio_X0_eta.pdf").c_str());
0062 
0063   TCanvas *c3 = new TCanvas("c3", "mat_X0/phi", 1200, 1200);
0064   c3->SetRightMargin(0.14);
0065   c3->SetTopMargin(0.14);
0066   c3->SetLeftMargin(0.14);
0067   c3->SetBottomMargin(0.14);
0068   Proj_phi_prop->Divide(Proj_phi_geant);
0069   Proj_phi_prop->GetYaxis()->SetTitle("X0 Val/Geant");
0070   Proj_phi_prop->SetMarkerStyle(7);
0071   Proj_phi_prop->Draw("HIST PC");
0072   c3->Print((name + "ratio_X0_phi.pdf").c_str());
0073 
0074   delete c1;
0075   delete c2;
0076   delete c3;
0077   delete vol;
0078   delete Unit_Map_prop;
0079   delete Unit_Map_geant;
0080 }
0081 
0082 
0083 /// Plot the material ratio between the geantino scan and the map validation for each detector.
0084 /// detectors : list of the ID of the volume constitutive of the detector/sub-detector
0085 /// nbprocess : number of parameter to be processed
0086 /// name : name of the output directory.
0087 /// name_prop : name of the output directory for the map valdation.
0088 /// name_geant : name of the output directory for the geantino scan.
0089 /// The map valdation and geantino scan plots are only saved if name_prop and name_geant are defined.
0090 
0091 void Mat_map_detector_plot_ratio(std::string input_file_prop = "", std::string input_file_geant = "", std::vector<int> detectors = vector<int>(), int nbprocess = -1, std::string name = "", std::string name_prop = "", std::string name_geant = ""){
0092 
0093   gStyle->SetOptStat(0);
0094   gStyle->SetOptTitle(0);
0095 
0096   std::vector<TH2F*> detector_hist_prop;
0097   std::vector<TH2F*> detector_hist_geant;
0098 
0099   Fill(detector_hist_prop, input_file_prop, detectors, nbprocess);
0100   Fill(detector_hist_geant, input_file_geant, detectors, nbprocess);
0101 
0102   plot(detector_hist_prop, detectors, name_prop);
0103   plot(detector_hist_geant, detectors, name_geant);
0104   plot_ratio(detector_hist_prop, detector_hist_geant, detectors, name);
0105 }