CMS 3D CMS Logo

StoreSecondary.cc
Go to the documentation of this file.
3 
6 
9 
11 
12 #include "G4Step.hh"
13 #include "G4Track.hh"
14 #include "G4HCofThisEvent.hh"
15 #include "CLHEP/Units/GlobalSystemOfUnits.h"
16 #include "CLHEP/Units/GlobalPhysicalConstants.h"
17 
18 #include <cmath>
19 #include <iostream>
20 #include <iomanip>
21 
23 
24  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("StoreSecondary");
25  treatSecondary = new TreatSecondary (m_p);
26 
27  produces<std::vector<math::XYZTLorentzVector> >("SecondaryMomenta");
28  produces<std::vector<int> >("SecondaryParticles");
29 
30  edm::LogInfo("CheckSecondary") << "Instantiate StoreSecondary to store "
31  << "secondaries after 1st hadronic inelastic"
32  << " interaction";
33 }
34 
36  delete treatSecondary;
37 }
38 
40 
41  std::unique_ptr<std::vector<math::XYZTLorentzVector> > secMom(new std::vector<math::XYZTLorentzVector>);
42  *secMom = secondaries;
43  e.put(std::move(secMom), "SecondaryMomenta");
44 
45  std::unique_ptr<std::vector<int> > secNumber(new std::vector<int>);
46  *secNumber = nsecs;
47  e.put(std::move(secNumber), "SecondaryParticles");
48 
49  LogDebug("CheckSecondary") << "StoreSecondary:: Event " << e.id() << " with "
50  << nsecs.size() << " hadronic collisions with "
51  << "secondaries produced in each step";
52  for (unsigned int i= 0; i < nsecs.size(); i++) {
53  LogDebug("CheckSecondary") << " " << nsecs[i] << " from " << procs[i];
54  }
55  LogDebug("CheckSecondary") << " and " << secondaries.size() << " secondaries"
56  << " produced in the first interactions:";
57  for (unsigned int i= 0; i < secondaries.size(); i++)
58  LogDebug("CheckSecondary") << "Secondary " << i << " " << secondaries[i];
59 }
60 
62 
63  nsecs.clear();
64  procs.clear();
65  secondaries.clear();
66 }
67 
69 
70  const G4Track * thTk = (*trk)();
72  if (nsecs.empty() && thTk->GetParentID() <= 0) storeIt = true;
73  else storeIt = false;
74  nHad = 0;
75 }
76 
77 void StoreSecondary::update(const G4Step * aStep) {
78 
80  int procID;
81  bool hadrInt;
82  double deltaE;
83  std::vector<int> charge;
84  std::vector<math::XYZTLorentzVector> tracks = treatSecondary->tracks(aStep,
85  name,
86  procID,
87  hadrInt,
88  deltaE,
89  charge);
90  if (hadrInt) {
91  ++nHad;
92  if (storeIt) {
93  int sec = (int)(tracks.size());
94  nsecs.push_back(sec);
95  procs.push_back(name);
96  if (nHad == 1) {
97  for (int i=0; i<sec; i++)
98  secondaries.push_back(tracks[i]);
99  }
100  }
101  }
102 }
#define LogDebug(id)
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
~StoreSecondary() override
std::vector< int > nsecs
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
StoreSecondary(const edm::ParameterSet &p)
std::vector< math::XYZTLorentzVector > tracks(const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges)
std::vector< math::XYZTLorentzVector > secondaries
std::vector< std::string > procs
void initTrack(const G4Track *trk)
TreatSecondary * treatSecondary
edm::EventID id() const
Definition: EventBase.h:60
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:510