CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SubsystemNeutronWriter.cc
Go to the documentation of this file.
12 
13 #include "CLHEP/Random/RandFlat.h"
14 
15 #include <algorithm>
16 
17 using namespace std;
18 
19 class SortByTime {
20 public:
21  bool operator()(const PSimHit & h1, const PSimHit & h2) {
22  return (h1.tof() < h2.tof());
23  }
24 };
25 
26 
28 : theHitWriter(0),
29  useRandFlat(false),
30  theInputTag(pset.getParameter<edm::InputTag>("input")),
31  theNeutronTimeCut(pset.getParameter<double>("neutronTimeCut")),
32  theTimeWindow(pset.getParameter<double>("timeWindow")),
33  theT0(pset.getParameter<double>("t0")),
34  theNEvents(0),
35  initialized(false),
36  useLocalDetId_(true)
37 {
38  string writer = pset.getParameter<string>("writer");
39  string output = pset.getParameter<string>("output");
40  if(writer == "ASCII")
41  {
42  theHitWriter = new AsciiNeutronWriter(output);
43  }
44  else if (writer == "ROOT")
45  {
46  theHitWriter = new RootNeutronWriter(output);
47  }
48  else if (writer == "EDM")
49  {
50  produces<edm::PSimHitContainer>();
52  // write out the real DetId, not the local one
53  useLocalDetId_ = false;
54  // smear the times
56  if ( ! rng.isAvailable()) {
57  throw cms::Exception("Configuration")
58  << "SubsystemNeutronWriter requires the RandomNumberGeneratorService\n"
59  "which is not present in the configuration file. You must add the service\n"
60  "in the configuration file or remove the modules that require it.";
61  }
62  useRandFlat = true;
63  }
64  else
65  {
66  throw cms::Exception("NeutronWriter") << "Bad writer: "
67  << writer;
68  }
69 }
70 
71 
73 {
74  printStats();
75  delete theHitWriter;
76 }
77 
78 
79 
81  edm::LogInfo("SubsystemNeutronWriter") << "SubsystemNeutronWriter Statistics:\n";
82  for(map<int,int>::iterator mapItr = theCountPerChamberType.begin();
83  mapItr != theCountPerChamberType.end(); ++mapItr) {
84  edm::LogInfo("SubsystemNeutronWriter") << " theEventOccupancy[" << mapItr->first << "] = "
85  << mapItr->second << " / " << theNEvents << " / NCT \n";
86  }
87 }
88 
89 
91 {
92  CLHEP::HepRandomEngine* engine = nullptr;
93  if(useRandFlat) {
95  engine = &rng->getEngine(e.streamID());
96  }
98  ++theNEvents;
100  e.getByLabel(theInputTag, hits);
101 
102  // sort hits by chamber
103  map<int, edm::PSimHitContainer> hitsByChamber;
104  for(edm::PSimHitContainer::const_iterator hitItr = hits->begin();
105  hitItr != hits->end(); ++hitItr)
106  {
107  int chamberIndex = chamberId(hitItr->detUnitId());
108  hitsByChamber[chamberIndex].push_back(*hitItr);
109  }
110 
111  // now write out each chamber's contents
112  for(map<int, edm::PSimHitContainer>::iterator hitsByChamberItr = hitsByChamber.begin();
113  hitsByChamberItr != hitsByChamber.end(); ++hitsByChamberItr)
114  {
115  int chambertype = chamberType(hitsByChamberItr->first);
116  writeHits(chambertype, hitsByChamberItr->second, engine);
117  }
119 }
120 
121 
123 {
124  // should instantiate one of every chamber type, just so
125  // ROOT knows what file to associate them with
126  theHitWriter->initialize(chamberType);
127 }
128 
129 
130 void SubsystemNeutronWriter::writeHits(int chamberType, edm::PSimHitContainer & chamberHits,
131  CLHEP::HepRandomEngine* engine)
132 {
133 
134  sort(chamberHits.begin(), chamberHits.end(), SortByTime());
135  edm::PSimHitContainer cluster;
136  float startTime = -1000.;
137  float smearing = 0.;
138  for(size_t i = 0; i < chamberHits.size(); ++i) {
139  PSimHit hit = chamberHits[i];
140  float tof = hit.tof();
141  LogDebug("SubsystemNeutronWriter") << "found hit from part type " << hit.particleType()
142  << " at tof " << tof << " p " << hit.pabs()
143  << " on det " << hit.detUnitId()
144  << " chamber type " << chamberType;
145  if(tof > theNeutronTimeCut) {
146  if(tof > (startTime + theTimeWindow) ) { // 1st in cluster
147  startTime = tof;
148  // set the time to be [t0, t0+25] at start of event
149  smearing = theT0;
150  if(useRandFlat) {
151  smearing += CLHEP::RandFlat::shoot(engine, 25.);
152  }
153  if(!cluster.empty()) {
154  LogDebug("SubsystemNeutronWriter") << "filling old cluster";
155  writeCluster(chamberType, cluster);
156  cluster.clear();
157  }
158  LogDebug("SubsystemNeutronWriter") << "starting neutron cluster at time " << startTime
159  << " on detType " << chamberType;
160  }
161  adjust(hit, -1.*startTime, smearing);
162  cluster.push_back( hit );
163  }
164  }
165  // get any leftover cluster
166  if(!cluster.empty()) {
167  writeCluster(chamberType, cluster);
168  }
169 }
170 
171 
172 void SubsystemNeutronWriter::writeCluster(int chamberType, const edm::PSimHitContainer & cluster)
173 {
174  if(accept(cluster))
175  {
176  theHitWriter->writeCluster(chamberType, cluster);
177  updateCount(chamberType);
178  }
179 }
180 
181 
182 void SubsystemNeutronWriter::adjust(PSimHit & h, float timeOffset, float smearing) {
183  unsigned int detId = useLocalDetId_ ? localDetId(h.detUnitId()) : h.detUnitId();
184  float htime = h.timeOfFlight() + timeOffset + smearing;
185  // prevent float precision loss
186  if (h.timeOfFlight() > 1.E+6) {
187  htime = smearing;
188  }
189  h = PSimHit( h.entryPoint(), h.exitPoint(), h.pabs(), htime,
190  h.energyLoss(), h.particleType(),
191  detId, h.trackId(),
192  h.momentumAtEntry().theta(),
193  h.momentumAtEntry().phi() );
194 }
195 
196 
197 void SubsystemNeutronWriter::updateCount(int chamberType) {
198  map<int,int>::iterator entry = theCountPerChamberType.find(chamberType);
199  if(entry == theCountPerChamberType.end()) {
200  theCountPerChamberType.insert( pair<int,int>(chamberType, 1) );
201  } else {
202  ++(entry->second);
203  }
204 }
205 
206 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:72
virtual int chamberType(int globalDetId) const =0
virtual ~SubsystemNeutronWriter()
destructor prints statistics on number of events written
void initialize(int chamberType)
good practice to do once for each chamber type
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:47
virtual bool accept(const edm::PSimHitContainer &cluster) const =0
decides whether this cluster is good enough to be included
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
virtual void beginEvent(edm::Event &e, const edm::EventSetup &es)
Definition: NeutronWriter.h:20
void adjust(PSimHit &h, float timeOffset, float smearing)
helper to add time offsets and local det ID
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
virtual int localDetId(int globalDetId) const =0
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
virtual int chamberId(int globalDetId) const =0
float timeOfFlight() const
Definition: PSimHit.h:69
virtual void endEvent()
Definition: NeutronWriter.h:21
virtual void writeHits(int chamberType, edm::PSimHitContainer &chamberHits, CLHEP::HepRandomEngine *)
bool operator()(const PSimHit &h1, const PSimHit &h2)
bool isAvailable() const
Definition: Service.h:46
std::map< int, int > theCountPerChamberType
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
virtual void writeCluster(int detType, const edm::PSimHitContainer &simHits)=0
writes out a list of SimHits.
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
virtual void initialize(int detType)
Definition: NeutronWriter.h:19
SubsystemNeutronWriter(edm::ParameterSet const &pset)
void updateCount(int chamberType)
updates the counter
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
void writeCluster(int chamberType, const edm::PSimHitContainer &cluster)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
int particleType() const
Definition: PSimHit.h:85
StreamID streamID() const
Definition: Event.h:74
unsigned int trackId() const
Definition: PSimHit.h:102
virtual void produce(edm::Event &e, edm::EventSetup const &c)
std::vector< PSimHit > PSimHitContainer
volatile std::atomic< bool > shutdown_flag false
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
unsigned int detUnitId() const
Definition: PSimHit.h:93