Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:43

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