hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
StdhepMCParticleProcessor.cxx
Go to the documentation of this file.
1
8#include "utilities.h"
9
11 : Processor(name, process) {
12 }
13
16
18
19 std::cout << "Configuring StdhepMCParticleProcessor" << std::endl;
20 try{
21 mcPartCollStdhep_ = parameters.getString("mcPartCollStdhep", mcPartCollStdhep_);
22 mcPartCollRoot_ = parameters.getString("mcPartCollRoot", mcPartCollRoot_);
23 maxEvent_ = parameters.getInteger("maxEvent",maxEvent_);
24 skipEvent_ = parameters.getInteger("skipEvent",skipEvent_);
25 }
26 catch (std::runtime_error& error)
27 {
28 std::cout << error.what() << std::endl;
29 }
30}
31
32void StdhepMCParticleProcessor::initialize(std::string inFilename, std::string outFilename) {
33 // Init Files
34 std::cout << "[StdhepMCParticleProcessor] initialize" << std::endl;
35 inFilename_ = inFilename;
36 outF_ = new TFile(outFilename.c_str(),"RECREATE");
37
38 std::cout << "Convert " << inFilename_ << " to ROOT " << outF_->GetName() << std::endl;
39
40 //TTree to store stdhep as root
41 tree_ = new TTree("HPS_Event","conversion");
42 tree_->Branch(mcPartCollRoot_.c_str(),&mc_particles_);
43}
44
46 std::cout << "[StdhepMCParticleProcessor] Starting process()" << std::endl;
47
48 std::cout << "opening file : " << inFilename_ << std::endl;
49 UTIL::LCStdHepRdr rdr(inFilename_.c_str());
50 rdr.printHeader();
51 std::stringstream description;
52
53 description << " file generated with LCIO stdhepjob from " << inFilename_;
54
55 int count = skipEvent_;
56
57 try {
58
59 while( maxEvent_ < 0 || count < maxEvent_ ){
60
61 if (mc_particles_.size() > 0){
62 for (std::vector<MCParticle*>::iterator it = mc_particles_.begin(); it != mc_particles_.end(); ++it){
63 delete *it;
64 }
65 mc_particles_.clear();
66 }
67
68 std::unique_ptr<IMPL::LCEventImpl> evt( new IMPL::LCEventImpl() ) ;
69 evt->setRunNumber(0) ;
70 evt->setEventNumber(count) ;
71
72 // read the next stdhep event and add an MCParticle collection to the event
73 rdr.updateNextEvent(evt.get(),mcPartCollStdhep_.c_str()) ;
74
75 //Get collection from event
76 EVENT::LCCollection* lc_particles{nullptr};
77 try
78 {
79 lc_particles = evt->getCollection(mcPartCollStdhep_.c_str());
80 //Loop through all particles in event
81 for (int iparticle = 0; iparticle < lc_particles->getNumberOfElements(); ++iparticle){
82 //Get particle from LCEvent
83 IMPL::MCParticleImpl* lc_particle = static_cast<IMPL::MCParticleImpl*>(lc_particles->getElementAt(iparticle));
84
85 //Make MCParticle to build
86 MCParticle* particle = new MCParticle();
87
88 //Set charge of HpsMCParticle
89 particle->setCharge(lc_particle->getCharge());
90
91 // Set the HpsMCParticle type
92 particle->setTime(lc_particle->getTime());
93
94 // Set the energy of the HpsMCParticle
95 particle->setEnergy(lc_particle->getEnergy());
96
97 // Set the momentum of the HpsMCParticle
98 particle->setMomentum(lc_particle->getMomentum());
99
100 // Set the momentum of HpsMCParticle at Endpoint
101 particle->setEndpointMomentum(lc_particle->getMomentumAtEndpoint());
102
103 // Set the mass of the HpsMCParticle
104 particle->setMass(lc_particle->getMass());
105
106 // Set the PDG of the particle
107 particle->setPDG(lc_particle->getPDG());
108
109 // Set the LCIO id of the particle
110 particle->setID(lc_particle->id());
111
112 // Set the PDG of the particle
113 std::vector<EVENT::MCParticle*> parentVec = lc_particle->getParents();
114 if(parentVec.size() > 0) particle->setMomPDG(parentVec.at(0)->getPDG());
115
116 // Set the generator status of the particle
117 particle->setGenStatus(lc_particle->getGeneratorStatus());
118
119 // Set the generator status of the particle
120 particle->setSimStatus(lc_particle->getSimulatorStatus());
121
122 // Set the vertex position of the particle
123 particle->setVertexPosition(lc_particle->getVertex());
124
125 // Set the vertex position of the particle
126 particle->setEndPoint(lc_particle->getEndpoint());
127
128 mc_particles_.push_back(particle);
129 }
130
131 tree_->Fill();
132 }
133 catch (EVENT::DataNotAvailableException e)
134 {
135 std::cout << e.what() << std::endl;
136 }
137
138 ++count ;
139
140 } // evt loop
141 }
142
143 catch( IO::EndOfDataException& e ) {
144 }
145
146 std::cout << " converted " << count << std::endl ;
147
148 std::cout << "==================================================== "
149 << std::endl << std::endl ;
150
151 return true;
152}
153
155 outF_->cd();
156 outF_->Write();
157 outF_->Close();
158}
159
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
Processor used to translate StdHep MCParticles to ROOT MCParticle objects.
void setVertexPosition(const double *vtx_pos)
void setEndpointMomentum(const double *momentum_ep)
void setPDG(const int pdg)
Definition MCParticle.h:60
void setMass(const double mass)
Definition MCParticle.h:102
void setID(const int id)
Definition MCParticle.h:67
void setGenStatus(const int gen)
Definition MCParticle.h:81
void setMomPDG(const int momPDG)
Definition MCParticle.h:74
void setCharge(const int charge)
Definition MCParticle.h:53
void setMomentum(const double *momentum)
void setEnergy(const double energy)
Definition MCParticle.h:95
void setTime(const double time)
Definition MCParticle.h:109
void setEndPoint(const double *ep_pos)
void setSimStatus(const int sim)
Definition MCParticle.h:88
description
Base class for all event processing components.
Definition Processor.h:34
Processor used to translate StdHep MCParticles to ROOT MCParticle objects. more details.
std::string inFilename_
stdhep input file
TFile * outF_
root tuple outfile
std::string mcPartCollStdhep_
name temporary lcio collection
virtual void configure(const ParameterSet &parameters)
Configure the Processor.
virtual void finalize()
Callback for the Processor to take any necessary action when the processing of events finishes.
int skipEvent_
skipped event numbers to convet
TTree * tree_
TTree holds converted MCParticles for each event.
int maxEvent_
max stdhep event number to convert
std::vector< MCParticle * > mc_particles_
list of converted MCParticles
virtual void initialize(std::string inFilename, std::string outFilename)
Callback for the Processor to take any necessary action when the processing of events starts.
StdhepMCParticleProcessor(const std::string &name, Process &process)
Class constructor.
std::string mcPartCollRoot_
name root collection
virtual bool process()
Process the event and put new data products into it.