CMS 3D CMS Logo

RHStopDump.cc
Go to the documentation of this file.
2 
6 
8  : mStream(parameters.getParameter<std::string>("stoppedFile").c_str()),
9  mProducer(parameters.getUntrackedParameter<std::string>("producer", "g4SimHits")),
10  tokNames_(consumes<std::vector<std::string> >(edm::InputTag(mProducer, "StoppedParticlesName"))),
11  tokenXs_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesX"))),
12  tokenYs_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesY"))),
13  tokenZs_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesZ"))),
14  tokenTs_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesTime"))),
15  tokenIds_(consumes<std::vector<int> >(edm::InputTag(mProducer, "StoppedParticlesPdgId"))),
16  tokenMasses_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesMass"))),
17  tokenCharges_(consumes<std::vector<float> >(edm::InputTag(mProducer, "StoppedParticlesCharge"))) {}
18 
21  const edm::Handle<std::vector<float> >& xs = fEvent.getHandle(tokenXs_);
22  const edm::Handle<std::vector<float> >& ys = fEvent.getHandle(tokenYs_);
23  const edm::Handle<std::vector<float> >& zs = fEvent.getHandle(tokenZs_);
24  const edm::Handle<std::vector<float> >& ts = fEvent.getHandle(tokenTs_);
25  const edm::Handle<std::vector<int> >& ids = fEvent.getHandle(tokenIds_);
26  const edm::Handle<std::vector<float> >& masses = fEvent.getHandle(tokenMasses_);
28 
29  if (names->size() != xs->size() || xs->size() != ys->size() || ys->size() != zs->size()) {
30  edm::LogError("RHStopDump") << "mismatch array sizes name/x/y/z:" << names->size() << '/' << xs->size() << '/'
31  << ys->size() << '/' << zs->size() << std::endl;
32  } else {
33  for (size_t i = 0; i < names->size(); ++i) {
34  mStream << (*names)[i] << ' ' << (*xs)[i] << ' ' << (*ys)[i] << ' ' << (*zs)[i] << ' ' << (*ts)[i] << std::endl;
35  mStream << (*ids)[i] << ' ' << (*masses)[i] << ' ' << (*charges)[i] << std::endl;
36  }
37  }
38 }
const edm::EDGetTokenT< std::vector< float > > tokenYs_
Definition: RHStopDump.h:23
std::ofstream mStream
Definition: RHStopDump.h:19
const edm::EDGetTokenT< std::vector< int > > tokenIds_
Definition: RHStopDump.h:26
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: RHStopDump.cc:19
Log< level::Error, false > LogError
const std::string names[nVars_]
RHStopDump(const edm::ParameterSet &)
Definition: RHStopDump.cc:7
const edm::EDGetTokenT< std::vector< float > > tokenXs_
Definition: RHStopDump.h:22
const edm::EDGetTokenT< std::vector< float > > tokenZs_
Definition: RHStopDump.h:24
const edm::EDGetTokenT< std::vector< float > > tokenCharges_
Definition: RHStopDump.h:28
HLT enums.
charges
only generated particles of these IDs are considered
const edm::EDGetTokenT< std::vector< std::string > > tokNames_
Definition: RHStopDump.h:21
const edm::EDGetTokenT< std::vector< float > > tokenTs_
Definition: RHStopDump.h:25
const edm::EDGetTokenT< std::vector< float > > tokenMasses_
Definition: RHStopDump.h:27