HighLAND
FGDECalSMRDMatchEffSystematics.cxx
1 #include "FGDECalSMRDMatchEffSystematics.hxx"
2 #include "CutUtils.hxx"
3 #include "SystematicUtils.hxx"
4 #include "SystId.hxx"
5 
6 //#define DEBUG
7 
8 
9 //********************************************************************
10 FGDECalSMRDMatchEffSystematics::FGDECalSMRDMatchEffSystematics(bool computecounters):
11  FGDECalMatchEffSystematics(computecounters, "FGDECalSMRDMatchEff", k2D_EFF_ASSYMMETRIC){
12  //********************************************************************
13 
14 }
15 
16 //********************************************************************
18  //********************************************************************
19 
20 
21  if(_computecounters)
23 
24  (void)event;
25 
26  // Get the SystBox for this event, and the appropriate selection and branch
27  SystBoxB* SystBox = GetSystBox(event, box.SelectionEnabledIndex, box.SuccessfulBranch);
28 
29  Weight_h eventWeight = 1.;
30 
31 #ifdef DEBUG
32  std::cout << "FGDECalSMRDMatchEffSystematics::Apply(): " << SystBox->nRelevantTrueObjects<< std::endl;
33 #endif
34 
35  // Loop over relevant true tracks
36  for (Int_t itrk=0;itrk<SystBox->nRelevantTrueObjects;itrk++){
37 
38  AnaTrueParticleB* truePart = static_cast<AnaTrueParticleB*>(SystBox->RelevantTrueObjects[itrk]);
39 
40  AnaTrackB* recoTrack = static_cast<AnaTrackB*>(SystBox->RelevantTrueObjectsReco[itrk]);
41 
42  if (!truePart || !recoTrack) continue;
43 
44  // Do fine-tuning of the track relevance via the selection
45  if (!sel.IsRelevantTrueObjectForSystematicInToy(event, box, truePart, SystId::kFgdECalSmrdMatchEff, box.SuccessfulBranch)) continue;
46 
47  SubDetId::SubDetEnum smrd_det[5];
48 
49  int nsmrd_det = anaUtils::GetSMRDDetCrossed(truePart, smrd_det);
50 
51  if(nsmrd_det==0) continue;
52 
53  // Get ecal segment (wont work if call for TECAL)
54  AnaParticleB* ecal = anaUtils::GetSegmentWithMostNodesInDet(*(recoTrack), SubDetId::kECAL);
55 
56  if(!ecal)
57  continue;
58 
59  Float_t* dir = anaUtils::GetSLineDir(ecal->PositionStart, ecal->PositionEnd);
60 
61  // Check the direction of the reco track
62  dir[2] *= ((recoTrack->DirectionStart[2]*dir[2] >= 0) ? 1. : -1.);
63 
64  BinnedParamsParams params;
65  int index;
66 
67  if (!GetBinValues(smrd_det[0], dir[2], params, index)) continue;
68 
69  // Found the correspondence, now check detector bits
70  bool found = SubDetId::GetDetectorUsed(recoTrack->Detector, smrd_det[0]);
71 
72  eventWeight *= systUtils::ComputeEffLikeWeight(found, toy, GetIndex(), index, params);
73 
74 #ifdef DEBUG
75  std::cout<<"fgd-ecal-smrd found "<< found<<" eventWeight "<<eventWeight<<std::endl;
76 #endif
77 
78  if(_computecounters)
79  UpdateEfficiencyCounter(index,found);
80 
81  }
82 
83  return eventWeight;
84 }
85 
86 //********************************************************************
88  //********************************************************************
89 
90  (void)event;
91 
92  //should cross SMRD so to be potentially reconstructable
93  if (anaUtils::TrueParticleCrossesSMRD(static_cast<const AnaTrueParticleB*>(&track)))
94  return true;
95 
96  return false;
97 
98 }
99 
100 //**************************************************
102  //**************************************************
103 
104  (void)event;
105  //should use Barrel ECal since ECal iso-tracking efficiency is treated elsewhere
106 
107  if (SubDetId::GetDetectorUsed( track.Detector, SubDetId::kFGD) &&
108  SubDetId::GetDetectorUsed( track.Detector, SubDetId::kTECAL) &&
109  !SubDetId::GetDetectorUsed( track.Detector, SubDetId::kDSECAL) &&
110  !SubDetId::GetDetectorUsed( track.Detector, SubDetId::kPECAL) &&
111  !SubDetId::GetDetectorUsed( track.Detector, SubDetId::kTPC))
112  return true;
113 
114  return false;
115 }
116 
117 
Float_t PositionStart[4]
The reconstructed start position of the particle.
unsigned long Detector
Int_t GetIndex() const
Return the index of this systematic.
Int_t SelectionEnabledIndex
The enabled index of this selection this ToyBox belongs to.
Definition: ToyBoxB.hxx:49
SystBoxB * GetSystBox(const AnaEventC &event, Int_t isel=0, Int_t ibranch=0) const
Get the SystBox corresponding to a selection, branch and event.
static bool GetDetectorUsed(unsigned long BitField, SubDetId::SubDetEnum det)
Method to see if a certain subdetector or subdetector system is used.
Definition: SubDetId.cxx:40
bool GetBinValues(Float_t value, Float_t &mean, Float_t &sigma)
Gets the bin values for a 1D source.
bool UpdateEfficiencyCounter(Int_t index, bool correct)
Update the efficiency variables _ncorrect and _nwrong.
Representation of a true Monte Carlo trajectory/particle.
int GetSMRDDetCrossed(const AnaTrueParticleB *track, SubDetId::SubDetEnum det[])
SMRD detectors crossed.
Definition: TruthUtils.cxx:367
bool TrueParticleCrossesSMRD(const AnaTrueParticleB *track, Float_t min_sep=48.)
Definition: TruthUtils.cxx:186
Int_t SuccessfulBranch
The branch that is successful for this toy in the selection this ToyBox belongs to.
Definition: ToyBoxB.hxx:46
SubDetEnum
Enumeration of all detector systems and subdetectors.
Definition: SubDetId.hxx:25
Representation of a global track.
Int_t nRelevantTrueObjects
Array of Relevant True RecObjects for each systematic.
Definition: SystBoxB.hxx:24
bool IsRelevantTrueObject(const AnaEventC &event, const AnaTrueObjectC &track) const
Is this track relevant for this systematic ?
bool IsRelevantRecObject(const AnaEventC &event, const AnaRecObjectC &track) const
Is this track relevant for this systematic ?
Float_t DirectionStart[3]
The reconstructed start direction of the particle.
AnaRecObjectC ** RelevantTrueObjectsReco
Definition: SystBoxB.hxx:29
Weight_h ComputeWeight(const ToyExperiment &, const AnaEventC &, const ToyBoxB &)
Apply the systematic.
Float_t * GetSLineDir(Float_t *start, Float_t *end)
Direction assuming straight line between given start and end points.
void InitializeEfficiencyCounter()
Initialize counters.
virtual bool IsRelevantTrueObjectForSystematicInToy(const AnaEventC &, const ToyBoxB &, AnaTrueObjectC *, SystId_h syst_index, Int_t branch=0) const
Is this true track relevant for a given systematic (after selection, called for each toy) ...
AnaParticleB * GetSegmentWithMostNodesInDet(const AnaTrackB &track, SubDetId::SubDetEnum det)
Method to get the subtrack with most nodes in a given detector.
Representation of a reconstructed particle (track or shower).
Float_t PositionEnd[4]
The reconstructed end position of the particle.