CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
StoreSecondary.cc
Go to the documentation of this file.
2 
5 
8 
10 
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 #include "G4VProcess.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  // produces<std::vector<std::string> >("SecondaryProcesses");
30 
31  edm::LogInfo("CheckSecondary") << "Instantiate StoreSecondary to store "
32  << "secondaries after 1st hadronic inelastic"
33  << " interaction";
34 }
35 
37  delete treatSecondary;
38 }
39 
41 
42  std::auto_ptr<std::vector<math::XYZTLorentzVector> > secMom(new std::vector<math::XYZTLorentzVector>);
43  *secMom = secondaries;
44  e.put(secMom, "SecondaryMomenta");
45 
46  std::auto_ptr<std::vector<int> > secNumber(new std::vector<int>);
47  *secNumber = nsecs;
48  e.put(secNumber, "SecondaryParticles");
49 
50  /*
51  std::auto_ptr<std::vector<std::string> > secProc(new std::vector<std::string>);
52  *secProc = procs;
53  e.put(secProc, "SecondaryProcesses");
54  */
55 
56  LogDebug("CheckSecondary") << "StoreSecondary:: Event " << e.id() << " with "
57  << nsecs.size() << " hadronic collisions with "
58  << "secondaries produced in each step";
59  for (unsigned int i= 0; i < nsecs.size(); i++)
60  LogDebug("CheckSecondary") << " " << nsecs[i] << " from " << procs[i];
61  LogDebug("CheckSecondary") << " and " << secondaries.size() << " secondaries"
62  << " produced in the first interactions:";
63  for (unsigned int i= 0; i < secondaries.size(); i++)
64  LogDebug("CheckSecondary") << "Secondary " << i << " " << secondaries[i];
65 }
66 
68 
69  nsecs.clear();
70  procs.clear();
71  secondaries.clear();
72 }
73 
75 
76  const G4Track * thTk = (*trk)();
78  if (nsecs.size() == 0 && thTk->GetParentID() <= 0) storeIt = true;
79  else storeIt = false;
80  nHad = 0;
81 }
82 
83 void StoreSecondary::update(const G4Step * aStep) {
84 
86  int procID;
87  bool hadrInt;
88  double deltaE;
89  std::vector<int> charge;
90  std::vector<math::XYZTLorentzVector> tracks = treatSecondary->tracks(aStep,
91  name,
92  procID,
93  hadrInt,
94  deltaE,
95  charge);
96  if (hadrInt) {
97  nHad++;
98  if (storeIt) {
99  int sec = (int)(tracks.size());
100  nsecs.push_back(sec);
101  procs.push_back(name);
102  if (nHad == 1) {
103  for (int i=0; i<sec; i++)
104  secondaries.push_back(tracks[i]);
105  }
106  }
107  }
108 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual ~StoreSecondary()
std::vector< int > nsecs
void update(const BeginOfEvent *evt)
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
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
std::vector< std::string > procs
void produce(edm::Event &, const edm::EventSetup &)
tuple tracks
Definition: testEve_cfg.py:39
void initTrack(const G4Track *trk)
TreatSecondary * treatSecondary
edm::EventID id() const
Definition: EventBase.h:60