** Warning **
Issuing rollback() due to DESTROY without explicit disconnect() of DBD::mysql::db handle dbname=lxr_eic at /usr/local/share/lxr/lxr-2.3.7/lib/LXR/Common.pm line 1161, <GEN5> line 1.
Last-Modified: Wed, 16 Sep 2025 09:15:53 GMT
Content-Type: text/html; charset=utf-8
/master/detector_benchmarks/benchmarks/zdc_neutron/analysis/fwd_neutrons_recon.C
File indexing completed on 2025-09-16 08:16:52
0001
0002 void fwd_neutrons_recon (std ::string inputfile , std ::string outputfile ){
0003
0004
0005 gStyle ->SetOptStat (0);
0006 gStyle ->SetPadBorderMode (0);
0007 gStyle ->SetFrameBorderMode (0);
0008 gStyle ->SetFrameLineWidth (2);
0009 gStyle ->SetLabelSize (0.035,"X" );
0010 gStyle ->SetLabelSize (0.035,"Y" );
0011
0012
0013 gStyle ->SetTitleXSize (0.04);
0014 gStyle ->SetTitleXOffset (0.9);
0015 gStyle ->SetTitleYSize (0.04);
0016 gStyle ->SetTitleYOffset (0.9);
0017
0018
0019 TH2 *h1_neut = new TH2D ("h1_neut" ,"Neutron true energy vs. polar angle" ,100,0.6,2.2,100,0,200);
0020 h1_neut ->GetXaxis ()->SetTitle ("#theta [deg]" ); h1_neut ->GetXaxis ()->CenterTitle ();
0021 h1_neut ->GetYaxis ()->SetTitle ("E [GeV]" ); h1_neut ->GetYaxis ()->CenterTitle ();
0022
0023 TH2 *h2_neut = new TH2D ("h2_neut" ,"Neutron true azimuthal angle vs. polar angle around p axis" ,100,-0.1,12,100,-200,200);
0024 h2_neut ->GetXaxis ()->SetTitle ("#theta^{*} [mRad]" ); h2_neut ->GetXaxis ()->CenterTitle ();
0025 h2_neut ->GetYaxis ()->SetTitle ("#phi^{*} [deg]" ); h2_neut ->GetYaxis ()->CenterTitle ();
0026
0027 TH2 *h2a_neut = new TH2D ("h2a_neut" ,"Neutron true azimuthal angle vs. cosine of polar angle around p axis" ,100,0.99996,1,100,-200,200);
0028 h2a_neut ->GetXaxis ()->SetTitle ("cos(#theta^{*})" ); h2a_neut ->GetXaxis ()->CenterTitle ();
0029 h2a_neut ->GetYaxis ()->SetTitle ("#phi^{*} [deg]" ); h2a_neut ->GetYaxis ()->CenterTitle ();
0030
0031 TH1 *h1_ecal_adc = new TH1D ("h1_ecal_adc" ,"ECal ADC amplitude spectrum" ,1000,-0.1,35000);
0032 h1_ecal_adc ->GetXaxis ()->SetTitle ("ADC Channel" ); h1_ecal_adc ->GetXaxis ()->CenterTitle ();
0033 h1_ecal_adc ->SetLineColor (kRed );h1_ecal_adc ->SetLineWidth (2);
0034
0035 TH1 *h1_hcal_adc = new TH1D ("h1_hcal_adc" ,"HCal ADC amplitude spectrum" ,1000,-0.1,35000);
0036 h1_hcal_adc ->GetXaxis ()->SetTitle ("ADC Channel" ); h1_hcal_adc ->GetXaxis ()->CenterTitle ();
0037 h1_hcal_adc ->SetLineColor (kBlue );h1_hcal_adc ->SetLineWidth (2);
0038
0039 TH1 *h1_ecal = new TH1D ("h1_ecal" ,"Total reconstructed hit energy sum" ,100,-0.1,2);
0040 h1_ecal ->GetXaxis ()->SetTitle ("E [GeV]" ); h1_ecal ->GetXaxis ()->CenterTitle ();
0041 h1_ecal ->SetLineColor (kRed );h1_ecal ->SetLineWidth (2);
0042
0043 TH1 *h1_hcal = new TH1D ("h1_hcal" ,"Total reconstructed hit energy sum" ,100,-0.1,2);
0044 h1_hcal ->GetXaxis ()->SetTitle ("E [GeV]" ); h1_hcal ->GetXaxis ()->CenterTitle ();
0045 h1_hcal ->SetLineColor (kBlue );h1_hcal ->SetLineWidth (2);
0046
0047
0048 TH1 *h2_ecal = new TH1D ("h2_ecal" ,"Total reconstructed hit energy sum" ,100,-0.1,2);
0049 h2_ecal ->GetXaxis ()->SetTitle ("E [GeV]" ); h2_ecal ->GetXaxis ()->CenterTitle ();
0050 h2_ecal ->SetLineColor (kRed );h2_ecal ->SetLineWidth (2);
0051
0052
0053 TH1 *h2_hcal = new TH1D ("h2_hcal" ,"Total reconstructed hit energy sum" ,100,-0.1,2);
0054 h2_hcal ->GetXaxis ()->SetTitle ("E [GeV]" ); h2_hcal ->GetXaxis ()->CenterTitle ();
0055 h2_hcal ->SetLineColor (kBlue );h2_hcal ->SetLineWidth (2);
0056
0057
0058 TH1 *h3_ecal = new TH1D ("h3_ecal" ,"Number of reconstructed clusters in ECal" ,5,0,5);
0059 h3_ecal ->GetXaxis ()->SetTitle ("Number of clusters" ); h3_ecal ->GetXaxis ()->CenterTitle ();
0060 h3_ecal ->GetXaxis ()->SetNdivisions (405);h3_ecal ->GetXaxis ()->CenterLabels ();
0061 h3_ecal ->SetLineColor (kRed );h3_ecal ->SetLineWidth (2);
0062
0063
0064 TH1 *h3_hcal = new TH1D ("h3_hcal" ,"Number of reconstructed clusters in HCal" ,5,0,5);
0065 h3_hcal ->GetXaxis ()->SetTitle ("Number of clusters" ); h3_hcal ->GetXaxis ()->CenterTitle ();
0066 h3_hcal ->GetXaxis ()->SetNdivisions (405);h3_hcal ->GetXaxis ()->CenterLabels ();
0067 h3_hcal ->SetLineColor (kBlue );h3_hcal ->SetLineWidth (2);
0068
0069
0070 TH1 *h4_hcal = new TH1D ("h4_hcal" ,"HCal cluster energy: Events w/ 1 HCal Cluster" ,100,-0.1,120);
0071 h4_hcal ->GetXaxis ()->SetTitle ("E [GeV]" ); h4_hcal ->GetXaxis ()->CenterTitle ();
0072 h4_hcal ->SetLineColor (kBlue );h4_hcal ->SetLineWidth (2);
0073
0074
0075 TH1 *h1_neut_rec = new TH1D ("h1_neut_rec" ,"Reconstructed Neutron Energy" ,100,-0.1,120);
0076 h1_neut_rec ->GetXaxis ()->SetTitle ("E [GeV]" ); h1_neut_rec ->GetXaxis ()->CenterTitle ();
0077 h1_neut_rec ->SetLineColor (kBlue );h1_neut_rec ->SetLineWidth (2);
0078
0079
0080 TH1 *h2_neut_rec = new TH2D ("h2_neut_rec" ,"Reconstructed Neutron local hit position at ZDC HCal front face" ,100,-200,200,100,-200,200);
0081 h2_neut_rec ->GetXaxis ()->SetTitle ("x [mm]" ); h2_neut_rec ->GetXaxis ()->CenterTitle ();
0082 h2_neut_rec ->GetYaxis ()->SetTitle ("y [mm]" ); h2_neut_rec ->GetYaxis ()->CenterTitle ();
0083
0084
0085 TFile * file = new TFile (inputfile .c_str ());
0086 TTree *tree = (TTree *) file ->Get ("events" );
0087
0088 cout <<"Total number of events to analyze is " <<tree ->GetEntries ()<<endl ;
0089
0090
0091 TTreeReader tr (tree );
0092
0093 TTreeReaderArray <int > gen_status (tr , "MCParticles.generatorStatus" );
0094 TTreeReaderArray <int > gen_pid (tr , "MCParticles.PDG" );
0095 TTreeReaderArray <double > gen_px (tr , "MCParticles.momentum.x" );
0096 TTreeReaderArray <double > gen_py (tr , "MCParticles.momentum.y" );
0097 TTreeReaderArray <double > gen_pz (tr , "MCParticles.momentum.z" );
0098 TTreeReaderArray <double > gen_mass (tr , "MCParticles.mass" );
0099 TTreeReaderArray <float > gen_charge (tr , "MCParticles.charge" );
0100 TTreeReaderArray <double > gen_vx (tr , "MCParticles.vertex.x" );
0101 TTreeReaderArray <double > gen_vy (tr , "MCParticles.vertex.y" );
0102 TTreeReaderArray <double > gen_vz (tr , "MCParticles.vertex.z" );
0103
0104
0105 TTreeReaderArray <int > ecal_adc (tr ,"EcalFarForwardZDCRawHits.amplitude" );
0106 TTreeReaderArray <float > ecal_hit_e (tr ,"EcalFarForwardZDCRecHits.energy" );
0107 TTreeReaderArray <float > ecal_cluster_energy (tr , "EcalFarForwardZDCClusters.energy" );
0108
0109
0110 TTreeReaderArray <int > hcal_adc (tr ,"HcalFarForwardZDCRawHits.amplitude" );
0111 TTreeReaderArray <float > hcal_hit_e (tr , "HcalFarForwardZDCRecHits.energy" );
0112 TTreeReaderArray <float > hcal_cluster_energy (tr , "HcalFarForwardZDCClusters.energy" );
0113
0114
0115 TTreeReaderArray <float > rec_neutron_energy (tr ,"ReconstructedFarForwardZDCNeutrals.energy" );
0116 TTreeReaderArray <float > rec_neutron_px (tr ,"ReconstructedFarForwardZDCNeutrals.momentum.x" );
0117 TTreeReaderArray <float > rec_neutron_py (tr ,"ReconstructedFarForwardZDCNeutrals.momentum.y" );
0118 TTreeReaderArray <float > rec_neutron_pz (tr ,"ReconstructedFarForwardZDCNeutrals.momentum.z" );
0119 TTreeReaderArray <int > rec_neutron_PDG (tr ,"ReconstructedFarForwardZDCNeutrals.PDG" );
0120
0121 int counter (0);
0122
0123 TLorentzVector neut_true ;
0124 TLorentzVector neut_true_rot ;
0125
0126 float ecal_e_tot (0);
0127 float hcal_e_tot (0);
0128
0129
0130 while (tr .Next ()) {
0131
0132 if (counter %100==0) cout <<"Analyzing event " <<counter <<endl ;
0133 counter ++;
0134
0135
0136 ecal_e_tot = 0;
0137 hcal_e_tot = 0;
0138
0139
0140 for (int igen =0;igen <gen_status .GetSize ();igen ++){
0141 if (gen_status [igen ]==1 && gen_pid [igen ]==2112){
0142
0143 neut_true .SetXYZM (gen_px [igen ],gen_py [igen ],gen_pz [igen ],gen_mass [igen ]);
0144 h1_neut ->Fill (neut_true .Theta ()*TMath ::RadToDeg (),neut_true .E ());
0145
0146
0147 neut_true_rot = neut_true ;
0148 neut_true_rot .RotateY (0.025);
0149 h2_neut ->Fill (neut_true_rot .Theta ()*1000.,neut_true_rot .Phi ()*TMath ::RadToDeg ());
0150 h2a_neut ->Fill (std ::cos (neut_true_rot .Theta ()),neut_true_rot .Phi ()*TMath ::RadToDeg ());
0151
0152 }
0153 }
0154
0155
0156 for (int iadc =0;iadc <ecal_adc .GetSize ();iadc ++){
0157 h1_ecal_adc ->Fill (ecal_adc [iadc ]);
0158 }
0159
0160
0161 for (int ihit =0;ihit <ecal_hit_e .GetSize ();ihit ++){
0162 ecal_e_tot += ecal_hit_e [ihit ];
0163 }
0164
0165
0166 for (int iadc =0;iadc <hcal_adc .GetSize ();iadc ++){
0167 h1_hcal_adc ->Fill (hcal_adc [iadc ]);
0168 }
0169
0170
0171 for (int ihit =0;ihit <hcal_hit_e .GetSize ();ihit ++){
0172 hcal_e_tot += hcal_hit_e [ihit ];
0173 }
0174
0175
0176 int ecal_clus_size = ecal_cluster_energy .GetSize ();
0177
0178
0179 int hcal_clus_size = hcal_cluster_energy .GetSize ();
0180
0181
0182 h1_ecal ->Fill (ecal_e_tot );
0183 h1_hcal ->Fill (hcal_e_tot );
0184
0185 if ( neut_true_rot .Theta ()*1000. < 3.5 ){
0186 h2_ecal ->Fill (ecal_e_tot );
0187 h2_hcal ->Fill (hcal_e_tot );
0188
0189 h3_ecal ->Fill (ecal_clus_size );
0190 h3_hcal ->Fill (hcal_clus_size );
0191
0192
0193 if (hcal_clus_size ==1) h4_hcal ->Fill (hcal_cluster_energy [0]);
0194 }
0195
0196
0197 for (int irec =0;irec <rec_neutron_energy .GetSize ();irec ++){
0198 if (rec_neutron_PDG [irec ] != 2112)
0199 continue ;
0200 if (neut_true_rot .Theta ()*1000. < 3.5 && hcal_clus_size ==1){
0201 h1_neut_rec ->Fill (rec_neutron_energy [irec ]);
0202
0203 TVector3 neut_rec (rec_neutron_px [irec ],rec_neutron_py [irec ],rec_neutron_pz [irec ]);
0204
0205
0206 TVector3 neut_rec_rot = neut_rec ;
0207 neut_rec_rot .RotateY (0.025);
0208
0209
0210 auto proj_x = 35.8 * 1000. * tan (neut_rec_rot .Theta ()) * cos (neut_rec_rot .Phi ()) ;
0211 auto proj_y = 35.8 * 1000. * tan (neut_rec_rot .Theta ()) * sin (neut_rec_rot .Phi ()) ;
0212 h2_neut_rec ->Fill (proj_x ,proj_y );
0213 }
0214 }
0215
0216 }
0217
0218
0219 TCanvas *c1 = new TCanvas ("c1" );
0220 h1_neut ->Draw ("colz" );
0221
0222 TCanvas *c2 = new TCanvas ("c2" );
0223 h2_neut ->Draw ("colz" );
0224
0225 TCanvas *c2a = new TCanvas ("c2a" );
0226 h2a_neut ->Draw ("colz" );
0227
0228 TCanvas *c3a = new TCanvas ("c3a" );
0229 c3a ->SetLogy ();
0230 h1_ecal_adc ->Draw ();
0231
0232 TCanvas *c3b = new TCanvas ("c3b" );
0233 c3b ->SetLogy ();
0234 h1_hcal_adc ->Draw ();
0235
0236 TCanvas *c4 = new TCanvas ("c4" );
0237 h1_ecal ->Draw ("" );
0238 h1_hcal ->Draw ("same" );
0239
0240 TLegend *leg4 = new TLegend (0.6,0.6,0.85,0.8);
0241 leg4 ->SetBorderSize (0);leg4 ->SetFillStyle (0);
0242 leg4 ->AddEntry (h1_ecal ,"Sum of digitized ZDC Ecal hit energies" ,"l" );
0243 leg4 ->AddEntry (h1_hcal ,"Sum of digitized ZDC Hcal hit energies" ,"l" );
0244 leg4 ->Draw ();
0245
0246 TCanvas *c5 = new TCanvas ("c5" );
0247 h2_ecal ->Draw ("" );
0248 h2_hcal ->Draw ("same" );
0249
0250 TLegend *leg5 = new TLegend (0.5,0.6,0.9,0.8);
0251 leg5 ->SetBorderSize (0);leg5 ->SetFillStyle (0);
0252 leg5 ->SetHeader ("Require neutron #theta^{*} (wrt proton beam) < 3.5 mRad" );
0253 leg5 ->AddEntry (h1_ecal ,"Sum of digitized ZDC Ecal hit energies" ,"l" );
0254 leg5 ->AddEntry (h1_hcal ,"Sum of digitized ZDC Hcal hit energies" ,"l" );
0255 leg5 ->Draw ();
0256
0257 TCanvas *c6a = new TCanvas ("c6a" );
0258 h3_ecal ->Draw ();
0259
0260 TLegend *leg6 = new TLegend (0.5,0.7,0.9,0.9);
0261 leg6 ->SetBorderSize (0);leg6 ->SetFillStyle (0);
0262 leg6 ->SetHeader ("Require neutron #theta^{*} (wrt proton beam) < 3.5 mRad" );
0263 leg6 ->Draw ();
0264
0265 TCanvas *c6b = new TCanvas ("c6b" );
0266 h3_hcal ->Draw ();
0267 leg6 ->Draw ();
0268
0269 TCanvas *c7 = new TCanvas ("c7" );
0270 h4_hcal ->Draw ();
0271
0272 TLegend *leg7 = new TLegend (0.15,0.7,0.5,0.9);
0273 leg7 ->SetBorderSize (0);leg7 ->SetFillStyle (0);
0274 leg7 ->SetHeader ("Require neutron #theta^{*} (wrt proton beam) < 3.5 mRad" );
0275 leg7 ->Draw ();
0276
0277 TCanvas *c8 = new TCanvas ("c8" );
0278 h1_neut_rec ->Draw ();
0279
0280 TLegend *leg8 = new TLegend (0.15,0.7,0.7,0.9);
0281 leg8 ->SetBorderSize (0);leg8 ->SetFillStyle (0);
0282 leg8 ->SetHeader ("Generated neutron #theta^{*} < 3.5 mRad + 1 HCal cluster" );
0283 leg8 ->Draw ();
0284
0285 TCanvas *c9 = new TCanvas ("c9" );
0286 h2_neut_rec ->Draw ("colz" );
0287
0288 TLatex *tex9 = new TLatex (-150,150,"Generated neutron #theta^{*} < 3.5 mRad + 1 HCal cluster" );
0289 tex9 ->SetTextSize (0.03);
0290 tex9 ->Draw ();
0291
0292
0293 c1 ->Print (Form ("%s[" ,outputfile .c_str ()));
0294 c1 ->Print (Form ("%s" ,outputfile .c_str ()));
0295 c2 ->Print (Form ("%s" ,outputfile .c_str ()));
0296 c2a ->Print (Form ("%s" ,outputfile .c_str ()));
0297 c3a ->Print (Form ("%s" ,outputfile .c_str ()));
0298 c3b ->Print (Form ("%s" ,outputfile .c_str ()));
0299 c4 ->Print (Form ("%s" ,outputfile .c_str ()));
0300 c5 ->Print (Form ("%s" ,outputfile .c_str ()));
0301 c6a ->Print (Form ("%s" ,outputfile .c_str ()));
0302 c6b ->Print (Form ("%s" ,outputfile .c_str ()));
0303 c7 ->Print (Form ("%s" ,outputfile .c_str ()));
0304 c8 ->Print (Form ("%s" ,outputfile .c_str ()));
0305 c9 ->Print (Form ("%s" ,outputfile .c_str ()));
0306 c9 ->Print (Form ("%s]" ,outputfile .c_str ()));
0307
0308 }