File indexing completed on 2025-08-05 08:16:05
0001 #include "HeadReco.h"
0002
0003 #include <ffaobjects/EventHeader.h>
0004 #include <ffaobjects/EventHeaderv1.h>
0005 #include <ffaobjects/RunHeader.h>
0006 #include <ffaobjects/RunHeaderv1.h>
0007
0008 #include <phhepmc/PHHepMCGenEvent.h>
0009 #include <phhepmc/PHHepMCGenEventMap.h>
0010
0011 #include <fun4all/Fun4AllReturnCodes.h>
0012 #include <fun4all/Fun4AllServer.h>
0013 #include <fun4all/SubsysReco.h> // for SubsysReco
0014
0015 #include <phool/PHCompositeNode.h>
0016 #include <phool/PHIODataNode.h> // for PHIODataNode
0017 #include <phool/PHNode.h> // for PHNode
0018 #include <phool/PHNodeIterator.h> // for PHNodeIterator
0019 #include <phool/PHObject.h> // for PHObject
0020 #include <phool/getClass.h>
0021 #include <phool/recoConsts.h>
0022
0023 #include <Event/Event.h>
0024
0025 #include <HepMC/GenEvent.h>
0026 #include <HepMC/HeavyIon.h> // for HeavyIon
0027
0028 #include <iterator> // for operator!=, reverse_iterator
0029 #include <map> // for _Rb_tree_iterator
0030 #include <ranges>
0031 #include <utility> // for pair
0032
0033 HeadReco::HeadReco(const std::string &name)
0034 : SubsysReco(name)
0035 {
0036 }
0037
0038
0039
0040 int HeadReco::Init(PHCompositeNode *topNode)
0041 {
0042 PHNodeIterator iter(topNode);
0043 PHCompositeNode *runNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "RUN"));
0044 RunHeader *runheader = new RunHeaderv1();
0045 PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(runheader, "RunHeader", "PHObject");
0046 runNode->addNode(newNode);
0047
0048 PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
0049 EventHeader *eventheader = new EventHeaderv1();
0050 newNode = new PHIODataNode<PHObject>(eventheader, "EventHeader", "PHObject");
0051 dstNode->addNode(newNode);
0052
0053 return Fun4AllReturnCodes::EVENT_OK;
0054 }
0055
0056 int HeadReco::InitRun(PHCompositeNode *topNode)
0057 {
0058 recoConsts *rc = recoConsts::instance();
0059 RunHeader *runheader = findNode::getClass<RunHeader>(topNode, "RunHeader");
0060 runheader->set_RunNumber(rc->get_IntFlag("RUNNUMBER"));
0061 return Fun4AllReturnCodes::EVENT_OK;
0062 }
0063
0064 int HeadReco::process_event(PHCompositeNode *topNode)
0065 {
0066 Fun4AllServer *se = Fun4AllServer::instance();
0067 EventHeader *evtheader = findNode::getClass<EventHeader>(topNode, "EventHeader");
0068 PHHepMCGenEventMap *genevtmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
0069
0070 if (genevtmap)
0071 {
0072 for (auto &iter : std::ranges::reverse_view(*genevtmap))
0073 {
0074 PHHepMCGenEvent *genevt = iter.second;
0075 int embed_flag = genevt->get_embedding_id();
0076 if (embed_flag == 0)
0077 {
0078 HepMC::GenEvent *hepmcevt = genevt->getEvent();
0079
0080 if (hepmcevt)
0081 {
0082 HepMC::HeavyIon *hi = hepmcevt->heavy_ion();
0083 if (hi)
0084 {
0085 evtheader->set_ImpactParameter(hi->impact_parameter());
0086 evtheader->set_EventPlaneAngle(hi->event_plane_angle());
0087 evtheader->set_eccentricity(hi->eccentricity());
0088 evtheader->set_ncoll(hi->Ncoll());
0089 evtheader->set_npart(hi->Npart_targ() + hi->Npart_proj());
0090 }
0091 }
0092 }
0093 }
0094 }
0095 else
0096 {
0097 Event *evt = findNode::getClass<Event>(topNode, "PRDF");
0098 if (evt)
0099 {
0100 evtheader->set_EvtType(evt->getEvtType());
0101 }
0102 }
0103 evtheader->set_RunNumber(se->RunNumber());
0104 evtheader->set_EvtSequence(se->EventNumber());
0105 if (Verbosity() > 0)
0106 {
0107 evtheader->identify();
0108 }
0109
0110 return Fun4AllReturnCodes::EVENT_OK;
0111 }