Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:17:45

0001 #include "ClusterJetInput.h"
0002 
0003 #include "Jet.h"
0004 #include "Jetv2.h"
0005 
0006 #include <calobase/RawCluster.h>
0007 #include <calobase/RawClusterContainer.h>
0008 #include <calobase/RawClusterUtility.h>
0009 #include <globalvertex/GlobalVertex.h>
0010 #include <globalvertex/GlobalVertexMap.h>
0011 #include <phool/getClass.h>
0012 
0013 #include <CLHEP/Vector/ThreeVector.h>  // for Hep3Vector
0014 
0015 // standard includes
0016 #include <cassert>
0017 #include <iostream>
0018 #include <map>      // for _Rb_tree_const_iterator
0019 #include <utility>  // for pair
0020 #include <vector>
0021 
0022 ClusterJetInput::ClusterJetInput(Jet::SRC input)
0023   : m_Input(input)
0024 {
0025 }
0026 
0027 void ClusterJetInput::identify(std::ostream &os)
0028 {
0029   os << "   ClusterJetInput: ";
0030   if (m_Input == Jet::CEMC_CLUSTER)
0031   {
0032     os << "CLUSTER_CEMC to Jet::CEMC_CLUSTER";
0033   }
0034   else if (m_Input == Jet::EEMC_CLUSTER)
0035   {
0036     os << "CLUSTER_EEMC to Jet::EEMC_CLUSTER";
0037   }
0038   else if (m_Input == Jet::HCALIN_CLUSTER)
0039   {
0040     os << "CLUSTER_HCALIN to Jet::HCALIN_CLUSTER";
0041   }
0042   else if (m_Input == Jet::HCALOUT_CLUSTER)
0043   {
0044     os << "CLUSTER_HCALOUT to Jet::HCALOUT_CLUSTER";
0045   }
0046   os << std::endl;
0047 }
0048 
0049 std::vector<Jet *> ClusterJetInput::get_input(PHCompositeNode *topNode)
0050 {
0051   if (m_Verbosity > 0)
0052   {
0053     std::cout << "ClusterJetInput::process_event -- entered" << std::endl;
0054   }
0055 
0056   CLHEP::Hep3Vector vertex(0, 0, 0);
0057   GlobalVertexMap *vertexmap = findNode::getClass<GlobalVertexMap>(topNode, "GlobalVertexMap");
0058   if (!vertexmap)
0059   {
0060     std::cout << "ClusterJetInput::get_input - Fatal Error - GlobalVertexMap node is missing. Please turn on the do_global flag in the main macro in order to reconstruct the global vertex." << std::endl;
0061     assert(vertexmap);  // force quit
0062 
0063     return std::vector<Jet *>();
0064   }
0065 
0066   if (vertexmap->empty())
0067   {
0068     std::cout << "ClusterJetInput::get_input - Warning - GlobalVertexMap node is empty. Continue as if vtxz = (0,0,0)." << std::endl;
0069   }
0070   else
0071   {
0072     GlobalVertex *vtx = vertexmap->begin()->second;
0073     if (vtx)
0074     {
0075       if (m_use_vertextype)
0076       {
0077         auto typeStartIter = vtx->find_vertexes(m_vertex_type);
0078         auto typeEndIter = vtx->end_vertexes();
0079         for (auto iter = typeStartIter; iter != typeEndIter; ++iter)
0080         {
0081           const auto &[type, vertexVec] = *iter;
0082           if (type != m_vertex_type)
0083           {
0084             continue;
0085           }
0086           for (const auto *v : vertexVec)
0087           {
0088             if (!v)
0089             {
0090               continue;
0091             }
0092             vertex.set(v->get_x(), v->get_y(), v->get_z());
0093           }
0094         }
0095       }
0096       else
0097       {
0098         vertex.set(vtx->get_x(), vtx->get_y(), vtx->get_z());
0099       }
0100     }
0101 
0102     if (std::isnan(vertex.z()))
0103     {
0104       vertex.set(0, 0, 0);
0105     }
0106   }
0107 
0108   RawClusterContainer *clusters = nullptr;
0109   if (m_Input == Jet::CEMC_CLUSTER)
0110   {
0111     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_CEMC");
0112     if (!clusters)
0113     {
0114       return std::vector<Jet *>();
0115     }
0116   }
0117   else if (m_Input == Jet::EEMC_CLUSTER)
0118   {
0119     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_EEMC");
0120     if (!clusters)
0121     {
0122       return std::vector<Jet *>();
0123     }
0124   }
0125   else if (m_Input == Jet::HCALIN_CLUSTER)
0126   {
0127     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_HCALIN");
0128     if (!clusters)
0129     {
0130       return std::vector<Jet *>();
0131     }
0132   }
0133   else if (m_Input == Jet::HCALOUT_CLUSTER)
0134   {
0135     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_HCALOUT");
0136     if (!clusters)
0137     {
0138       return std::vector<Jet *>();
0139     }
0140   }
0141   else if (m_Input == Jet::HCAL_TOPO_CLUSTER)
0142   {
0143     clusters = findNode::getClass<RawClusterContainer>(topNode, "TOPOCLUSTER_HCAL");
0144     if (!clusters)
0145     {
0146       return std::vector<Jet *>();
0147     }
0148   }
0149   else if (m_Input == Jet::ECAL_TOPO_CLUSTER)
0150   {
0151     clusters = findNode::getClass<RawClusterContainer>(topNode, "TOPOCLUSTER_EMCAL");
0152     if (!clusters)
0153     {
0154       return std::vector<Jet *>();
0155     }
0156   }
0157   else if (m_Input == Jet::ECAL_HCAL_TOPO_CLUSTER)
0158   {
0159     clusters = findNode::getClass<RawClusterContainer>(topNode, "TOPOCLUSTER_ALLCALO");
0160     if (!clusters)
0161     {
0162       return std::vector<Jet *>();
0163     }
0164   }
0165   else if (m_Input == Jet::FEMC_CLUSTER)
0166   {
0167     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_FEMC");
0168     if (!clusters)
0169     {
0170       return std::vector<Jet *>();
0171     }
0172   }
0173   else if (m_Input == Jet::FHCAL_CLUSTER)
0174   {
0175     clusters = findNode::getClass<RawClusterContainer>(topNode, "CLUSTER_FHCAL");
0176     if (!clusters)
0177     {
0178       return std::vector<Jet *>();
0179     }
0180   }
0181   else
0182   {
0183     return std::vector<Jet *>();
0184   }
0185 
0186   std::vector<Jet *> pseudojets;
0187   RawClusterContainer::ConstRange begin_end = clusters->getClusters();
0188   RawClusterContainer::ConstIterator rtiter;
0189   for (rtiter = begin_end.first; rtiter != begin_end.second; ++rtiter)
0190   {
0191     RawCluster *cluster = rtiter->second;
0192 
0193     CLHEP::Hep3Vector E_vec_cluster = RawClusterUtility::GetEVec(*cluster, vertex);
0194 
0195     Jet *jet = new Jetv2();
0196     jet->set_px(E_vec_cluster.x());
0197     jet->set_py(E_vec_cluster.y());
0198     jet->set_pz(E_vec_cluster.z());
0199     jet->set_e(cluster->get_energy());
0200     jet->insert_comp(m_Input, cluster->get_id());
0201     pseudojets.push_back(jet);
0202   }
0203 
0204   if (m_Verbosity > 0)
0205   {
0206     std::cout << "ClusterJetInput::process_event -- exited" << std::endl;
0207   }
0208 
0209   return pseudojets;
0210 }