3 #include "TEfficiency.h"
9 #include "StMuDSTMaker/COMMON/StMuDst.h"
10 #include "StMuDSTMaker/COMMON/StMuMcTrack.h"
11 #include "StMuDSTMaker/COMMON/StMuMcVertex.h"
12 #include "StMuDSTMaker/COMMON/StMuPrimaryVertex.h"
18 tvx::HistContainer(name, motherDir, option)
20 Add(
new TH1I(
"hNumVertices",
"; Number of Vertices; Counts; ", 10, 0, 10));
21 Add(
new TH1I(
"hNumTracksPerEvent",
"; Number of Tracks per Event; Counts", 100, 0, 100));
22 Add(
new TH1I(
"hNumTracksPerEventPrimary",
"; Number of Primary Tracks per Event; Counts", 100, 0, 100));
25 Add(
new TH1I(
"McRecMulT",
"; Reconstructable multiplicity for trigger Mc Vertex; Num. of Reco Vertices; ", 50, 0, 500) );
26 Add(
new TH1I(
"McRecMulAny",
"; MC Vertex Multiplicity; Num. of Reco Vertices matching MC (idTruth); ", 50, 0, 500) );
27 Add(
new TH1I(
"McRecMulGood",
"; MC Vertex Multiplicity; Num. of Max Rank Reco Vertices matching MC (idTruth); ", 50, 0, 500) );
29 Add(
new TH2I(
"hMcMultVsQATruth",
"; MC Vertex Multiplicity; QA Truth; Counts", 50, 0, 50, 50, 0, 1) );
35 h(
"hNumVertices")->Fill( event.primaryVertices()->GetEntriesFast() );
36 h(
"hNumTracksPerEvent")->Fill( event.globalTracks()->GetEntriesFast() );
37 h(
"hNumTracksPerEventPrimary")->Fill( event.primaryTracks()->GetEntriesFast() );
42 const StMuPrimaryVertex *recoVertex,
const StMuPrimaryVertex *recoVertexMaxRank)
47 auto isReconstructible = [](
const StMuMcTrack* mcTrack) {
return mcTrack->No_tpc_hit() >= 15; };
48 int nRecoTracks = std::count_if(mcTracks.begin(), mcTracks.end(), isReconstructible);
51 h(
"McRecMulT")->Fill(nRecoTracks);
54 h(
"McRecMulAny")->Fill(nRecoTracks);
55 h(
"hMcMultVsQATruth")->Fill(nRecoTracks, recoVertex->qaTruth()/100.);
57 if (recoVertex == recoVertexMaxRank)
58 h(
"McRecMulGood")->Fill(nRecoTracks);
66 TEfficiency *eff_any =
new TEfficiency(*h(
"McRecMulAny"), *h(
"McRecMulT"));
68 std::string eff_title = std::string(
";") + h(
"McRecMulAny")->GetXaxis()->GetTitle() +
"; Efficiency";
70 eff_any->SetTitle(eff_title.c_str() );
71 eff_any->SetMarkerColor(kBlack);
72 eff_any->SetMarkerStyle(kFullCross);
74 TH1 *eff_hist =
static_cast<TH1*
>( h(
"McRecMulAny")->Clone() );
75 eff_hist->Reset(
"ICES");
76 eff_hist->SetName( (std::string(eff_hist->GetName()) +
"_eff").c_str() );
77 eff_hist->SetYTitle(
"Efficiency");
78 eff_hist->GetListOfFunctions()->Add(eff_any);
84 TEfficiency *eff_good =
new TEfficiency(*h(
"McRecMulGood"), *h(
"McRecMulT"));
86 eff_title = std::string(
";") + h(
"McRecMulGood")->GetXaxis()->GetTitle() +
"; Efficiency";
88 eff_good->SetTitle(eff_title.c_str() );
89 eff_good->SetMarkerColor(kRed);
90 eff_good->SetMarkerStyle(kFullCross);
92 eff_hist =
static_cast<TH1*
>( h(
"McRecMulGood")->Clone() );
93 eff_hist->Reset(
"ICES");
94 eff_hist->SetName( (std::string(eff_hist->GetName()) +
"_eff").c_str() );
95 eff_hist->SetYTitle(
"Efficiency");
96 eff_hist->GetListOfFunctions()->Add(eff_good);
vtxeval::VectorMcTracks getMcTracksMatchingMcVertex(const StMuDst &stMuDst, const StMuMcVertex *mcVertex)
StarEventHistContainer(const std::string name, TDirectory *motherDir=nullptr, const std::string option="")
std::vector< const StMuMcTrack * > VectorMcTracks
void FillEfficyHists(const StMuDst &event, const StMuMcVertex &mcVertex, const StMuPrimaryVertex *recoVertex=nullptr, const StMuPrimaryVertex *recoVertexMaxRank=nullptr)
void FillHists(const StMuDst &event)