Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DCustomAction_p3pi_hists.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DCustomAction_p3pi_hists.cc
4 // Created: Wed Jan 21 16:53:41 EST 2015
5 // Creator: jrsteven (on Linux ifarm1401 2.6.32-431.el6.x86_64 x86_64)
6 //
7 
9 
10 void DCustomAction_p3pi_hists::Initialize(JEventLoop* locEventLoop)
11 {
12  //CREATE THE HISTOGRAMS
13  //Since we are creating histograms, the contents of gDirectory will be modified: must use JANA-wide ROOT lock
14  japp->RootWriteLock(); //ACQUIRE ROOT LOCK!!
15  {
16  //Required: Create a folder in the ROOT output file that will contain all of the output ROOT objects (if any) for this action.
17  //If another thread has already created the folder, it just changes to it.
19 
20  dEgamma = GetOrCreate_Histogram<TH1I>("Egamma", "TAGGER photon energy; E_{#gamma}", 400, 0., 12.);
21 
22  dMM2_M3pi = GetOrCreate_Histogram<TH2I>("MM2_M3pi", "MM^{2} off #pi^{+}#pi^{-}#pi^{0} vs M_{#pi^{+}#pi^{-}#pi^{0}}; M_{#pi^{+}#pi^{-}#pi^{0}}; MM^{2}", 200, 0.0, 2.0, 200, -1., 1.);
23  dProton_dEdx_P = GetOrCreate_Histogram<TH2I>("Proton_dEdx_P","dE/dx vs p; p; dE/dx",200,0,2,500,0,5);
24  dProton_P_Theta = GetOrCreate_Histogram<TH2I>("Proton_P_Theta","p vs #theta; #theta; p (GeV/c)",180,0,180,120,0,12);
25 
26  dDeltaE_M3pi = GetOrCreate_Histogram<TH2I>("dDeltaE_M3pi", "#Delta E vs M_{#pi^{+}#pi^{-}#pi^{0}}; M_{#pi^{+}#pi^{-}#pi^{0}}; #Delta E (tagger - tracks)", 200, 0.0, 2.0, 200, -5., 5.);
27  dPhi3pi_PhiP = GetOrCreate_Histogram<TH2I>("Phi3pi_PhiP", "#phi_{#pi^{+}#pi^{-}#pi^{0}} vs. #phi_{p}; #phi_{p}; #phi_{#pi^{+}#pi^{-}#pi^{0}}", 360, -180.0, 180.0, 360, -180., 180.);
28  dDeltaPhi_M3pi = GetOrCreate_Histogram<TH2I>("DeltaPhi_M3pi", "#Delta#phi vs M_{#pi^{+}#pi^{-}#pi^{0}}; M_{#pi^{+}#pi^{-}#pi^{0}}; #Delta#phi", 200, 0.0, 2.0, 360, 0.0, 360.0);
29 
30  dEgamma_M3pi_ProtonTag = GetOrCreate_Histogram<TH2I>("dEgamma_M3pi_ProtonTag", "E_{#gamma} vs M_{#pi^{+}#pi^{-}#pi^{0}}; M_{#pi^{+}#pi^{-}#pi^{0}}; E_{#gamma}", 200, 0.0, 2.0, 240, 0., 6.);
31  dMppizero_M2pi = GetOrCreate_Histogram<TH2I>("Mppizero_M2pi", "Psuedo-Dalitz p#pi^{+}#pi^{-}#pi^{0}; M_{#pi^{+}#pi^{-}}; M_{p#pi^{0}}", 200, 0.0, 1.0, 200, 1.0, 4.0);
32  dMppiplus_M2pi = GetOrCreate_Histogram<TH2I>("Mppiplus_M2pi", "Psuedo-Dalitz p#pi^{+}#pi^{-}#pi^{0}; M_{#pi^{0}#pi^{-}}; M_{p#pi^{+}}", 200, 0.0, 1.0, 200, 1.0, 4.0);
33  dMppiminus_M2pi = GetOrCreate_Histogram<TH2I>("Mppiminus_M2pi", "Psuedo-Dalitz p#pi^{+}#pi^{-}#pi^{0}; M_{#pi^{0}#pi^{+}}; M_{p#pi^{-}}", 200, 0.0, 1.0, 200, 1.0, 4.0);
34 
35  dMM2_M3pi_CoplanarTag = GetOrCreate_Histogram<TH2I>("MM2_M3pi_CoplanarTag", "MM^{2} off p #pi^{+}#pi^{-}#pi^{0} vs M_{#pi^{+}#pi^{-}#pi^{0}}: Coplanar Tag; M_{#pi^{+}#pi^{-}#pi^{0}}; MM^{2}", 200, 0.0, 2.0, 200, -1., 1.);
36  dDeltaE_M3pi_CoplanarTag = GetOrCreate_Histogram<TH2I>("dDeltaE_M3pi_CoplanarTag", "#Delta E vs M_{#pi^{+}#pi^{-}#pi^{0}}: Coplanar Tag; M_{#pi^{+}#pi^{-}#pi^{0}}; #Delta E (tagger - p#pi^{+}#pi^{-}#pi^{0})", 200, 0.0, 2.0, 200, -5., 5.);
37  dMM2_DeltaE_CoplanarTag = GetOrCreate_Histogram<TH2I>("dMM2_DeltaE_CoplanarTag", "MM^{2} vs #Delta E: Coplanar Tag; #Delta E (tagger - p#pi^{+}#pi^{-}#pi^{0}); MM^{2}", 200, -5., 5., 200, -1., 1.);
38 
39  dMM2_M3pi_ProtonTag = GetOrCreate_Histogram<TH2I>("MM2_M3pi_ProtonTag", "MM^{2} off p #pi^{+}#pi^{-}#pi^{0} vs M_{#pi^{+}#pi^{-}#pi^{0}}: Proton Tag; M_{#pi^{+}#pi^{-}#pi^{0}}; MM^{2}", 200, 0.0, 2.0, 200, -1., 1.);
40  dDeltaE_M3pi_ProtonTag = GetOrCreate_Histogram<TH2I>("dDeltaE_M3pi_ProtonTag", "#Delta E vs M_{#pi^{+}#pi^{-}#pi^{0}}: Proton Tag; M_{#pi^{+}#pi^{-}#pi^{0}}; #Delta E (tagger - p#pi^{+}#pi^{-}#pi^{0})", 200, 0.0, 2.0, 200, -5., 5.);
41  dMM2_DeltaE_ProtonTag = GetOrCreate_Histogram<TH2I>("dMM2_DeltaE_ProtonTag", "MM^{2} vs #Delta E: Proton Tag; #Delta E (tagger - p#pi^{+}#pi^{-}#pi^{0}); MM^{2}", 200, -5., 5., 200, -1., 1.);
42 
43  //Return to the base directory
45  }
46  japp->RootUnLock(); //RELEASE ROOT LOCK!!
47 }
48 
49 bool DCustomAction_p3pi_hists::Perform_Action(JEventLoop* locEventLoop, const DParticleCombo* locParticleCombo)
50 {
51  const DParticleComboStep* locParticleComboStep = locParticleCombo->Get_ParticleComboStep(0);
52 
53  // get beam photon energy and final state particles
54  auto locBeamPhoton = Get_UseKinFitResultsFlag() ? locParticleComboStep->Get_InitialParticle() : locParticleComboStep->Get_InitialParticle_Measured();
55  auto locParticles = Get_UseKinFitResultsFlag() ? locParticleComboStep->Get_FinalParticles() : locParticleComboStep->Get_FinalParticles_Measured();
56  double locBeamPhotonEnergy = locBeamPhoton->energy();
57 
58  DLorentzVector locSumInitP4;
59  locSumInitP4.SetXYZM(0, 0, 0, 0.938);
60  locSumInitP4 += locBeamPhoton->lorentzMomentum();
61 
62  // calculate missing mass
64 
65  // calculate 3pi mass
68  DLorentzVector locPiPlusP4 = locParticleCombo->Get_ParticleComboStep(1)->Get_FinalParticle_Measured(0)->lorentzMomentum();
69  DLorentzVector locPiMinusP4 = locParticleCombo->Get_ParticleComboStep(1)->Get_FinalParticle_Measured(1)->lorentzMomentum();
70 
71  // reconstructed proton variables
72  double dEdx = 0.;
73  DLorentzVector locProtonP4;
74  const DChargedTrack* locChargedTrack = static_cast<const DChargedTrack*>(locParticleComboStep->Get_FinalParticle_SourceObject(1));
75  const DChargedTrackHypothesis* locChargedTrackHypothesis = locChargedTrack->Get_Hypothesis(Proton);
76  dEdx = locChargedTrackHypothesis->Get_TrackTimeBased()->dEdx()*1e6;
77  locProtonP4 = locChargedTrackHypothesis->lorentzMomentum();
78 
79  double dEdxCut = 2.2;
80 
81  //FILL HISTOGRAMS
82  //Since we are filling histograms local to this action, it will not interfere with other ROOT operations: can use action-wide ROOT lock
83  //Note, the mutex is unique to this DReaction + action_string combo: actions of same class with different hists will have a different mutex
84  Lock_Action(); //ACQUIRE ROOT LOCK!!
85  {
86  // Fill histograms here
87  dEgamma->Fill(locBeamPhotonEnergy);
88 
89  dMM2_M3pi->Fill(locOmegaP4.M(), locMissingP4.M2());
90  dDeltaE_M3pi->Fill(locOmegaP4.M(),locMissingP4.E());
91 
92  double locDeltaPhi = (locProtonP4.Phi() - locOmegaP4.Phi())*180./TMath::Pi();
93  if(locDeltaPhi > 360.) locDeltaPhi -= 360.;
94  if(locDeltaPhi < 0.) locDeltaPhi += 360.;
95  dPhi3pi_PhiP->Fill(locProtonP4.Phi()*180./TMath::Pi(), locOmegaP4.Phi()*180./TMath::Pi());
96  dDeltaPhi_M3pi->Fill(locOmegaP4.M(), locDeltaPhi);
97 
98  // require proton and omega are back-to-back
99  if(locDeltaPhi < 175. || locDeltaPhi > 185.)
100  {
101  Unlock_Action(); //RELEASE ROOT LOCK!!
102  return true;
103  }
104  dMM2_M3pi_CoplanarTag->Fill(locOmegaP4.M(), locMissingP4.M2());
105  dDeltaE_M3pi_CoplanarTag->Fill(locOmegaP4.M(),locMissingP4.E());
106  if(locOmegaP4.M() > 0.7 && locOmegaP4.M() < 0.9)
107  dMM2_DeltaE_CoplanarTag->Fill(locMissingP4.E(), locMissingP4.M2());
108 
109  // tag proton with dE/dx
110  if(dEdx > dEdxCut) {
111  dMM2_M3pi_ProtonTag->Fill(locOmegaP4.M(), locMissingP4.M2());
112  dDeltaE_M3pi_ProtonTag->Fill(locOmegaP4.M(),locMissingP4.E());
113  if(locOmegaP4.M() > 0.7 && locOmegaP4.M() < 0.9)
114  dMM2_DeltaE_ProtonTag->Fill(locMissingP4.E(), locMissingP4.M2());
115 
116  // 3pi invariant mass for exclusive
117  if(fabs(locMissingP4.M2()) < 0.05 && fabs(locMissingP4.E()) < 0.5) {
118  dEgamma_M3pi_ProtonTag->Fill(locOmegaP4.M(),locBeamPhotonEnergy);
119 
120  if(locBeamPhotonEnergy > 2.5 && locBeamPhotonEnergy < 3.0){
121  dMppizero_M2pi->Fill((locPiPlusP4+locPiMinusP4).M(), (locProtonP4+locPi0P4).M());
122  dMppiplus_M2pi->Fill((locPi0P4+locPiMinusP4).M(), (locProtonP4+locPiPlusP4).M());
123  dMppiminus_M2pi->Fill((locPiPlusP4+locPi0P4).M(), (locProtonP4+locPiMinusP4).M());
124  }
125  }
126  }
127 
128  if(fabs(locMissingP4.M2()) < 0.05 && fabs(locMissingP4.E()) < 0.5) {
129 
130  dProton_dEdx_P->Fill(locProtonP4.Vect().Mag(), dEdx);
131  dProton_P_Theta->Fill(locProtonP4.Vect().Theta()*180/TMath::Pi(), locProtonP4.Vect().Mag());
132 
133  }
134  }
135  Unlock_Action(); //RELEASE ROOT LOCK!!
136 
137  return true; //return false if you want to use this action to apply a cut (and it fails the cut!)
138 }
vector< const DKinematicData * > Get_FinalParticles_Measured(void) const
TDirectoryFile * ChangeTo_BaseDirectory(void)
const DKinematicData * Get_FinalParticle_Measured(size_t locFinalParticleIndex) const
const DReaction * Get_Reaction(void) const
void Initialize(JEventLoop *locEventLoop)
const DTrackTimeBased * Get_TrackTimeBased(void) const
TDirectoryFile * CreateAndChangeTo_ActionDirectory(void)
bool Perform_Action(JEventLoop *locEventLoop, const DParticleCombo *locParticleCombo)
const DKinematicData * Get_InitialParticle_Measured(void) const
TLorentzVector DLorentzVector
const DChargedTrackHypothesis * Get_Hypothesis(Particle_t locPID) const
Definition: DChargedTrack.h:59
bool Get_UseKinFitResultsFlag(void) const
JApplication * japp
const DAnalysisUtilities * dAnalysisUtilities
DLorentzVector Calc_FinalStateP4(const DReaction *locReaction, const DParticleCombo *locParticleCombo, size_t locStepIndex, bool locUseKinFitDataFlag) const
DLorentzVector Calc_MissingP4(const DReaction *locReaction, const DParticleCombo *locParticleCombo, bool locUseKinFitDataFlag) const
vector< const DKinematicData * > Get_FinalParticles(void) const
DLorentzVector lorentzMomentum(void) const
const JObject * Get_FinalParticle_SourceObject(size_t locFinalParticleIndex) const
const DKinematicData * Get_InitialParticle(void) const
double dEdx(void) const
const DParticleComboStep * Get_ParticleComboStep(size_t locStepIndex) const