Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:13:17

0001 //____________________________________________________________________________..
0002 
0003 
0004 #include "JetKinematicCheck.h"
0005 #include <fun4all/Fun4AllReturnCodes.h>
0006 #include <phool/PHCompositeNode.h>
0007 
0008 #include <fun4all/PHTFileServer.h>
0009 #include <phool/getClass.h>
0010 
0011 #include <jetbase/JetContainer.h>
0012 #include <jetbase/Jetv2.h>
0013 
0014 
0015 #include <TFile.h>
0016 #include <TH3D.h>
0017 #include <TH2D.h>
0018 #include <TH1D.h>
0019 #include <TPad.h>
0020 #include <TLegend.h>
0021 #include <cmath>
0022 #include <string>
0023 #include <vector>
0024 
0025 
0026 //____________________________________________________________________________..
0027 
0028 JetKinematicCheck::JetKinematicCheck(const std::string& recojetnameR02, const std::string& recojetnameR03, const std::string& recojetnameR04, const std::string& outputfilename):
0029 SubsysReco("JetKinematicCheck")
0030   , m_recoJetNameR02(recojetnameR02)
0031   , m_recoJetNameR03(recojetnameR03)
0032   , m_recoJetNameR04(recojetnameR04)
0033   , m_outputFileName(outputfilename)
0034   , m_etaRange(-1.1, 1.1)
0035   , m_ptRange(10, 100)
0036 
0037 
0038 
0039 {
0040   std::cout << "JetKinematicCheck::JetKinematicCheck(const std::string &name) Calling ctor" << std::endl;
0041 }
0042 
0043 //____________________________________________________________________________..
0044 JetKinematicCheck::~JetKinematicCheck()
0045 {
0046   std::cout << "JetKinematicCheck::~JetKinematicCheck() Calling dtor" << std::endl;
0047 }
0048 
0049 //____________________________________________________________________________..
0050 int JetKinematicCheck::Init(PHCompositeNode *topNode)
0051 {
0052 
0053 
0054   PHTFileServer::get().open(m_outputFileName, "RECREATE");
0055 
0056   // initialize histograms
0057 
0058   jet_spectra_r02 = new TH1D("h_spectra_r02", "", 19, 5, 100);
0059   jet_spectra_r02->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0060 
0061   jet_spectra_r03 = new TH1D("h_spectra_r03", "", 19, 5, 100);
0062   jet_spectra_r03->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0063 
0064   jet_spectra_r04 = new TH1D("h_spectra_r04", "", 19, 5, 100);
0065   jet_spectra_r04->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0066 
0067   jet_eta_phi_r02 = new TH2D("h_eta_phi_r02","", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0068   jet_eta_phi_r02->GetXaxis()->SetTitle("#eta");
0069   jet_eta_phi_r02->GetYaxis()->SetTitle("#Phi");
0070 
0071   jet_eta_phi_r03 = new TH2D("h_eta_phi_r03","", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0072   jet_eta_phi_r03->GetXaxis()->SetTitle("#eta");
0073   jet_eta_phi_r03->GetYaxis()->SetTitle("#Phi");
0074  
0075   jet_eta_phi_r04 = new TH2D("h_eta_phi_r04","", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0076   jet_eta_phi_r04->GetXaxis()->SetTitle("#eta");
0077   jet_eta_phi_r04->GetYaxis()->SetTitle("#Phi");
0078  
0079   jet_mass_pt_r02 = new TH2D("h_jet_mass_pt_r02","", 19, 5, 100, 15, 0, 15);
0080   jet_mass_pt_r02->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0081   jet_mass_pt_r02->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}]");
0082 
0083   jet_mass_pt_r03 = new TH2D("h_jet_mass_pt_r03","", 19, 5, 100, 15, 0, 15);
0084   jet_mass_pt_r03->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0085   jet_mass_pt_r03->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}]");
0086 
0087   jet_mass_pt_r04 = new TH2D("h_jet_mass_pt_r04","", 19, 5, 100, 15, 0, 15);
0088   jet_mass_pt_r04->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0089   jet_mass_pt_r04->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}]");
0090 
0091   jet_mass_eta_r02 = new TH2D("h_jet_mass_eta_r02","", 24, -1.1, 1.1, 15, 0, 15);
0092   jet_mass_eta_r02->GetXaxis()->SetTitle("#eta");
0093   jet_mass_eta_r02->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}]");
0094 
0095   jet_mass_eta_r03 = new TH2D("h_jet_mass_eta_r03","", 24, -1.1, 1.1, 15, 0, 15);
0096   jet_mass_eta_r03->GetXaxis()->SetTitle("#eta");
0097   jet_mass_eta_r03->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}");
0098 
0099   jet_mass_eta_r04 = new TH2D("h_jet_mass_eta_r04","", 24, -1.1, 1.1, 15, 0, 15);
0100   jet_mass_eta_r04->GetXaxis()->SetTitle("#eta");
0101   jet_mass_eta_r04->GetYaxis()->SetTitle("Jet Mass [GeV/c^{2}]");
0102 
0103 
0104 
0105   std::cout << "JetKinematicCheck::Init(PHCompositeNode *topNode) Initializing" << std::endl;
0106   return Fun4AllReturnCodes::EVENT_OK;
0107 }
0108 
0109 
0110 //____________________________________________________________________________..
0111 int JetKinematicCheck::InitRun(PHCompositeNode *topNode)
0112 {
0113   std::cout << "JetKinematicCheck::InitRun(PHCompositeNode *topNode) Initializing for Run XXX" << std::endl;
0114   return Fun4AllReturnCodes::EVENT_OK;
0115 }
0116 
0117 
0118 
0119 
0120 //____________________________________________________________________________..
0121 int JetKinematicCheck::process_event(PHCompositeNode *topNode)
0122 {
0123 
0124   std::vector<std::string> m_recoJetName_array = {m_recoJetNameR02, m_recoJetNameR03, m_recoJetNameR04};
0125   m_radii = {0.2, 0.3, 0.4};
0126   int n_radii = m_radii.size();
0127 
0128   // Loop over each reco jet radii from array
0129   for(int i = 0; i < n_radii; i++){
0130 
0131     std::string recoJetName = m_recoJetName_array[i];
0132 
0133     JetContainer* jets = findNode::getClass<JetContainer>(topNode, recoJetName);
0134     if (!jets)
0135       {
0136     std::cout
0137       << "MyJetAnalysis::process_event - Error can not find DST Reco JetContainer node "
0138       << recoJetName << std::endl;
0139     exit(-1);
0140       }
0141 
0142 
0143     //loop over jets
0144     for (auto jet : *jets)
0145       {
0146     bool eta_cut = (jet->get_eta() >= m_etaRange.first) and (jet->get_eta() <= m_etaRange.second);
0147     bool pt_cut = (jet->get_pt() >= m_ptRange.first) and (jet->get_pt() <= m_ptRange.second);
0148     if ((not eta_cut) or (not pt_cut)) continue;
0149     if(jet->get_pt() < 1) continue; // to remove noise jets
0150 
0151     if(i == 0){
0152       jet_spectra_r02->Fill(jet->get_pt());
0153       jet_eta_phi_r02->Fill(jet->get_eta(), jet->get_phi());
0154       jet_mass_pt_r02->Fill(jet->get_pt(), jet->get_mass());
0155       jet_mass_eta_r02->Fill(jet->get_eta(), jet->get_mass());
0156     }
0157 
0158     else if(i == 1){
0159       jet_spectra_r03->Fill(jet->get_pt());
0160       jet_eta_phi_r03->Fill(jet->get_eta(), jet->get_phi());
0161       jet_mass_pt_r03->Fill(jet->get_pt(), jet->get_mass());
0162       jet_mass_eta_r03->Fill(jet->get_eta(), jet->get_mass());
0163     
0164     }
0165 
0166     else if(i == 2){
0167       jet_spectra_r04->Fill(jet->get_pt());
0168       jet_eta_phi_r04->Fill(jet->get_eta(), jet->get_phi());
0169       jet_mass_pt_r04->Fill(jet->get_pt(), jet->get_mass());
0170       jet_mass_eta_r04->Fill(jet->get_eta(), jet->get_mass());
0171     }
0172       }
0173   }
0174 
0175 
0176   std::cout << "JetKinematicCheck::process_event(PHCompositeNode *topNode) Processing Event" << std::endl;
0177   return Fun4AllReturnCodes::EVENT_OK;
0178 }
0179 
0180 int count = 0;
0181 
0182 //____________________________________________________________________________..
0183 int JetKinematicCheck::ResetEvent(PHCompositeNode *topNode)
0184 {
0185   std::cout << "JetKinematicCheck::ResetEvent(PHCompositeNode *topNode) Resetting internal structures, prepare for next event" << std::endl;
0186 
0187   std:: cout << count << std::endl;
0188   count++;
0189 
0190   return Fun4AllReturnCodes::EVENT_OK;
0191 }
0192 
0193 //____________________________________________________________________________..
0194 int JetKinematicCheck::EndRun(const int runnumber)
0195 {
0196   std::cout << "JetKinematicCheck::EndRun(const int runnumber) Ending Run for Run " << runnumber << std::endl;
0197   return Fun4AllReturnCodes::EVENT_OK;
0198 }
0199 
0200 //____________________________________________________________________________..
0201 int JetKinematicCheck::End(PHCompositeNode *topNode)
0202 {
0203 
0204   std::cout << "JetKinematicCheck::End - Output to " << m_outputFileName << std::endl;
0205   PHTFileServer::get().cd(m_outputFileName);
0206 
0207 
0208   //for jet spectra [R02]
0209   TLegend *leg1 = new TLegend(.7,.9,.9,1);
0210   leg1->SetFillStyle(0);
0211   leg1->SetBorderSize(0);
0212   leg1->SetTextSize(0.06);
0213   leg1->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0214   leg1->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0215   jet_spectra_r02->SetMarkerStyle(kFullCircle);
0216   jet_spectra_r02->SetMarkerColor(kBlack);
0217   jet_spectra_r02->SetLineColor(kBlack);
0218   jet_spectra_r02->SetTitle("Jet Spectra [R = 0.2]");
0219   jet_spectra_r02->GetYaxis()->SetTitle("Counts");
0220   jet_spectra_r02->GetListOfFunctions()->Add(leg1);
0221   jet_spectra_r02->SetStats(0);
0222   jet_spectra_r02->Write("Jet_Spectra_R02_pp");
0223 
0224 
0225   //for jet spectra [R03]
0226   TLegend *leg2 = new TLegend(.7,.9,.9,1);
0227   leg2->SetFillStyle(0);
0228   leg2->SetBorderSize(0);
0229   leg2->SetTextSize(0.06);
0230   leg2->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0231   leg2->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0232   jet_spectra_r03->SetMarkerStyle(kFullCircle);
0233   jet_spectra_r03->SetMarkerColor(kBlack);
0234   jet_spectra_r03->SetLineColor(kBlack);
0235   jet_spectra_r03->SetTitle("Jet Spectra [R = 0.3]");
0236   jet_spectra_r03->GetYaxis()->SetTitle("Counts");
0237   jet_spectra_r03->GetListOfFunctions()->Add(leg2);
0238   jet_spectra_r03->SetStats(0);
0239   jet_spectra_r03->Write("Jet_Spectra_R03_pp");
0240  
0241   //for jet spectra [R04]
0242   TLegend *leg3 = new TLegend(.7,.9,.9,1);
0243   leg3->SetFillStyle(0);
0244   leg3->SetBorderSize(0);
0245   leg3->SetTextSize(0.06);
0246   leg3->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0247   leg3->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0248   jet_spectra_r04->SetMarkerStyle(kFullCircle);
0249   jet_spectra_r04->SetMarkerColor(kBlack);
0250   jet_spectra_r04->SetLineColor(kBlack);
0251   jet_spectra_r04->SetTitle("Jet Spectra [R = 0.4]");
0252   jet_spectra_r04->GetYaxis()->SetTitle("Counts");
0253   jet_spectra_r04->GetListOfFunctions()->Add(leg3);
0254   jet_spectra_r04->SetStats(0);
0255   jet_spectra_r04->Write("Jet_Spectra_R04_pp");
0256  
0257    
0258   //for jet eta-phi [R02]
0259   TLegend *leg4 = new TLegend(.7,.9,.9,1);
0260   leg4->SetFillStyle(0);
0261   leg4->SetBorderSize(0);
0262   leg4->SetTextSize(0.06);
0263   leg4->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0264   leg4->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),""); 
0265   jet_eta_phi_r02->SetStats(0);
0266   jet_eta_phi_r02->SetTitle("Jet Eta-Phi [R = 0.2]");
0267   jet_eta_phi_r02->GetListOfFunctions()->Add(leg4);
0268   jet_eta_phi_r02->Write("Jet_Eta_Phi_R02_pp");
0269 
0270  
0271   //for jet eta-phi [R03]
0272   TLegend *leg5 = new TLegend(.7,.9,.9,1);
0273   leg5->SetFillStyle(0);
0274   leg5->SetBorderSize(0);
0275   leg5->SetTextSize(0.06);
0276   leg5->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0277   leg5->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),""); 
0278   jet_eta_phi_r03->SetStats(0);
0279   jet_eta_phi_r03->SetTitle("Jet Eta-Phi [R = 0.3]");
0280   jet_eta_phi_r03->GetListOfFunctions()->Add(leg5);
0281   jet_eta_phi_r03->Write("Jet_Eta_Phi_R03_pp");
0282   
0283   //for jet eta-phi [R04]
0284   TLegend *leg6 = new TLegend(.7,.9,.9,1);
0285   leg6->SetFillStyle(0);
0286   leg6->SetBorderSize(0);
0287   leg6->SetTextSize(0.06);
0288   leg6->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0289   leg6->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0290   jet_eta_phi_r04->SetStats(0);
0291   jet_eta_phi_r04->SetTitle("Jet Eta-Phi [R = 0.4]");
0292   jet_eta_phi_r04->GetListOfFunctions()->Add(leg6);
0293   jet_eta_phi_r04->Write("Jet_Eta_Phi_R04_pp");
0294 
0295   //for jet mass vs pt [R02]
0296   TLegend *leg7 = new TLegend(.7,.9,.9,1);
0297   leg7->SetFillStyle(0);
0298   leg7->SetBorderSize(0);
0299   leg7->SetTextSize(0.06);
0300   leg7->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0301   leg7->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0302   jet_mass_pt_r02->SetStats(0);
0303   jet_mass_pt_r02->SetTitle("Jet Mass vs p_{T} [R = 0.2]");
0304   jet_mass_pt_r02->GetListOfFunctions()->Add(leg7);
0305   jet_mass_pt_r02->Write("Jet_Mass_pt_R02_pp");
0306 
0307   //for average jet mass vs pt [R02]
0308   TLegend *leg8 = new TLegend(.7,.9,.9,1);
0309   leg8->SetFillStyle(0);
0310   leg8->SetBorderSize(0);
0311   leg8->SetTextSize(0.06);
0312   leg8->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0313   leg8->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0314   jet_mass_pt_1D_r02 = (TH1D*)jet_mass_pt_r02->ProfileX();
0315   jet_mass_pt_1D_r02->SetStats(0);
0316   jet_mass_pt_1D_r02->SetTitle("Average Jet Mass vs p_{T} [R = 0.2]");
0317   jet_mass_pt_1D_r02->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0318   jet_mass_pt_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0319   jet_mass_pt_1D_r02->SetMarkerStyle(8);
0320   jet_mass_pt_1D_r02->SetMarkerColor(1);
0321   jet_mass_pt_1D_r02->SetLineColor(1);
0322   jet_mass_pt_1D_r02->GetListOfFunctions()->Add(leg8);
0323   jet_mass_pt_1D_r02->Write("Average_Jet_Mass_pt_R02_pp");
0324   
0325 
0326   //for jet mass vs pt [R03]
0327   TLegend *leg9 = new TLegend(.7,.9,.9,1);
0328   leg9->SetFillStyle(0);
0329   leg9->SetBorderSize(0);
0330   leg9->SetTextSize(0.06);
0331   leg9->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0332   leg9->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0333   jet_mass_pt_r03->SetStats(0);
0334   jet_mass_pt_r03->SetTitle("Jet Mass vs p_{T} [R = 0.3]");
0335   jet_mass_pt_r03->GetListOfFunctions()->Add(leg9);
0336   jet_mass_pt_r03->Write("Jet_Mass_pt_R03_pp");
0337 
0338 
0339   //for average jet mass vs pt [R03]
0340   TLegend *leg10 = new TLegend(.7,.9,.9,1);
0341   leg10->SetFillStyle(0);
0342   leg10->SetBorderSize(0);
0343   leg10->SetTextSize(0.06);
0344   leg10->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0345   leg10->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0346   jet_mass_pt_1D_r03 = (TH1D*)jet_mass_pt_r03->ProfileX();
0347   jet_mass_pt_1D_r03->SetStats(0);
0348   jet_mass_pt_1D_r03->SetTitle("Average Jet Mass vs p_{T} [R = 0.3]");
0349   jet_mass_pt_1D_r03->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0350   jet_mass_pt_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0351   jet_mass_pt_1D_r03->SetMarkerStyle(8);
0352   jet_mass_pt_1D_r03->SetMarkerColor(1);
0353   jet_mass_pt_1D_r03->SetLineColor(1);
0354   jet_mass_pt_1D_r03->GetListOfFunctions()->Add(leg10);
0355   jet_mass_pt_1D_r03->Write("Average_Jet_Mass_pt_R03_pp");
0356 
0357 
0358   //for jet mass vs pt [R04]
0359   TLegend *leg11 = new TLegend(.7,.9,.9,1);
0360   leg11->SetFillStyle(0);
0361   leg11->SetBorderSize(0);
0362   leg11->SetTextSize(0.06);
0363   leg11->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0364   leg11->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0365   jet_mass_pt_r04->SetStats(0);
0366   jet_mass_pt_r04->SetTitle("Jet Mass vs p_{T} [R = 0.4]");
0367   jet_mass_pt_r04->GetListOfFunctions()->Add(leg11);
0368   jet_mass_pt_r04->Write("Jet_Mass_pt_R04_pp");
0369 
0370   //for average jet mass vs pt [R04]
0371   TLegend *leg12 = new TLegend(.7,.9,.9,1);
0372   leg12->SetFillStyle(0);
0373   leg12->SetBorderSize(0);
0374   leg12->SetTextSize(0.06);
0375   leg12->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0376   leg12->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0377   jet_mass_pt_1D_r04 = (TH1D*)jet_mass_pt_r04->ProfileX();
0378   jet_mass_pt_1D_r04->SetStats(0);
0379   jet_mass_pt_1D_r04->SetTitle("Average Jet Mass vs p_{T} [R = 0.4]");
0380   jet_mass_pt_1D_r04->GetXaxis()->SetTitle("p_{T} [GeV/c]");
0381   jet_mass_pt_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0382   jet_mass_pt_1D_r04->SetMarkerStyle(8);
0383   jet_mass_pt_1D_r04->SetMarkerColor(1);
0384   jet_mass_pt_1D_r04->SetLineColor(1);
0385   jet_mass_pt_1D_r04->GetListOfFunctions()->Add(leg12);
0386   jet_mass_pt_1D_r04->Write("Average_Jet_Mass_pt_R04_pp");
0387 
0388 
0389   //for jet mass vs eta [R02]
0390   TLegend *leg13 = new TLegend(.7,.9,.9,1);
0391   leg13->SetFillStyle(0);
0392   leg13->SetBorderSize(0);
0393   leg13->SetTextSize(0.06);
0394   leg13->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0395   leg13->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0396   jet_mass_eta_r02->SetStats(0);
0397   jet_mass_eta_r02->SetTitle("Jet Mass vs #eta [R = 0.2]");
0398   jet_mass_eta_r02->GetListOfFunctions()->Add(leg13);
0399   jet_mass_eta_r02->Write("Jet_Mass_Eta_R02_pp");
0400 
0401   //for average jet mass vs eta [R02]
0402   TLegend *leg14 = new TLegend(.7,.9,.9,1);
0403   leg14->SetFillStyle(0);
0404   leg14->SetBorderSize(0);
0405   leg14->SetTextSize(0.06);
0406   leg14->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0407   leg14->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0408   jet_mass_eta_1D_r02 = (TH1D*)jet_mass_eta_r02->ProfileX();
0409   jet_mass_eta_1D_r02->SetStats(0);
0410   jet_mass_eta_1D_r02->SetTitle("Average Jet Mass vs #eta [R = 0.2]");
0411   jet_mass_eta_1D_r02->GetXaxis()->SetTitle("#eta");
0412   jet_mass_eta_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0413   jet_mass_eta_1D_r02->SetMarkerStyle(8);
0414   jet_mass_eta_1D_r02->SetMarkerColor(1);
0415   jet_mass_eta_1D_r02->SetLineColor(1);
0416   jet_mass_eta_1D_r02->GetListOfFunctions()->Add(leg14);
0417   jet_mass_eta_1D_r02->Write("Average_Jet_Mass_Eta_R02_pp");
0418 
0419   //for jet mass vs eta [R03]
0420   TLegend *leg15 = new TLegend(.7,.9,.9,1);
0421   leg15->SetFillStyle(0);
0422   leg15->SetBorderSize(0);
0423   leg15->SetTextSize(0.06);
0424   leg15->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0425   leg15->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0426   jet_mass_eta_r03->SetStats(0);
0427   jet_mass_eta_r03->SetTitle("Jet Mass vs #eta [R = 0.3]");
0428   jet_mass_eta_r03->GetListOfFunctions()->Add(leg15);
0429   jet_mass_eta_r03->Write("Jet_Mass_Eta_R03_pp");
0430 
0431   //for average jet mass vs eta [R03]
0432   TLegend *leg16 = new TLegend(.7,.9,.9,1);
0433   leg16->SetFillStyle(0);
0434   leg16->SetBorderSize(0);
0435   leg16->SetTextSize(0.06);
0436   leg16->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0437   leg16->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0438   jet_mass_eta_1D_r03 = (TH1D*)jet_mass_eta_r03->ProfileX();
0439   jet_mass_eta_1D_r03->SetStats(0);
0440   jet_mass_eta_1D_r03->SetTitle("Average Jet Mass vs #eta [R = 0.3]");
0441   jet_mass_eta_1D_r03->GetXaxis()->SetTitle("#eta");
0442   jet_mass_eta_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0443   jet_mass_eta_1D_r03->SetMarkerStyle(8);
0444   jet_mass_eta_1D_r03->SetMarkerColor(1);
0445   jet_mass_eta_1D_r03->SetLineColor(1);
0446   jet_mass_eta_1D_r03->GetListOfFunctions()->Add(leg16);
0447   jet_mass_eta_1D_r03->Write("Average_Jet_Mass_Eta_R03_pp");
0448 
0449 
0450   //for jet mass vs eta [R04]
0451   TLegend *leg17 = new TLegend(.7,.9,.9,1);
0452   leg17->SetFillStyle(0);
0453   leg17->SetBorderSize(0);
0454   leg17->SetTextSize(0.06);
0455   leg17->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0456   leg17->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0457   jet_mass_eta_r04->SetStats(0);
0458   jet_mass_eta_r04->SetTitle("Jet Mass vs #eta [R = 0.4]");
0459   jet_mass_eta_r04->GetListOfFunctions()->Add(leg17);
0460   jet_mass_eta_r04->Write("Jet_Mass_Eta_R04_pp");
0461 
0462   //for average jet mass vs eta [R04]
0463   TLegend *leg18 = new TLegend(.7,.9,.9,1);
0464   leg18->SetFillStyle(0);
0465   leg18->SetBorderSize(0);
0466   leg18->SetTextSize(0.06);
0467   leg18->AddEntry((TObject*)0, Form("%2.0f < p_{T} < %2.0f [GeV/c]", m_ptRange.first, m_ptRange.second),"");
0468   leg18->AddEntry((TObject*)0, Form("%1.1f < #eta < %1.1f", m_etaRange.first, m_etaRange.second),"");
0469   jet_mass_eta_1D_r04 = (TH1D*)jet_mass_eta_r04->ProfileX();
0470   jet_mass_eta_1D_r04->SetStats(0);
0471   jet_mass_eta_1D_r04->SetTitle("Average Jet Mass vs #eta [R = 0.4]");
0472   jet_mass_eta_1D_r04->GetXaxis()->SetTitle("#eta");
0473   jet_mass_eta_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV/c^{2}]");
0474   jet_mass_eta_1D_r04->SetMarkerStyle(8);
0475   jet_mass_eta_1D_r04->SetMarkerColor(1);
0476   jet_mass_eta_1D_r04->SetLineColor(1);
0477   jet_mass_eta_1D_r04->GetListOfFunctions()->Add(leg18);
0478   jet_mass_eta_1D_r04->Write("Average_Jet_Mass_Eta_R04_pp");
0479 
0480 
0481 
0482 
0483   std::cout << "JetKinematicCheck::End(PHCompositeNode *topNode) This is the End..." << std::endl;
0484   return Fun4AllReturnCodes::EVENT_OK;
0485 }
0486 
0487 //____________________________________________________________________________..
0488 int JetKinematicCheck::Reset(PHCompositeNode *topNode)
0489 {
0490  std::cout << "JetKinematicCheck::Reset(PHCompositeNode *topNode) being Reset" << std::endl;
0491   return Fun4AllReturnCodes::EVENT_OK;
0492 }
0493 
0494 //____________________________________________________________________________..
0495 void JetKinematicCheck::Print(const std::string &what) const
0496 {
0497   std::cout << "JetKinematicCheck::Print(const std::string &what) const Printing info for " << what << std::endl;
0498 }