Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DCustomAction_p2gamma_hists.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DCustomAction_p2gamma_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_p2gamma_hists::Initialize(JEventLoop* locEventLoop)
11 {
12  DApplication* dapp=dynamic_cast<DApplication*>(locEventLoop->GetJApplication());
13  JCalibration *jcalib = dapp->GetJCalibration((locEventLoop->GetJEvent()).GetRunNumber());
14 
15  // Parameters for event selection to fill histograms
16  endpoint_energy = 12.;
17  map<string, double> photon_endpoint_energy;
18  if(jcalib->Get("/PHOTON_BEAM/endpoint_energy", photon_endpoint_energy) == false) {
19  endpoint_energy = photon_endpoint_energy["PHOTON_BEAM_ENDPOINT_ENERGY"];
20  }
22 
23  cohmin_energy = 0.;
24  cohedge_energy = 12.;
25  map<string, double> photon_beam_param;
26  if(jcalib->Get("test/PHOTON_BEAM/coherent_energy", photon_beam_param) == false) {
27  cohmin_energy = photon_beam_param["cohmin_energy"];
28  cohedge_energy = photon_beam_param["cohedge_energy"];
29  }
30 
31  dEdxCut = 2.2;
32  minMM2Cut = -0.05;
33  maxMM2Cut = 0.05;
34  missingEnergyCut = 1.0;
35  min2gMassCut = 0.10;
36  max2gMassCut = 0.16;
37 
38  //CREATE THE HISTOGRAMS
39  //Since we are creating histograms, the contents of gDirectory will be modified: must use JANA-wide ROOT lock
40  japp->RootWriteLock(); //ACQUIRE ROOT LOCK!!
41  {
42  //Required: Create a folder in the ROOT output file that will contain all of the output ROOT objects (if any) for this action.
43  //If another thread has already created the folder, it just changes to it.
45 
46  dEgamma = GetOrCreate_Histogram<TH1I>("Egamma", "TAGGER photon energy; E_{#gamma}", endpoint_energy_bins, 0., endpoint_energy);
47 
48  dMM2_M2g = GetOrCreate_Histogram<TH2I>("MM2_M2g", "MM^{2} off p #gamma#gamma vs M_{#gamma#gamma}; M_{#gamma#gamma}; MM^{2}", 500, 0.0, 1.0, 200, -1., 1.);
49  dProton_dEdx_P = GetOrCreate_Histogram<TH2I>("Proton_dEdx_P","dE/dx vs p; p; dE/dx",200,0,2,500,0,5);
50  dProton_P_Theta = GetOrCreate_Histogram<TH2I>("Proton_P_Theta","p vs #theta; #theta; p (GeV/c)",180,0,180,120,0,12);
51  dProtonPhi_Egamma = GetOrCreate_Histogram<TH2I>("ProtonPhi_Egamma","#phi vs E_{#gamma}; E_{#gamma}; proton #phi",endpoint_energy_bins, 0., endpoint_energy,360,-180,180);
52  dProtonPhi_Theta = GetOrCreate_Histogram<TH2I>("ProtonPhi_Theta","#phi vs #theta_{CM}; #pi^{0} #theta_{CM}; proton #phi",180,0,180,360,-180,180);
53  dProtonPhi_t = GetOrCreate_Histogram<TH2I>("ProtonPhi_t","#psi vs |t|; |t|; proton #phi",1000,0.,5.,360,-180,180);
54 
55  dPi0Phi_Egamma = GetOrCreate_Histogram<TH2I>("Pi0Phi_Egamma","#phi vs E_{#gamma}; E_{#gamma}; #pi^{0} #phi",endpoint_energy_bins, 0., endpoint_energy,360,-180,180);
56  dPi0Phi_Theta = GetOrCreate_Histogram<TH2I>("Pi0Phi_Theta","#phi vs #theta_{CM}; #pi^{0} #theta_{CM}; #pi^{0} #phi",180,0,180,360,-180,180);
57  dPi0EgammaCorr = GetOrCreate_Histogram<TH2I>("Pi0EgammaCorr","E_{#gamma 1} vs E_{#gamma 2}; E_{#gamma 2}; E_{#gamma 1}",100,0.,5., 100.,0.,5.);
58 
59  dDeltaE_M2g = GetOrCreate_Histogram<TH2I>("dDeltaE_M2g", "#Delta E vs M_{#gamma#gamma}; M_{#gamma#gamma}; #Delta E (tagger - p#gamma#gamma)", 500, 0.0, 1.0, 200, -5., 5.);
60  dPhi2g_PhiP = GetOrCreate_Histogram<TH2I>("Phi2g_PhiP", "#phi_{#gamma#gamma} vs. #phi_{p}; #phi_{p}; #phi_{#gamma#gamma}", 360, -180.0, 180.0, 360, -180., 180.);
61  dDeltaPhi_M2g = GetOrCreate_Histogram<TH2I>("DeltaPhi_M2g", "#Delta#phi vs M_{#gamma#gamma}; M_{#gamma#gamma}; #Delta#phi", 500, 0.0, 1.0, 360, 0.0, 360.0);
62 
63  dEgamma_M2g_ProtonTag = GetOrCreate_Histogram<TH2I>("dEgamma_M2g_ProtonTag", "E_{#gamma} vs M_{#gamma#gamma}; M_{#gamma#gamma}; E_{#gamma}", 500, 0.0, 1.0, 240, 0., 6.);
64 
65  dMM2_M2g_CoplanarTag = GetOrCreate_Histogram<TH2I>("MM2_M2g_CoplanarTag", "MM^{2} off p #gamma#gamma vs M_{#gamma#gamma}: Coplanar Tag; M_{#gamma#gamma}; MM^{2}", 500, 0.0, 1.0, 200, -1., 1.);
66  dDeltaE_M2g_CoplanarTag = GetOrCreate_Histogram<TH2I>("dDeltaE_M2g_CoplanarTag", "#Delta E vs M_{#gamma#gamma}: Coplanar Tag; M_{#gamma#gamma}; #Delta E (tagger - p#gamma#gamma)", 500, 0.0, 1.0, 200, -5., 5.);
67  dMM2_DeltaE_CoplanarTag = GetOrCreate_Histogram<TH2I>("dMM2_DeltaE_CoplanarTag", "MM^{2} vs #Delta E: Coplanar Tag; #Delta E (tagger - p#gamma#gamma); MM^{2}", 200, -5., 5., 200, -1., 1.);
68 
69  dMM2_M2g_ProtonTag = GetOrCreate_Histogram<TH2I>("MM2_M2g_ProtonTag", "MM^{2} off p #gamma#gamma vs M_{#gamma#gamma}: Proton Tag; M_{#gamma#gamma}; MM^{2}", 500, 0.0, 1.0, 200, -1., 1.);
70  dDeltaE_M2g_ProtonTag = GetOrCreate_Histogram<TH2I>("dDeltaE_M2g_ProtonTag", "#Delta E vs M_{#gamma#gamma}: Proton Tag; M_{#gamma#gamma}; #Delta E (tagger - p#gamma#gamma)", 500, 0.0, 1.0, 200, -5., 5.);
71  dMM2_DeltaE_ProtonTag = GetOrCreate_Histogram<TH2I>("dMM2_DeltaE_ProtonTag", "MM^{2} vs #Delta E: Proton Tag; #Delta E (tagger - p#gamma#gamma); MM^{2}", 200, -5., 5., 200, -1., 1.);
72 
73  dMinEgamma_M2g = GetOrCreate_Histogram<TH2I>("MinEgamma_M2g", "Minimum E_{#gamma} vs M_{#gamma#gamma}; M_{#gamma#gamma}; Minimum E_{#gamma}", 500, 0., 1., 200, 0., 2.);
74 
75  //Return to the base directory
77  }
78  japp->RootUnLock(); //RELEASE ROOT LOCK!!
79 }
80 
81 bool DCustomAction_p2gamma_hists::Perform_Action(JEventLoop* locEventLoop, const DParticleCombo* locParticleCombo)
82 {
83  const DParticleComboStep* locParticleComboStep = locParticleCombo->Get_ParticleComboStep(0);
84 
85  // get beam photon energy and final state particles
86  auto locBeamPhoton = Get_UseKinFitResultsFlag() ? locParticleComboStep->Get_InitialParticle() : locParticleComboStep->Get_InitialParticle_Measured();
87  auto locParticles = Get_UseKinFitResultsFlag() ? locParticleComboStep->Get_FinalParticles() : locParticleComboStep->Get_FinalParticles_Measured();
88  double locBeamPhotonEnergy = locBeamPhoton->energy();
89 
90  // calculate missing mass
91  DLorentzVector locMissingP4;
92  DLorentzVector locProtonP4Init(0,0,0,0.938);
93  locMissingP4 += locProtonP4Init;
94  locMissingP4 += locBeamPhoton->lorentzMomentum();
95  DLorentzVector locSumInitP4 = locMissingP4;
96 
97  // reconstructed proton variables
98  const DChargedTrack* locChargedTrack = static_cast<const DChargedTrack*>(locParticleComboStep->Get_FinalParticle_SourceObject(2));
99  const DChargedTrackHypothesis* locChargedTrackHypothesis = locChargedTrack->Get_Hypothesis(Proton);
100  auto locTrackTimeBased = locChargedTrackHypothesis->Get_TrackTimeBased();
101  double dEdx = locTrackTimeBased->dEdx()*1e6; // convert to keV
102  DLorentzVector locProtonP4 = locChargedTrackHypothesis->lorentzMomentum();
103 
104  // calculate missing mass
105  DLorentzVector loc2g_P4;
106  set<pair<const JObject*, Particle_t> > locSourceObjects;
107  for(size_t loc_i = 0; loc_i < 3; ++loc_i) {
108  locMissingP4 -= locParticles[loc_i]->lorentzMomentum();
109 
110  if(locParticles[loc_i]->PID() == Gamma) {
111  locSourceObjects.insert(pair<const JObject*, Particle_t>(locParticleComboStep->Get_FinalParticle_SourceObject(loc_i), locParticles[loc_i]->PID()));
112  loc2g_P4 += locParticles[loc_i]->lorentzMomentum();
113  }
114  }
115 
116  bool fillHist = true;
117  if(dPreviousSourceObjects.find(locSourceObjects) != dPreviousSourceObjects.end())
118  fillHist = false; //dupe: already histed!
119  dPreviousSourceObjects.insert(locSourceObjects);
120 
121  DLorentzVector locGamma1P4 = locParticleCombo->Get_ParticleComboStep(0)->Get_FinalParticle_Measured(0)->lorentzMomentum();
122  DLorentzVector locGamma2P4 = locParticleCombo->Get_ParticleComboStep(0)->Get_FinalParticle_Measured(1)->lorentzMomentum();
123  double locMinEgamma = min(locGamma1P4.E(), locGamma2P4.E());
124 
125  // production kinematics
126  TVector3 locBoostVector = -1.*locSumInitP4.BoostVector();
127  TLorentzVector locGamma_P4CMFrame = locBeamPhoton->lorentzMomentum();
128  TLorentzVector locPi0_P4CMFrame = loc2g_P4;
129  locGamma_P4CMFrame.Boost(locBoostVector);
130  locPi0_P4CMFrame.Boost(locBoostVector);
131  double locThetaCM = locGamma_P4CMFrame.Vect().Angle(locPi0_P4CMFrame.Vect());
132  TLorentzVector locDelta = (locProtonP4 - locProtonP4Init);
133  double t = locDelta.M2();
134 
135  //FILL HISTOGRAMS
136  //Since we are filling histograms local to this action, it will not interfere with other ROOT operations: can use action-wide ROOT lock
137  //Note, the mutex is unique to this DReaction + action_string combo: actions of same class with different hists will have a different mutex
138  Lock_Action(); //ACQUIRE ROOT LOCK!!
139  {
140  // Fill histograms here
141  dEgamma->Fill(locBeamPhotonEnergy);
142 
143  if(fillHist)
144  dMinEgamma_M2g->Fill(loc2g_P4.M(), locMinEgamma);
145 
146  dMM2_M2g->Fill(loc2g_P4.M(), locMissingP4.M2());
147  dDeltaE_M2g->Fill(loc2g_P4.M(),locMissingP4.E());
148 
149  double locDeltaPhi = (locProtonP4.Phi() - loc2g_P4.Phi())*180./TMath::Pi();
150  if(locDeltaPhi > 360.) locDeltaPhi -= 360.;
151  if(locDeltaPhi < 0.) locDeltaPhi += 360.;
152  dPhi2g_PhiP->Fill(locProtonP4.Phi()*180./TMath::Pi(), loc2g_P4.Phi()*180./TMath::Pi());
153  dDeltaPhi_M2g->Fill(loc2g_P4.M(), locDeltaPhi);
154 
155  // require proton and pi0 are back-to-back
156  if(locDeltaPhi > 175. && locDeltaPhi < 185.) {
157 
158  dMM2_M2g_CoplanarTag->Fill(loc2g_P4.M(), locMissingP4.M2());
159  dDeltaE_M2g_CoplanarTag->Fill(loc2g_P4.M(), locMissingP4.E());
160  if(loc2g_P4.M() > min2gMassCut && loc2g_P4.M() < max2gMassCut)
161  dMM2_DeltaE_CoplanarTag->Fill(locMissingP4.E(), locMissingP4.M2());
162 
163  // tag proton with dE/dx
164  if(dEdx > dEdxCut) {
165  dMM2_M2g_ProtonTag->Fill(loc2g_P4.M(), locMissingP4.M2());
166  dDeltaE_M2g_ProtonTag->Fill(loc2g_P4.M(), locMissingP4.E());
167  if(loc2g_P4.M() > min2gMassCut && loc2g_P4.M() < max2gMassCut)
168  dMM2_DeltaE_ProtonTag->Fill(locMissingP4.E(), locMissingP4.M2());
169 
170  // di-photon invariant mass for exclusive g+p -> p + 2g
171  if(locMissingP4.M2() > minMM2Cut && locMissingP4.M2() < maxMM2Cut && fabs(locMissingP4.E()) < missingEnergyCut) {
172  dEgamma_M2g_ProtonTag->Fill(loc2g_P4.M(),locBeamPhotonEnergy);
173  }
174  }
175 
176  // for pi0 candidates require recoil proton
177  if(loc2g_P4.M() > min2gMassCut && loc2g_P4.M() < max2gMassCut && locMissingP4.M2() > minMM2Cut && locMissingP4.M2() < maxMM2Cut && fabs(locMissingP4.E()) < missingEnergyCut) {
178  dProton_dEdx_P->Fill(locProtonP4.Vect().Mag(), dEdx);
179  dProton_P_Theta->Fill(locProtonP4.Vect().Theta()*180/TMath::Pi(), locProtonP4.Vect().Mag());
180  dProtonPhi_Egamma->Fill(locBeamPhotonEnergy, locProtonP4.Phi()*180/TMath::Pi());
181  dPi0Phi_Egamma->Fill(locBeamPhotonEnergy, loc2g_P4.Phi()*180/TMath::Pi());
182 
183  if(locBeamPhotonEnergy > cohmin_energy && locBeamPhotonEnergy < cohedge_energy){
184  dProtonPhi_Theta->Fill(locThetaCM*180/TMath::Pi(), locProtonP4.Phi()*180/TMath::Pi());
185  dProtonPhi_t->Fill(fabs(t), locProtonP4.Phi()*180/TMath::Pi());
186  dPi0Phi_Theta->Fill(locThetaCM*180/TMath::Pi(), loc2g_P4.Phi()*180/TMath::Pi());
187  dPi0EgammaCorr->Fill(locParticles[0]->energy(),locParticles[1]->energy());
188  }
189  }
190  }
191  }
192  Unlock_Action(); //RELEASE ROOT LOCK!!
193 
194  return true; //return false if you want to use this action to apply a cut (and it fails the cut!)
195 }
vector< const DKinematicData * > Get_FinalParticles_Measured(void) const
TDirectoryFile * ChangeTo_BaseDirectory(void)
DApplication * dapp
const DKinematicData * Get_FinalParticle_Measured(size_t locFinalParticleIndex) const
const DTrackTimeBased * Get_TrackTimeBased(void) const
TDirectoryFile * CreateAndChangeTo_ActionDirectory(void)
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
vector< const DKinematicData * > Get_FinalParticles(void) const
bool Perform_Action(JEventLoop *locEventLoop, const DParticleCombo *locParticleCombo)
DLorentzVector lorentzMomentum(void) const
const JObject * Get_FinalParticle_SourceObject(size_t locFinalParticleIndex) const
const DKinematicData * Get_InitialParticle(void) const
double dEdx(void) const
set< set< pair< const JObject *, Particle_t > > > dPreviousSourceObjects
const DParticleComboStep * Get_ParticleComboStep(size_t locStepIndex) const
void Initialize(JEventLoop *locEventLoop)