Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
p2pi_trees/DCustomAction_p2pi_hists.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DCustomAction_p2pi_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_p2pi_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  }
21  else {
22  jout<<"No /PHOTON_BEAM/endpoint_energy for this run number: using default of 12 GeV"<<endl;
23  }
25 
26  cohmin_energy = 0.;
27  cohedge_energy = 12.;
28  map<string, double> photon_beam_param;
29  if(jcalib->Get("/ANALYSIS/beam_asymmetry/coherent_energy", photon_beam_param) == false) {
30  cohmin_energy = photon_beam_param["cohmin_energy"];
31  cohedge_energy = photon_beam_param["cohedge_energy"];
32  }
33  else {
34  jout<<"No /ANALYSIS/beam_asymmetry/coherent_energy for this run number: using default range of 0-12 GeV"<<endl;
35  }
36 
37  dEdxCut = 2.2;
38  minMMCut = 0.8;
39  maxMMCut = 1.05;
40  minMM2Cut = -0.006;
41  maxMM2Cut = 0.004;
42  missingEnergyCut = 1.0;
43  minRhoMassCut = 0.6;
44  maxRhoMassCut = 0.88;
45 
46  // get PID algos
47  const DParticleID* locParticleID = NULL;
48  locEventLoop->GetSingle(locParticleID);
49  dParticleID = locParticleID;
50 
51  locEventLoop->GetSingle(dAnalysisUtilities);
52 
53  // check if a particle is missing
54  Get_Reaction()->Get_MissingPID(dMissingPID);
55 
56  //CREATE THE HISTOGRAMS
57  //Since we are creating histograms, the contents of gDirectory will be modified: must use JANA-wide ROOT lock
58  japp->RootWriteLock(); //ACQUIRE ROOT LOCK!!
59  {
60  //Required: Create a folder in the ROOT output file that will contain all of the output ROOT objects (if any) for this action.
61  //If another thread has already created the folder, it just changes to it.
63 
64  dEgamma = GetOrCreate_Histogram<TH1I>("Egamma", "TAGGER photon energy; E_{#gamma}", endpoint_energy_bins, 0., endpoint_energy);
65 
66  if(dMissingPID == Proton) {
67  dMM_M2pi = GetOrCreate_Histogram<TH2I>("MM_M2pi", "MM off #pi^{+}#pi^{-} vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; MM", 200, 0.0, 2.0, 200, 0., 4.);
68  dMM_M2pi_noEle = GetOrCreate_Histogram<TH2I>("MM_M2pi_noEle", "MM off #pi^{+}#pi^{-} vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; MM", 200, 0.0, 2.0, 200, 0., 4.);
69  dt_M2pi_noEle = GetOrCreate_Histogram<TH2I>("t_M2pi_noEle", "|t| vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; t", 200, 0.0, 2.0, 1000, 0., 5.);
70  dPiPlusPsi_t = GetOrCreate_Histogram<TH2I>("PiPlusPsi_t","#psi vs |t|; |t|; #pi^{+} #psi",1000,0.,5.,360,-180,180);
71  }
72  else {
73  dMM2_M2pi = GetOrCreate_Histogram<TH2I>("MM2_M2pi", "MM^{2} off #pi^{+}#pi^{-} vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; MM^{2}", 200, 0.0, 2.0, 200, -0.1, 0.1);
74  dProton_dEdx_P = GetOrCreate_Histogram<TH2I>("Proton_dEdx_P","dE/dx vs p; p; dE/dx",200,0,2,500,0,25);
75  dProton_P_Theta = GetOrCreate_Histogram<TH2I>("Proton_P_Theta","p vs #theta; #theta; p (GeV/c)",180,0,180,120,0,12);
76  dDeltaE_M2pi = GetOrCreate_Histogram<TH2I>("DeltaE_M2pi", "#DeltaE vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; #DeltaE (tagger - tracks)", 200, 0.0, 2.0, 200, -5., 5.);
77  dDeltaE_M2pi_ProtonTag = GetOrCreate_Histogram<TH2I>("DeltaE_M2pi_ProtonTag", "#DeltaE vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; #DeltaE (tagger - tracks)", 200, 0.0, 2.0, 200, -5., 5.);
78 
79  dDalitz_p2pi = GetOrCreate_Histogram<TH2I>("Dalitz_p2pi", "Dalitz p#pi^{+}#pi^{-}; M_{#pi^{+}#pi^{-}}^{2}; M_{p#pi^{+}}^{2}", 200, 0.0, 5.0, 200, 0., 5.);
80  dMppiplus_M2pi = GetOrCreate_Histogram<TH2I>("Mppiplus_M2pi", "Psuedo-Dalitz p#pi^{+}#pi^{-}; M_{#pi^{+}#pi^{-}}; M_{p#pi^{+}}", 200, 0.0, 2.0, 200, 1.0, 3.0);
81  dMppiminus_M2pi = GetOrCreate_Histogram<TH2I>("Mppiminus_M2pi", "Psuedo-Dalitz p#pi^{+}#pi^{-}; M_{#pi^{+}#pi^{-}}; M_{p#pi^{-}}", 200, 0.0, 2.0, 200, 1.0, 3.0);
82 
83  dProtonPhi_Egamma = GetOrCreate_Histogram<TH2I>("ProtonPhi_Egamma","#phi vs E_{#gamma}; E_{#gamma}; proton #phi",endpoint_energy_bins,0,endpoint_energy,360,-180,180);
84  dPiPlusPsi_Egamma = GetOrCreate_Histogram<TH2I>("PiPlusPsi_Egamma","#psi vs E_{#gamma}; E_{#gamma}; #pi^{+} #phi",endpoint_energy_bins,0,endpoint_energy,360,-180,180);
85  dPiPlusPsi_t = GetOrCreate_Histogram<TH2I>("PiPlusPsi_t","#psi vs |t|; |t|; #pi^{+} #psi",1000,0.,5.,360,-180,180);
86 
87  dEgamma_M2pi = GetOrCreate_Histogram<TH2I>("Egamma_M2pi", "E_{#gamma} vs M_{#pi^{+}#pi^{-}}; M_{#pi^{+}#pi^{-}}; E_{#gamma}", 200, 0.0, 2.0, endpoint_energy_bins,0,endpoint_energy);
88 
89  dBaryonM_CosTheta_Egamma1 = GetOrCreate_Histogram<TH2I>("BaryonM_CosTheta_Egamma1", "Baryon M_{p#pi^{+}} vs cos#theta: E_{#gamma} < 2.5; cos#theta; M_{p#pi^{+}}", 100, -1, 1, 100, 1.0, 2.5);
90  dBaryonM_CosTheta_Egamma2 = GetOrCreate_Histogram<TH2I>("BaryonM_CosTheta_Egamma2", "Baryon M_{p#pi^{+}} vs cos#theta: 2.5 < E_{#gamma} < 3.0; cos#theta; M_{p#pi^{+}}", 100, -1, 1, 100, 1.0, 2.5);
91  dBaryonM_CosTheta_Egamma3 = GetOrCreate_Histogram<TH2I>("BaryonM_CosTheta_Egamma3", "Baryon M_{p#pi^{+}} vs cos#theta: E_{#gamma} > 3.0; cos#theta; M_{p#pi^{+}}", 100, -1, 1, 100, 1.0, 2.5);
92  }
93 
94  //Return to the base directory
96  }
97  japp->RootUnLock(); //RELEASE ROOT LOCK!!
98 }
99 
100 bool DCustomAction_p2pi_hists::Perform_Action(JEventLoop* locEventLoop, const DParticleCombo* locParticleCombo)
101 {
102 
103  const DDetectorMatches* locDetectorMatches = NULL;
104  locEventLoop->GetSingle(locDetectorMatches);
105 
106  // should only have one reaction step
107  const DParticleComboStep* locParticleComboStep = locParticleCombo->Get_ParticleComboStep(0);
108 
109  // get beam photon energy and final state particles
110  const DKinematicData* locBeamPhoton = NULL;
111  deque<const DKinematicData*> locParticles;
112  if(!Get_UseKinFitResultsFlag()) { //measured
113  locBeamPhoton = locParticleComboStep->Get_InitialParticle_Measured();
114  locParticleComboStep->Get_FinalParticles_Measured(locParticles);
115  }
116  else {
117  locBeamPhoton = locParticleComboStep->Get_InitialParticle();
118  locParticleComboStep->Get_FinalParticles(locParticles);
119  }
120  double locBeamPhotonEnergy = locBeamPhoton->energy();
121 
122  // calculate 2pi P4
123  DLorentzVector locP4_2pi;
124  for(size_t loc_i = 0; loc_i < 3; ++loc_i) {
125  if(locParticles[loc_i] == NULL) continue; // missing particle
126  if(locParticles[loc_i]->PID() == PiPlus || locParticles[loc_i]->PID() == PiMinus)
127  locP4_2pi += locParticles[loc_i]->lorentzMomentum();
128  }
129 
130  // calculate missing P4
131  DLorentzVector locMissingP4 = dAnalysisUtilities->Calc_MissingP4(locParticleCombo,Get_UseKinFitResultsFlag());
132 
133  // reconstructed proton variables
134  double dEdx = 0.;
135  DLorentzVector locProtonP4;
136  if(dMissingPID != Proton) {
137  const DChargedTrack* locChargedTrack = static_cast<const DChargedTrack*>(locParticleComboStep->Get_FinalParticle_SourceObject(0));
138  const DChargedTrackHypothesis* locChargedTrackHypothesis = locChargedTrack->Get_Hypothesis(Proton);
139  dEdx = locChargedTrackHypothesis->dEdx()*1e6; // convert to keV
140  locProtonP4 = locChargedTrackHypothesis->lorentzMomentum();
141  }
142  else locProtonP4 = locMissingP4;
143 
144  // production kinematics
145  DLorentzVector locSumInitP4;
146  DLorentzVector locProtonP4Init(0,0,0,0.938);
147  locSumInitP4 += locProtonP4Init;
148  locSumInitP4 += locBeamPhoton->lorentzMomentum();
149  TLorentzVector locDelta = (locProtonP4 - locProtonP4Init);
150  double t = locDelta.M2();
151 
152  TLorentzVector locPiPlus_P4 = locParticles[1]->lorentzMomentum();
153 
154  // boost to resonance frame for angular distributions
155  TLorentzRotation resonanceBoost( -locP4_2pi.BoostVector() );
156  TLorentzVector recoil_res = resonanceBoost * locProtonP4;
157  TLorentzVector p1_res = resonanceBoost * locPiPlus_P4;
158 
159  // normal to the production plane
160  TVector3 y = (locBeamPhoton->lorentzMomentum().Vect().Unit().Cross(-locProtonP4.Vect().Unit())).Unit();
161 
162  // choose helicity frame: z-axis opposite recoil proton in rho rest frame
163  TVector3 z = -1. * recoil_res.Vect().Unit();
164  TVector3 x = y.Cross(z).Unit();
165  TVector3 angles( (p1_res.Vect()).Dot(x),
166  (p1_res.Vect()).Dot(y),
167  (p1_res.Vect()).Dot(z) );
168 
169  double cosTheta = angles.CosTheta();
170  double phi = angles.Phi();
171 
172  TVector3 eps(1.0, 0.0, 0.0); // beam polarization vector
173  double Phi = atan2(y.Dot(eps), locBeamPhoton->lorentzMomentum().Vect().Unit().Dot(eps.Cross(y)));
174 
175  double locPsi = phi - Phi;
176  if(locPsi < -1*TMath::Pi()) locPsi += 2*TMath::Pi();
177  if(locPsi > TMath::Pi()) locPsi -= 2*TMath::Pi();
178 
179  //FILL HISTOGRAMS
180  //Since we are filling histograms local to this action, it will not interfere with other ROOT operations: can use action-wide ROOT lock
181  //Note, the mutex is unique to this DReaction + action_string combo: actions of same class with different hists will have a different mutex
182  Lock_Action(); //ACQUIRE ROOT LOCK!!
183  {
184  // Fill histograms here
185  dEgamma->Fill(locBeamPhotonEnergy);
186 
187  if(dMissingPID == Proton) {
188  if(locBeamPhotonEnergy > cohmin_energy && locBeamPhotonEnergy < cohedge_energy) {
189  dMM_M2pi->Fill(locP4_2pi.M(), locMissingP4.M());
190 
191  if(locParticles[1]->lorentzMomentum().Theta()*180/TMath::Pi() > 2. && locParticles[2]->lorentzMomentum().Theta()*180/TMath::Pi() > 2.) {
192  dMM_M2pi_noEle->Fill(locP4_2pi.M(), locMissingP4.M());
193 
194  if(locMissingP4.M() > minMMCut && locMissingP4.M() < maxMMCut) {
195  dt_M2pi_noEle->Fill(locP4_2pi.M(), fabs(t));
196 
197  if(locP4_2pi.M() > minRhoMassCut && locP4_2pi.M() < maxRhoMassCut){
198  dPiPlusPsi_t->Fill(fabs(t), locPsi*180/TMath::Pi());
199  }
200  }
201  }
202  }
203  }
204  else {
205  dMM2_M2pi->Fill(locP4_2pi.M(), locMissingP4.M2());
206  dDeltaE_M2pi->Fill(locP4_2pi.M(),locMissingP4.E());
207 
208  if(locMissingP4.M2() > minMM2Cut && locMissingP4.M2() < maxMM2Cut && fabs(locMissingP4.E()) < missingEnergyCut){
209  dProton_dEdx_P->Fill(locProtonP4.Vect().Mag(), dEdx);
210  dProton_P_Theta->Fill(locProtonP4.Vect().Theta()*180/TMath::Pi(), locProtonP4.Vect().Mag());
211 
212  if(dEdx > dEdxCut) {
213 
214  dEgamma_M2pi->Fill(locP4_2pi.M(), locBeamPhotonEnergy);
215  dDeltaE_M2pi_ProtonTag->Fill(locP4_2pi.M(),locMissingP4.E());
216 
217  DLorentzVector locP4_ppiplus = locProtonP4;
218  DLorentzVector locP4_ppiminus = locProtonP4;
219  locP4_ppiplus += locParticles[1]->lorentzMomentum();
220  locP4_ppiminus += locParticles[2]->lorentzMomentum();
221  if(locBeamPhotonEnergy > cohmin_energy && locBeamPhotonEnergy < cohedge_energy){
222  dDalitz_p2pi->Fill(locP4_2pi.M2(), locP4_ppiplus.M2());
223  dMppiplus_M2pi->Fill(locP4_2pi.M(), locP4_ppiplus.M());
224  dMppiminus_M2pi->Fill(locP4_2pi.M(), locP4_ppiminus.M());
225  }
226 
227  // some p pi+ mass distributions (hunt for Delta++...)
228  if(locBeamPhotonEnergy < cohmin_energy)
229  dBaryonM_CosTheta_Egamma1->Fill(cosTheta,locP4_ppiplus.M());
230  else if(locBeamPhotonEnergy < cohedge_energy)
231  dBaryonM_CosTheta_Egamma2->Fill(cosTheta,locP4_ppiplus.M());
232  else
233  dBaryonM_CosTheta_Egamma3->Fill(cosTheta,locP4_ppiplus.M());
234 
235 
236  if(locP4_2pi.M() > minRhoMassCut && locP4_2pi.M() < maxRhoMassCut){
237  dProtonPhi_Egamma->Fill(locBeamPhotonEnergy, locProtonP4.Phi()*180/TMath::Pi());
238  dPiPlusPsi_Egamma->Fill(locBeamPhotonEnergy, locPsi*180/TMath::Pi());
239 
240  if(locBeamPhotonEnergy > cohmin_energy && locBeamPhotonEnergy < cohedge_energy){
241  dPiPlusPsi_t->Fill(fabs(t), locPsi*180/TMath::Pi());
242  }
243  }
244  }
245  }
246  }
247  }
248  Unlock_Action(); //RELEASE ROOT LOCK!!
249 
250  return true; //return false if you want to use this action to apply a cut (and it fails the cut!)
251 }
vector< const DKinematicData * > Get_FinalParticles_Measured(void) const
const DAnalysisUtilities * dAnalysisUtilities
TDirectoryFile * ChangeTo_BaseDirectory(void)
DApplication * dapp
double energy(void) const
const DReaction * Get_Reaction(void) const
Double_t x[NCHANNELS]
Definition: st_tw_resols.C:39
#define y
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
DLorentzVector Calc_MissingP4(const DReaction *locReaction, const DParticleCombo *locParticleCombo, bool locUseKinFitDataFlag) const
void Initialize(JEventLoop *locEventLoop)
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
const DParticleComboStep * Get_ParticleComboStep(size_t locStepIndex) const
bool Perform_Action(JEventLoop *locEventLoop, const DParticleCombo *locParticleCombo)