Skip to content

Commit

Permalink
Merge pull request #137 from JeffersonLab/aaust_AnalyzeCutActions
Browse files Browse the repository at this point in the history
Fixed uniqueness tracking for invariant mass histograms in AnalyzeCutActions.
  • Loading branch information
niwgit authored Nov 23, 2020
2 parents e22984f + b7fc276 commit c4ad7e0
Show file tree
Hide file tree
Showing 2 changed files with 8 additions and 10 deletions.
12 changes: 5 additions & 7 deletions libraries/DSelector/DHistogramActions.cc
Original file line number Diff line number Diff line change
Expand Up @@ -58,7 +58,6 @@ bool DHistogramAction_AnalyzeCutActions::Perform_ActionWeight(double weight = 1.

for (auto const &cut_iter : dHistsWithoutCuts)
{
dPreviouslyHistogrammed.clear();
bool locFill = true;
for (auto const &action_iter : dAllAnalysisActions)
{
Expand All @@ -78,28 +77,27 @@ bool DHistogramAction_AnalyzeCutActions::Perform_ActionWeight(double weight = 1.
}

if (locFill)
Fill_Hists(cut_iter.second, locIndexCombos, weight);
Fill_Hists(cut_iter.second, locIndexCombos, weight, cut_iter.first);
}

if (!locComboCut)
Fill_Hists(dHist_InvariantMass_allcuts, locIndexCombos, weight);
Fill_Hists(dHist_InvariantMass_allcuts, locIndexCombos, weight, "allcuts");

return true;
}

bool DHistogramAction_AnalyzeCutActions::Fill_Hists(TH1I* locHist, set<set<size_t> > locIndexCombos, double weight)
bool DHistogramAction_AnalyzeCutActions::Fill_Hists(TH1I* locHist, set<set<size_t> > locIndexCombos, double weight, string locActionName)
{
dPreviouslyHistogrammed.clear();
double locMass = 0.0;
set<set<size_t> >::iterator locComboIterator = locIndexCombos.begin();
for(; locComboIterator != locIndexCombos.end(); ++locComboIterator)
{
map<unsigned int, set<Int_t> > locSourceObjects;
TLorentzVector locFinalStateP4 = dAnalysisUtilities.Calc_FinalStateP4(dParticleComboWrapper, dStepIndex, *locComboIterator, locSourceObjects, dUseKinFitFlag);
locMass = locFinalStateP4.M();
if(dPreviouslyHistogrammed.find(locSourceObjects) == dPreviouslyHistogrammed.end())
if(dPreviouslyHistogrammed[locActionName].find(locSourceObjects) == dPreviouslyHistogrammed[locActionName].end())
{
dPreviouslyHistogrammed.insert(locSourceObjects);
dPreviouslyHistogrammed[locActionName].insert(locSourceObjects);
locHist->Fill(locMass, weight);
}
}
Expand Down
6 changes: 3 additions & 3 deletions libraries/DSelector/DHistogramActions.h
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@ class DHistogramAction_AnalyzeCutActions : public DAnalysisAction
bool Perform_ActionWeight(double weight);

private:
bool Fill_Hists(TH1I* locHist, set<set<size_t>> locIndexCombos, double weight);
bool Fill_Hists(TH1I* locHist, set<set<size_t>> locIndexCombos, double weight, string locActionName);
vector<DAnalysisAction*> dAllAnalysisActions;
Particle_t dInitialPID;
int dStepIndex;
Expand All @@ -57,8 +57,8 @@ class DHistogramAction_AnalyzeCutActions : public DAnalysisAction
// string comes from Get_ActionName() and TH1I* is the histogram without that cut
map<string, TH1I*> dHistsWithoutCuts;

// uniqueness tracking
set<map<unsigned int, set<Int_t> > > dPreviouslyHistogrammed;
// uniqueness tracking for each CutAction separately
map<string, set<map<unsigned int, set<Int_t> > > > dPreviouslyHistogrammed;
};

class DHistogramAction_ParticleComboKinematics : public DAnalysisAction
Expand Down

0 comments on commit c4ad7e0

Please sign in to comment.