Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:19

0001 //____________________________________________________________________________..
0002 
0003 #include "JetKinematicCheck.h"
0004 
0005 #include <calotrigger/TriggerAnalyzer.h>
0006 
0007 #include <fun4all/Fun4AllHistoManager.h>
0008 #include <fun4all/Fun4AllReturnCodes.h>
0009 
0010 #include <jetbase/JetContainer.h>
0011 #include <jetbase/Jetv2.h>
0012 
0013 #include <phool/PHCompositeNode.h>
0014 #include <phool/getClass.h>
0015 
0016 #include <TH1.h>
0017 #include <TH2.h>
0018 #include <TH3.h>
0019 #include <TLegend.h>
0020 #include <TPad.h>
0021 #include <TStyle.h>
0022 
0023 #include <algorithm>
0024 #include <cmath>
0025 #include <cstdlib>
0026 #include <string>
0027 #include <vector>
0028 //____________________________________________________________________________.. 
0029 JetKinematicCheck::JetKinematicCheck(const std::string &moduleName,
0030                                      const std::string &recojetnameR02,
0031                                      const std::string &recojetnameR03,
0032                                      const std::string &recojetnameR04,
0033                                      const std::string &recojetnameR05)
0034   : SubsysReco(moduleName)
0035   , m_moduleName(moduleName)
0036   , m_recoJetNameR02(recojetnameR02)
0037   , m_recoJetNameR03(recojetnameR03)
0038   , m_recoJetNameR04(recojetnameR04)
0039   , m_recoJetNameR05(recojetnameR05)
0040 {
0041   if (Verbosity() > 1)
0042   {
0043     std::cout << "JetKinematicCheck::JetKinematicCheck(const std::string &name) Calling ctor" << std::endl;
0044   }
0045 }
0046 
0047 //____________________________________________________________________________..
0048 JetKinematicCheck::~JetKinematicCheck()
0049 {
0050   if (Verbosity() > 1)
0051   {
0052     std::cout << "JetKinematicCheck::~JetKinematicCheck() Calling dtor" << std::endl;
0053   }
0054   delete m_analyzer;
0055 }
0056 
0057 //____________________________________________________________________________..
0058 int JetKinematicCheck::Init(PHCompositeNode * /*unused*/)
0059 {
0060   delete m_analyzer;
0061   m_analyzer = new TriggerAnalyzer();
0062   hm = QAHistManagerDef::getHistoManager();
0063   assert(hm);
0064 
0065   // make sure module name is lower case
0066   std::string smallModuleName = m_moduleName;
0067   std::transform(
0068       smallModuleName.begin(),
0069       smallModuleName.end(),
0070       smallModuleName.begin(),
0071       ::tolower);
0072 
0073   // construct histogram names
0074   std::vector<std::string> vecHistNames = {
0075       "r02_spectra",
0076       "r03_spectra",
0077       "r04_spectra",
0078       "r05_spectra",
0079       "r02_etavsphi",
0080       "r03_etavsphi",
0081       "r04_etavsphi",
0082       "r05_etavsphi",
0083       "r02_jetmassvspt",
0084       "r03_jetmassvspt",
0085       "r04_jetmassvspt",
0086       "r05_jetmassvspt",
0087       "r02_jetmassvseta",
0088       "r03_jetmassvseta",
0089       "r04_jetmassvseta",
0090       "r05_jetmassvseta",
0091       "r02_jetmassvsptprofile",
0092       "r03_jetmassvsptprofile",
0093       "r04_jetmassvsptprofile",
0094       "r05_jetmassvsptprofile",
0095       "r02_jetmassvsetaprofile",
0096       "r03_jetmassvsetaprofile",
0097       "r04_jetmassvsetaprofile",
0098       "r05_jetmassvsetaprofile",
0099       "r02_jetptvseta",
0100       "r02_jetptvsphi",
0101       "r03_jetptvseta",
0102       "r03_jetptvsphi",
0103       "r04_jetptvseta",
0104       "r04_jetptvsphi",
0105       "r05_jetptvseta",
0106       "r05_jetptvsphi"};
0107   
0108   for (auto &vecHistName : vecHistNames)
0109   {
0110     vecHistName.insert(0, "h_" + smallModuleName + "_");
0111     if (!m_histTag.empty())
0112     {
0113       vecHistName.append("_" + m_histTag);
0114     }
0115   }
0116 
0117   // initialize histograms
0118 
0119   h_jet_pt_spectra_r02 = new TH1D(vecHistNames[0].data(), "", 19, 5, 100);
0120   h_jet_pt_spectra_r02->GetXaxis()->SetTitle("p_{T} [GeV]");
0121 
0122   h_jet_pt_spectra_r03 = new TH1D(vecHistNames[1].data(), "", 19, 5, 100);
0123   h_jet_pt_spectra_r03->GetXaxis()->SetTitle("p_{T} [GeV]");
0124 
0125   h_jet_pt_spectra_r04 = new TH1D(vecHistNames[2].data(), "", 19, 5, 100);
0126   h_jet_pt_spectra_r04->GetXaxis()->SetTitle("p_{T} [GeV]");
0127 
0128   h_jet_pt_spectra_r05 = new TH1D(vecHistNames[3].data(), "", 19, 5, 100);
0129   h_jet_pt_spectra_r05->GetXaxis()->SetTitle("p_{T} [GeV]");
0130 
0131   h_jet_eta_phi_r02 = new TH2D(vecHistNames[4].data(), "", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0132   h_jet_eta_phi_r02->GetXaxis()->SetTitle("#eta");
0133   h_jet_eta_phi_r02->GetYaxis()->SetTitle("#Phi");
0134 
0135   h_jet_eta_phi_r03 = new TH2D(vecHistNames[5].data(), "", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0136   h_jet_eta_phi_r03->GetXaxis()->SetTitle("#eta");
0137   h_jet_eta_phi_r03->GetYaxis()->SetTitle("#Phi");
0138 
0139   h_jet_eta_phi_r04 = new TH2D(vecHistNames[6].data(), "", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0140   h_jet_eta_phi_r04->GetXaxis()->SetTitle("#eta");
0141   h_jet_eta_phi_r04->GetYaxis()->SetTitle("#Phi");
0142 
0143   h_jet_eta_phi_r05 = new TH2D(vecHistNames[7].data(), "", 24, -1.1, 1.1, 64, -M_PI, M_PI);
0144   h_jet_eta_phi_r05->GetXaxis()->SetTitle("#eta");
0145   h_jet_eta_phi_r05->GetYaxis()->SetTitle("#Phi");
0146 
0147   h_jet_mass_pt_r02 = new TH2D(vecHistNames[8].data(), "", 95, 5, 100, 72, -3, 15);
0148   h_jet_mass_pt_r02->GetXaxis()->SetTitle("p_{T} [GeV]");
0149   h_jet_mass_pt_r02->GetYaxis()->SetTitle("Jet Mass [GeV]");
0150 
0151   h_jet_mass_pt_r03 = new TH2D(vecHistNames[9].data(), "", 95, 5, 100, 72, -3, 15);
0152   h_jet_mass_pt_r03->GetXaxis()->SetTitle("p_{T} [GeV]");
0153   h_jet_mass_pt_r03->GetYaxis()->SetTitle("Jet Mass [GeV]");
0154 
0155   h_jet_mass_pt_r04 = new TH2D(vecHistNames[10].data(), "", 95, 5, 100, 72, -3, 15);
0156   h_jet_mass_pt_r04->GetXaxis()->SetTitle("p_{T} [GeV]");
0157   h_jet_mass_pt_r04->GetYaxis()->SetTitle("Jet Mass [GeV]");
0158 
0159   h_jet_mass_pt_r05 = new TH2D(vecHistNames[11].data(), "", 95, 5, 100, 72, -3, 15);
0160   h_jet_mass_pt_r05->GetXaxis()->SetTitle("p_{T} [GeV]");
0161   h_jet_mass_pt_r05->GetYaxis()->SetTitle("Jet Mass [GeV]");
0162 
0163   h_jet_mass_eta_r02 = new TH2D(vecHistNames[12].data(), "", 24, -1.1, 1.1, 72, -3, 15);
0164   h_jet_mass_eta_r02->GetXaxis()->SetTitle("#eta");
0165   h_jet_mass_eta_r02->GetYaxis()->SetTitle("Jet Mass [GeV]");
0166 
0167   h_jet_mass_eta_r03 = new TH2D(vecHistNames[13].data(), "", 24, -1.1, 1.1, 72, -3, 15);
0168   h_jet_mass_eta_r03->GetXaxis()->SetTitle("#eta");
0169   h_jet_mass_eta_r03->GetYaxis()->SetTitle("Jet Mass [GeV]");
0170 
0171   h_jet_mass_eta_r04 = new TH2D(vecHistNames[14].data(), "", 24, -1.1, 1.1, 72, -3, 15);
0172   h_jet_mass_eta_r04->GetXaxis()->SetTitle("#eta");
0173   h_jet_mass_eta_r04->GetYaxis()->SetTitle("Jet Mass [GeV]");
0174 
0175   h_jet_mass_eta_r05 = new TH2D(vecHistNames[15].data(), "", 24, -1.1, 1.1, 72, -3, 15);
0176   h_jet_mass_eta_r05->GetXaxis()->SetTitle("#eta");
0177   h_jet_mass_eta_r05->GetYaxis()->SetTitle("Jet Mass [GeV]");
0178 
0179   h_jet_average_mass_pt_1D_r02 = new TH1D(vecHistNames[16].data(), "", 19, 5, 100);
0180   h_jet_average_mass_pt_1D_r02->GetXaxis()->SetTitle("p_{T} [GeV]");
0181   h_jet_average_mass_pt_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0182 
0183   h_jet_average_mass_pt_1D_r03 = new TH1D(vecHistNames[17].data(), "", 19, 5, 100);
0184   h_jet_average_mass_pt_1D_r03->GetXaxis()->SetTitle("p_{T} [GeV]");
0185   h_jet_average_mass_pt_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0186 
0187   h_jet_average_mass_pt_1D_r04 = new TH1D(vecHistNames[18].data(), "", 19, 5, 100);
0188   h_jet_average_mass_pt_1D_r04->GetXaxis()->SetTitle("p_{T} [GeV]");
0189   h_jet_average_mass_pt_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0190 
0191   h_jet_average_mass_pt_1D_r05 = new TH1D(vecHistNames[19].data(), "", 19, 5, 100);
0192   h_jet_average_mass_pt_1D_r05->GetXaxis()->SetTitle("p_{T} [GeV]");
0193   h_jet_average_mass_pt_1D_r05->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0194 
0195   h_jet_average_mass_eta_1D_r02 = new TH2D(vecHistNames[20].data(), "", 24, -1.1, 1.1, 15, 0, 15);
0196   h_jet_average_mass_eta_1D_r02->GetXaxis()->SetTitle("#eta");
0197   h_jet_average_mass_eta_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0198 
0199   h_jet_average_mass_eta_1D_r03 = new TH2D(vecHistNames[21].data(), "", 24, -1.1, 1.1, 15, 0, 15);
0200   h_jet_average_mass_eta_1D_r03->GetXaxis()->SetTitle("#eta");
0201   h_jet_average_mass_eta_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0202 
0203   h_jet_average_mass_eta_1D_r04 = new TH2D(vecHistNames[22].data(), "", 24, -1.1, 1.1, 15, 0, 15);
0204   h_jet_average_mass_eta_1D_r04->GetXaxis()->SetTitle("#eta");
0205   h_jet_average_mass_eta_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0206 
0207   h_jet_average_mass_eta_1D_r05 = new TH2D(vecHistNames[23].data(), "", 24, -1.1, 1.1, 15, 0, 15);
0208   h_jet_average_mass_eta_1D_r05->GetXaxis()->SetTitle("#eta");
0209   h_jet_average_mass_eta_1D_r05->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0210 
0211   h_jet_pt_eta_r02 = new TH2D(vecHistNames[24].data(), "", 24, -1.1, 1.1, 19, 5, 100);
0212   h_jet_pt_eta_r02->GetXaxis()->SetTitle("eta");
0213   h_jet_pt_eta_r02->GetYaxis()->SetTitle("p_{T} [GeV]");
0214 
0215   h_jet_pt_phi_r02 = new TH2D(vecHistNames[25].data(), "", 64, -M_PI, M_PI, 19, 5, 100);
0216   h_jet_pt_phi_r02->GetXaxis()->SetTitle("phi");
0217   h_jet_pt_phi_r02->GetYaxis()->SetTitle("p_{T} [GeV]");
0218 
0219   h_jet_pt_eta_r03 = new TH2D(vecHistNames[26].data(), "", 24, -1.1, 1.1, 19, 5, 100);
0220   h_jet_pt_eta_r03->GetXaxis()->SetTitle("eta");
0221   h_jet_pt_eta_r03->GetYaxis()->SetTitle("p_{T} [GeV]");
0222 
0223   h_jet_pt_phi_r03 = new TH2D(vecHistNames[27].data(), "", 64, -M_PI, M_PI, 19, 5, 100);
0224   h_jet_pt_phi_r03->GetXaxis()->SetTitle("phi");
0225   h_jet_pt_phi_r03->GetYaxis()->SetTitle("p_{T} [GeV]");
0226 
0227   h_jet_pt_eta_r04 = new TH2D(vecHistNames[28].data(), "", 24, -1.1, 1.1, 19, 5, 100);
0228   h_jet_pt_eta_r04->GetXaxis()->SetTitle("eta");
0229   h_jet_pt_eta_r04->GetYaxis()->SetTitle("p_{T} [GeV]");
0230 
0231   h_jet_pt_phi_r04 = new TH2D(vecHistNames[29].data(), "", 64, -M_PI, M_PI, 19, 5, 100);
0232   h_jet_pt_phi_r04->GetXaxis()->SetTitle("phi");
0233   h_jet_pt_phi_r04->GetYaxis()->SetTitle("p_{T} [GeV]");
0234 
0235   h_jet_pt_eta_r05 = new TH2D(vecHistNames[30].data(), "", 24, -1.1, 1.1, 19, 5, 100);
0236   h_jet_pt_eta_r05->GetXaxis()->SetTitle("eta");
0237   h_jet_pt_eta_r05->GetYaxis()->SetTitle("p_{T} [GeV]");
0238 
0239   h_jet_pt_phi_r05 = new TH2D(vecHistNames[31].data(), "", 64, -M_PI, M_PI, 19, 5, 100);
0240   h_jet_pt_phi_r05->GetXaxis()->SetTitle("phi");
0241   h_jet_pt_phi_r05->GetYaxis()->SetTitle("p_{T} [GeV]");
0242   
0243   if (Verbosity() > 1)
0244   {
0245     std::cout << "JetKinematicCheck::Init(PHCompositeNode *topNode) Initializing" << std::endl;
0246   }
0247   return Fun4AllReturnCodes::EVENT_OK;
0248 }
0249 
0250 //____________________________________________________________________________..
0251 int JetKinematicCheck::InitRun(PHCompositeNode * /*unused*/)
0252 {
0253   if (Verbosity() > 1)
0254   {
0255     std::cout << "JetKinematicCheck::InitRun(PHCompositeNode *topNode) Initializing for Run XXX" << std::endl;
0256   }
0257   return Fun4AllReturnCodes::EVENT_OK;
0258 }
0259 
0260 //____________________________________________________________________________..
0261 int JetKinematicCheck::process_event(PHCompositeNode *topNode)
0262 {
0263   // if needed, check if selected trigger fired
0264   if (m_doTrgSelect)
0265   {
0266     m_analyzer->decodeTriggers(topNode);
0267     bool hasTrigger = JetQADefs::DidTriggerFire(m_trgToSelect, m_analyzer);
0268     if (!hasTrigger)
0269     {
0270       return Fun4AllReturnCodes::EVENT_OK;
0271     }
0272   }
0273 
0274   std::vector<std::string> m_recoJetName_array = {m_recoJetNameR02, m_recoJetNameR03, m_recoJetNameR04, m_recoJetNameR05};
0275   m_radii = {0.2, 0.3, 0.4, 0.5};
0276   int n_radii = m_radii.size();
0277 
0278   // Loop over each reco jet radii from array
0279   for (int i = 0; i < n_radii; i++)
0280   {
0281 
0282     // update eta range based on resolution parameter
0283     std::pair<float, float> etaRangeUse;
0284     if (m_restrictEtaRange)
0285     {
0286       etaRangeUse = {m_etaRange.first + m_radii[i], m_etaRange.second - m_radii[i]};
0287     }
0288     else
0289     {
0290       etaRangeUse = {m_etaRange.first, m_etaRange.second};
0291     }
0292 
0293     const std::string& recoJetName = m_recoJetName_array[i];
0294 
0295     JetContainer *jets = findNode::getClass<JetContainer>(topNode, recoJetName);
0296     if (!jets)
0297     {
0298       std::cout
0299           << "JetKinematicCheck::process_event - Error can not find DST Reco JetContainer node "
0300           << recoJetName << std::endl;
0301       return Fun4AllReturnCodes::EVENT_OK;
0302     }
0303 
0304     // loop over jets
0305     for (auto *jet : *jets)
0306     {
0307       bool eta_cut = (jet->get_eta() >= etaRangeUse.first) && (jet->get_eta() <= etaRangeUse.second);
0308       bool pt_cut = (jet->get_pt() >= m_ptRange.first) && (jet->get_pt() <= m_ptRange.second);
0309       if ((! eta_cut) || (! pt_cut))
0310       {
0311         continue;
0312       }
0313       if (jet->get_pt() < m_ptRange.first)
0314       {
0315         continue;  // to remove noise jets
0316       }
0317 
0318       if (i == 0)
0319       {
0320         h_jet_pt_spectra_r02->Fill(jet->get_pt());
0321         h_jet_eta_phi_r02->Fill(jet->get_eta(), jet->get_phi());
0322         h_jet_mass_pt_r02->Fill(jet->get_pt(), jet->get_mass());
0323         h_jet_mass_eta_r02->Fill(jet->get_eta(), jet->get_mass());
0324     h_jet_pt_eta_r02->Fill(jet->get_eta(), jet->get_pt());
0325     h_jet_pt_phi_r02->Fill(jet->get_phi(), jet->get_pt());
0326       }
0327 
0328       else if (i == 1)
0329       {
0330         h_jet_pt_spectra_r03->Fill(jet->get_pt());
0331         h_jet_eta_phi_r03->Fill(jet->get_eta(), jet->get_phi());
0332         h_jet_mass_pt_r03->Fill(jet->get_pt(), jet->get_mass());
0333         h_jet_mass_eta_r03->Fill(jet->get_eta(), jet->get_mass());
0334     h_jet_pt_eta_r03->Fill(jet->get_eta(), jet->get_pt());
0335     h_jet_pt_phi_r03->Fill(jet->get_phi(), jet->get_pt());
0336       }
0337 
0338       else if (i == 2)
0339       {
0340         h_jet_pt_spectra_r04->Fill(jet->get_pt());
0341         h_jet_eta_phi_r04->Fill(jet->get_eta(), jet->get_phi());
0342         h_jet_mass_pt_r04->Fill(jet->get_pt(), jet->get_mass());
0343         h_jet_mass_eta_r04->Fill(jet->get_eta(), jet->get_mass());
0344     h_jet_pt_eta_r04->Fill(jet->get_eta(), jet->get_pt());
0345     h_jet_pt_phi_r04->Fill(jet->get_phi(), jet->get_pt());
0346       }
0347 
0348       else if (i == 3)
0349       {
0350         h_jet_pt_spectra_r05->Fill(jet->get_pt());
0351         h_jet_eta_phi_r05->Fill(jet->get_eta(), jet->get_phi());
0352         h_jet_mass_pt_r05->Fill(jet->get_pt(), jet->get_mass());
0353         h_jet_mass_eta_r05->Fill(jet->get_eta(), jet->get_mass());
0354     h_jet_pt_eta_r05->Fill(jet->get_eta(), jet->get_pt());
0355     h_jet_pt_phi_r05->Fill(jet->get_phi(), jet->get_pt());
0356       }
0357     }
0358   }
0359 
0360   if (Verbosity() > 1)
0361   {
0362     std::cout << "JetKinematicCheck::process_event(PHCompositeNode *topNode) Processing Event" << std::endl;
0363   }
0364   return Fun4AllReturnCodes::EVENT_OK;
0365 }
0366 
0367 //____________________________________________________________________________..
0368 int JetKinematicCheck::EndRun(const int runnumber)
0369 {
0370   if (Verbosity() > 1)
0371   {
0372     std::cout << "JetKinematicCheck::EndRun(const int runnumber) Ending Run for Run " << runnumber << std::endl;
0373   }
0374   return Fun4AllReturnCodes::EVENT_OK;
0375 }
0376 
0377 //____________________________________________________________________________..
0378 int JetKinematicCheck::End(PHCompositeNode * /*unused*/)
0379 {
0380   if (Verbosity() > 1)
0381   {
0382     std::cout << "JetKinematicCheck::End(PHCompositeNode *topNode) Entering the end" << std::endl;
0383   }
0384 
0385   // for jet spectra [R02]
0386   h_jet_pt_spectra_r02->SetMarkerStyle(8);
0387   h_jet_pt_spectra_r02->SetMarkerColor(1);
0388   h_jet_pt_spectra_r02->SetLineColor(1);
0389   h_jet_pt_spectra_r02->SetTitle("Jet p_{T} Spectra for R02");
0390   h_jet_pt_spectra_r02->GetYaxis()->SetTitle("Counts");
0391   h_jet_pt_spectra_r02->SetStats(false);
0392 
0393   // for jet spectra [R03]
0394   h_jet_pt_spectra_r03->SetMarkerStyle(8);
0395   h_jet_pt_spectra_r03->SetMarkerColor(1);
0396   h_jet_pt_spectra_r03->SetLineColor(1);
0397   h_jet_pt_spectra_r03->SetTitle("Jet p_{T} Spectra for R03");
0398   h_jet_pt_spectra_r03->GetYaxis()->SetTitle("Counts");
0399   h_jet_pt_spectra_r03->SetStats(false);
0400 
0401   // for jet spectra [R04]
0402   h_jet_pt_spectra_r04->SetMarkerStyle(8);
0403   h_jet_pt_spectra_r04->SetMarkerColor(1);
0404   h_jet_pt_spectra_r04->SetLineColor(1);
0405   h_jet_pt_spectra_r04->SetTitle("Jet p_{T} Spectra for R04");
0406   h_jet_pt_spectra_r04->GetYaxis()->SetTitle("Counts");
0407   h_jet_pt_spectra_r04->SetStats(false);
0408 
0409   // for jet spectra [R05]
0410   h_jet_pt_spectra_r05->SetMarkerStyle(8);
0411   h_jet_pt_spectra_r05->SetMarkerColor(1);
0412   h_jet_pt_spectra_r05->SetLineColor(1);
0413   h_jet_pt_spectra_r05->SetTitle("Jet p_{T} Spectra for R05");
0414   h_jet_pt_spectra_r05->GetYaxis()->SetTitle("Counts");
0415   h_jet_pt_spectra_r05->SetStats(false);
0416 
0417   // for jet eta-phi [R02]
0418   h_jet_eta_phi_r02->SetStats(false);
0419   h_jet_eta_phi_r02->SetTitle("Jet #eta vs #phi for R02");
0420 
0421   // for jet eta-phi [R03]
0422   h_jet_eta_phi_r03->SetStats(false);
0423   h_jet_eta_phi_r03->SetTitle("Jet #eta vs #phi for R03");
0424 
0425   // for jet eta-phi [R04]
0426   h_jet_eta_phi_r04->SetStats(false);
0427   h_jet_eta_phi_r04->SetTitle("Jet #eta vs #phi for R04");
0428 
0429   // for jet eta-phi [R05]
0430   h_jet_eta_phi_r05->SetStats(false);
0431   h_jet_eta_phi_r05->SetTitle("Jet #eta vs #phi for R05");
0432 
0433   // for jet mass vs pt [R02]
0434   h_jet_mass_pt_r02->SetStats(false);
0435   h_jet_mass_pt_r02->SetTitle("Jet Mass vs p_{T} for R02");
0436 
0437   // jet pt vs eta
0438   h_jet_pt_eta_r02->SetStats(false);
0439   h_jet_pt_eta_r02->SetTitle("Jet p_{T} vs #eta for R02");
0440 
0441   h_jet_pt_eta_r03->SetStats(false);
0442   h_jet_pt_eta_r03->SetTitle("Jet p_{T} vs #eta for R03");
0443 
0444   h_jet_pt_eta_r04->SetStats(false);
0445   h_jet_pt_eta_r04->SetTitle("Jet p_{T} vs #eta for R04");
0446 
0447   h_jet_pt_eta_r05->SetStats(false);
0448   h_jet_pt_eta_r05->SetTitle("Jet p_{T} vs #eta for R05");
0449 
0450   // jet pt vs phi
0451   h_jet_pt_phi_r02->SetStats(false);
0452   h_jet_pt_phi_r02->SetTitle("Jet p_{T} vs #phi for R02");
0453 
0454   h_jet_pt_phi_r03->SetStats(false);
0455   h_jet_pt_phi_r03->SetTitle("Jet p_{T} vs #phi for R03");
0456 
0457   h_jet_pt_phi_r04->SetStats(false);
0458   h_jet_pt_phi_r04->SetTitle("Jet p_{T} vs #phi for R04");
0459 
0460   h_jet_pt_phi_r05->SetStats(false);
0461   h_jet_pt_phi_r05->SetTitle("Jet p_{T} vs #phi for R05");
0462   
0463   // for average jet mass vs pt [R02]
0464   h_jet_average_mass_pt_1D_r02 = (TH1D *) h_jet_mass_pt_r02->ProfileX();
0465   h_jet_average_mass_pt_1D_r02->SetStats(false);
0466   h_jet_average_mass_pt_1D_r02->SetTitle("Average Jet Mass vs p_{T} [GeV] R02");
0467   h_jet_average_mass_pt_1D_r02->GetXaxis()->SetTitle("p_{T} [GeV]");
0468   h_jet_average_mass_pt_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0469   h_jet_average_mass_pt_1D_r02->SetMarkerStyle(8);
0470   h_jet_average_mass_pt_1D_r02->SetMarkerColor(1);
0471   h_jet_average_mass_pt_1D_r02->SetLineColor(1);
0472 
0473   // for jet mass vs pt [R03]
0474   h_jet_mass_pt_r03->SetStats(false);
0475   h_jet_mass_pt_r03->SetTitle("Jet Mass vs p_{T}[GeV] R03");
0476 
0477   // for average jet mass vs pt [R03]
0478   h_jet_average_mass_pt_1D_r03 = (TH1D *) h_jet_mass_pt_r03->ProfileX();
0479   h_jet_average_mass_pt_1D_r03->SetStats(false);
0480   h_jet_average_mass_pt_1D_r03->SetTitle("Average Jet Mass vs p_{T} [GeV] R03");
0481   h_jet_average_mass_pt_1D_r03->GetXaxis()->SetTitle("p_{T} [GeV]");
0482   h_jet_average_mass_pt_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0483   h_jet_average_mass_pt_1D_r03->SetMarkerStyle(8);
0484   h_jet_average_mass_pt_1D_r03->SetMarkerColor(1);
0485   h_jet_average_mass_pt_1D_r03->SetLineColor(1);
0486 
0487   // for jet mass vs pt [R04]
0488   h_jet_mass_pt_r04->SetStats(false);
0489   h_jet_mass_pt_r04->SetTitle("Jet Mass vs p_{T}  [GeV] R04");
0490 
0491   // for average jet mass vs pt [R04]
0492   h_jet_average_mass_pt_1D_r04 = (TH1D *) h_jet_mass_pt_r04->ProfileX();
0493   h_jet_average_mass_pt_1D_r04->SetStats(false);
0494   h_jet_average_mass_pt_1D_r04->SetTitle("Average Jet Mass vs p_{T} [GeV] R04");
0495   h_jet_average_mass_pt_1D_r04->GetXaxis()->SetTitle("p_{T} [GeV]");
0496   h_jet_average_mass_pt_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0497   h_jet_average_mass_pt_1D_r04->SetMarkerStyle(8);
0498   h_jet_average_mass_pt_1D_r04->SetMarkerColor(1);
0499   h_jet_average_mass_pt_1D_r04->SetLineColor(1);
0500 
0501   // for jet mass vs pt [R05]
0502   h_jet_mass_pt_r05->SetStats(false);
0503   h_jet_mass_pt_r05->SetTitle("Jet Mass vs p_{T} [GeV] R05");
0504 
0505   // for average jet mass vs pt [R05]
0506   h_jet_average_mass_pt_1D_r05 = (TH1D *) h_jet_mass_pt_r05->ProfileX();
0507   h_jet_average_mass_pt_1D_r05->SetStats(false);
0508   h_jet_average_mass_pt_1D_r05->SetTitle("Average Jet Mass vs p_{T} [GeV] R05");
0509   h_jet_average_mass_pt_1D_r05->GetXaxis()->SetTitle("p_{T} [GeV]");
0510   h_jet_average_mass_pt_1D_r05->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0511   h_jet_average_mass_pt_1D_r05->SetMarkerStyle(8);
0512   h_jet_average_mass_pt_1D_r05->SetMarkerColor(1);
0513   h_jet_average_mass_pt_1D_r05->SetLineColor(1);
0514 
0515   // for jet mass vs eta [R02]
0516   h_jet_mass_eta_r02->SetStats(false);
0517   h_jet_mass_eta_r02->SetTitle("Jet Mass vs p_{T} [GeV] R02");
0518 
0519   // for average jet mass vs eta [R02]
0520   h_jet_average_mass_eta_1D_r02 = (TH1D *) h_jet_mass_eta_r02->ProfileX();
0521   h_jet_average_mass_eta_1D_r02->SetStats(false);
0522   h_jet_average_mass_eta_1D_r02->SetTitle("Average Jet Mass vs #eta R02");
0523   h_jet_average_mass_eta_1D_r02->GetXaxis()->SetTitle("#eta");
0524   h_jet_average_mass_eta_1D_r02->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0525   h_jet_average_mass_eta_1D_r02->SetMarkerStyle(8);
0526   h_jet_average_mass_eta_1D_r02->SetMarkerColor(1);
0527   h_jet_average_mass_eta_1D_r02->SetLineColor(1);
0528 
0529   // for jet mass vs eta [R03]
0530   h_jet_mass_eta_r03->SetStats(false);
0531   h_jet_mass_eta_r03->SetTitle("Jet Mass vs #eta [GeV] R03");
0532 
0533   // for average jet mass vs eta [R03]
0534   h_jet_average_mass_eta_1D_r03 = (TH1D *) h_jet_mass_eta_r03->ProfileX();
0535   h_jet_average_mass_eta_1D_r03->SetStats(false);
0536   h_jet_average_mass_eta_1D_r03->SetTitle("Average Jet Mass vs #eta R03");
0537   h_jet_average_mass_eta_1D_r03->GetXaxis()->SetTitle("#eta");
0538   h_jet_average_mass_eta_1D_r03->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0539   h_jet_average_mass_eta_1D_r03->SetMarkerStyle(8);
0540   h_jet_average_mass_eta_1D_r03->SetMarkerColor(1);
0541   h_jet_average_mass_eta_1D_r03->SetLineColor(1);
0542 
0543   // for jet mass vs eta [R04]
0544   h_jet_mass_eta_r04->SetStats(false);
0545   h_jet_mass_eta_r04->SetTitle("Jet Mass vs #eta [GeV] R04");
0546 
0547   // for average jet mass vs eta [R04]
0548   h_jet_average_mass_eta_1D_r04 = (TH1D *) h_jet_mass_eta_r04->ProfileX();
0549   h_jet_average_mass_eta_1D_r04->SetStats(false);
0550   h_jet_average_mass_eta_1D_r04->SetTitle("Average Jet Mass vs #eta R04");
0551   h_jet_average_mass_eta_1D_r04->GetXaxis()->SetTitle("#eta");
0552   h_jet_average_mass_eta_1D_r04->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0553   h_jet_average_mass_eta_1D_r04->SetMarkerStyle(8);
0554   h_jet_average_mass_eta_1D_r04->SetMarkerColor(1);
0555   h_jet_average_mass_eta_1D_r04->SetLineColor(1);
0556 
0557   // for jet mass vs eta [R05]
0558   h_jet_mass_eta_r05->SetStats(false);
0559   h_jet_mass_eta_r05->SetTitle("Jet Mass vs #eta [GeV] R05");
0560 
0561   // for average jet mass vs eta [R05]
0562   h_jet_average_mass_eta_1D_r05 = (TH1D *) h_jet_mass_eta_r05->ProfileX();
0563   h_jet_average_mass_eta_1D_r05->SetStats(false);
0564   h_jet_average_mass_eta_1D_r05->SetTitle("Average Jet Mass vs #eta R05");
0565   h_jet_average_mass_eta_1D_r05->GetXaxis()->SetTitle("#eta");
0566   h_jet_average_mass_eta_1D_r05->GetYaxis()->SetTitle("Average Jet Mass [GeV]");
0567   h_jet_average_mass_eta_1D_r05->SetMarkerStyle(8);
0568   h_jet_average_mass_eta_1D_r05->SetMarkerColor(1);
0569   h_jet_average_mass_eta_1D_r05->SetLineColor(1);
0570 
0571 
0572   
0573   hm->registerHisto(h_jet_pt_spectra_r02);
0574   hm->registerHisto(h_jet_pt_spectra_r03);
0575   hm->registerHisto(h_jet_pt_spectra_r04);
0576   hm->registerHisto(h_jet_pt_spectra_r05);
0577   hm->registerHisto(h_jet_eta_phi_r02);
0578   hm->registerHisto(h_jet_eta_phi_r03);
0579   hm->registerHisto(h_jet_eta_phi_r04);
0580   hm->registerHisto(h_jet_eta_phi_r05);
0581   hm->registerHisto(h_jet_mass_pt_r02);
0582   hm->registerHisto(h_jet_mass_pt_r03);
0583   hm->registerHisto(h_jet_mass_pt_r04);
0584   hm->registerHisto(h_jet_mass_pt_r05);
0585   hm->registerHisto(h_jet_pt_eta_r02);
0586   hm->registerHisto(h_jet_pt_eta_r03);
0587   hm->registerHisto(h_jet_pt_eta_r04);
0588   hm->registerHisto(h_jet_pt_eta_r05);
0589   hm->registerHisto(h_jet_pt_phi_r02);
0590   hm->registerHisto(h_jet_pt_phi_r03);
0591   hm->registerHisto(h_jet_pt_phi_r04);
0592   hm->registerHisto(h_jet_pt_phi_r05);
0593   hm->registerHisto(h_jet_average_mass_pt_1D_r02);
0594   hm->registerHisto(h_jet_average_mass_pt_1D_r03);
0595   hm->registerHisto(h_jet_average_mass_pt_1D_r04);
0596   hm->registerHisto(h_jet_average_mass_pt_1D_r05);
0597   hm->registerHisto(h_jet_mass_eta_r02);
0598   hm->registerHisto(h_jet_mass_eta_r03);
0599   hm->registerHisto(h_jet_mass_eta_r04);
0600   hm->registerHisto(h_jet_mass_eta_r05);
0601   hm->registerHisto(h_jet_average_mass_eta_1D_r02);
0602   hm->registerHisto(h_jet_average_mass_eta_1D_r03);
0603   hm->registerHisto(h_jet_average_mass_eta_1D_r04);
0604   hm->registerHisto(h_jet_average_mass_eta_1D_r05);
0605 
0606   if (Verbosity() > 1)
0607   {
0608     std::cout << "JetKinematicCheck::End(PHCompositeNode *topNode) This is the End..." << std::endl;
0609   }
0610   return Fun4AllReturnCodes::EVENT_OK;
0611 }