CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Protected Member Functions | Private Attributes
SubsystemNeutronWriter Class Referenceabstract

doesn't have to be a producer. Can act as an analyzer, too. More...

#include <SubsystemNeutronWriter.h>

Inheritance diagram for SubsystemNeutronWriter:
edm::stream::EDProducer<> CSCNeutronWriter DTNeutronWriter RPCNeutronWriter

Public Member Functions

virtual bool accept (const edm::PSimHitContainer &cluster) const =0
 decides whether this cluster is good enough to be included More...
 
virtual int chamberId (int globalDetId) const =0
 
virtual int chamberType (int globalDetId) const =0
 
void initialize (int chamberType)
 good practice to do once for each chamber type More...
 
virtual int localDetId (int globalDetId) const =0
 
void printStats ()
 
void produce (edm::Event &e, edm::EventSetup const &c) override
 
 SubsystemNeutronWriter (edm::ParameterSet const &pset)
 
 ~SubsystemNeutronWriter () override
 destructor prints statistics on number of events written More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void adjust (PSimHit &h, float timeOffset, float smearing)
 helper to add time offsets and local det ID More...
 
void updateCount (int chamberType)
 updates the counter More...
 
void writeCluster (int chamberType, const edm::PSimHitContainer &cluster)
 
virtual void writeHits (int chamberType, edm::PSimHitContainer &chamberHits, CLHEP::HepRandomEngine *)
 

Private Attributes

bool initialized
 
std::map< int, int > theCountPerChamberType
 
NeutronWritertheHitWriter
 
edm::InputTag theInputTag
 
double theNeutronTimeCut
 
int theNEvents
 
double theT0
 
double theTimeWindow
 
bool useLocalDetId_
 
bool useRandFlat
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T...>
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T...>
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

doesn't have to be a producer. Can act as an analyzer, too.

Definition at line 28 of file SubsystemNeutronWriter.h.

Constructor & Destructor Documentation

SubsystemNeutronWriter::SubsystemNeutronWriter ( edm::ParameterSet const &  pset)
explicit

Definition at line 24 of file SubsystemNeutronWriter.cc.

References Exception, edm::ParameterSet::getParameter(), edm::Service< T >::isAvailable(), convertSQLitetoXML_cfg::output, theHitWriter, useLocalDetId_, and useRandFlat.

25  : theHitWriter(nullptr),
26  useRandFlat(false),
27  theInputTag(pset.getParameter<edm::InputTag>("input")),
28  theNeutronTimeCut(pset.getParameter<double>("neutronTimeCut")),
29  theTimeWindow(pset.getParameter<double>("timeWindow")),
30  theT0(pset.getParameter<double>("t0")),
31  theNEvents(0),
32  initialized(false),
33  useLocalDetId_(true) {
34  string writer = pset.getParameter<string>("writer");
35  string output = pset.getParameter<string>("output");
36  if (writer == "ASCII") {
37  theHitWriter = new AsciiNeutronWriter(output);
38  } else if (writer == "ROOT") {
39  theHitWriter = new RootNeutronWriter(output);
40  } else if (writer == "EDM") {
41  produces<edm::PSimHitContainer>();
43  // write out the real DetId, not the local one
44  useLocalDetId_ = false;
45  // smear the times
47  if (!rng.isAvailable()) {
48  throw cms::Exception("Configuration")
49  << "SubsystemNeutronWriter requires the RandomNumberGeneratorService\n"
50  "which is not present in the configuration file. You must add the service\n"
51  "in the configuration file or remove the modules that require it.";
52  }
53  useRandFlat = true;
54  } else {
55  throw cms::Exception("NeutronWriter") << "Bad writer: " << writer;
56  }
57 }
bool isAvailable() const
Definition: Service.h:40
SubsystemNeutronWriter::~SubsystemNeutronWriter ( )
override

destructor prints statistics on number of events written

Definition at line 59 of file SubsystemNeutronWriter.cc.

References printStats(), and theHitWriter.

59  {
60  printStats();
61  delete theHitWriter;
62 }

Member Function Documentation

virtual bool SubsystemNeutronWriter::accept ( const edm::PSimHitContainer cluster) const
pure virtual

decides whether this cluster is good enough to be included

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by esMonitoring.FDJsonServer::handle_accept(), and writeCluster().

void SubsystemNeutronWriter::adjust ( PSimHit h,
float  timeOffset,
float  smearing 
)
protected

helper to add time offsets and local det ID

Definition at line 153 of file SubsystemNeutronWriter.cc.

References PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), localDetId(), PSimHit::momentumAtEntry(), PSimHit::pabs(), PSimHit::particleType(), PV3DBase< T, PVType, FrameType >::phi(), PV3DBase< T, PVType, FrameType >::theta(), PSimHit::timeOfFlight(), PSimHit::trackId(), and useLocalDetId_.

Referenced by writeHits().

153  {
154  unsigned int detId = useLocalDetId_ ? localDetId(h.detUnitId()) : h.detUnitId();
155  float htime = h.timeOfFlight() + timeOffset + smearing;
156  // prevent float precision loss
157  if (h.timeOfFlight() > 1.E+6) {
158  htime = smearing;
159  }
160  h = PSimHit(h.entryPoint(),
161  h.exitPoint(),
162  h.pabs(),
163  htime,
164  h.energyLoss(),
165  h.particleType(),
166  detId,
167  h.trackId(),
168  h.momentumAtEntry().theta(),
169  h.momentumAtEntry().phi());
170 }
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:55
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
virtual int localDetId(int globalDetId) const =0
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
float timeOfFlight() const
Definition: PSimHit.h:73
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
int particleType() const
Definition: PSimHit.h:89
unsigned int trackId() const
Definition: PSimHit.h:106
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
unsigned int detUnitId() const
Definition: PSimHit.h:97
virtual int SubsystemNeutronWriter::chamberId ( int  globalDetId) const
pure virtual

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by produce().

virtual int SubsystemNeutronWriter::chamberType ( int  globalDetId) const
pure virtual

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by produce(), and writeHits().

void SubsystemNeutronWriter::initialize ( int  chamberType)

good practice to do once for each chamber type

Definition at line 101 of file SubsystemNeutronWriter.cc.

References NeutronWriter::initialize(), and theHitWriter.

Referenced by CSCNeutronWriter::CSCNeutronWriter().

101  {
102  // should instantiate one of every chamber type, just so
103  // ROOT knows what file to associate them with
105 }
virtual int chamberType(int globalDetId) const =0
virtual void initialize(int detType)
Definition: NeutronWriter.h:19
virtual int SubsystemNeutronWriter::localDetId ( int  globalDetId) const
pure virtual

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by adjust().

void SubsystemNeutronWriter::printStats ( )

Definition at line 64 of file SubsystemNeutronWriter.cc.

References theCountPerChamberType, and theNEvents.

Referenced by ~SubsystemNeutronWriter().

64  {
65  edm::LogInfo("SubsystemNeutronWriter") << "SubsystemNeutronWriter Statistics:\n";
66  for (map<int, int>::iterator mapItr = theCountPerChamberType.begin(); mapItr != theCountPerChamberType.end();
67  ++mapItr) {
68  edm::LogInfo("SubsystemNeutronWriter")
69  << " theEventOccupancy[" << mapItr->first << "] = " << mapItr->second << " / " << theNEvents << " / NCT \n";
70  }
71 }
std::map< int, int > theCountPerChamberType
Log< level::Info, false > LogInfo
void SubsystemNeutronWriter::produce ( edm::Event e,
edm::EventSetup const &  c 
)
override

Definition at line 73 of file SubsystemNeutronWriter.cc.

References NeutronWriter::beginEvent(), chamberId(), chamberType(), NeutronWriter::endEvent(), edm::Event::getByLabel(), edm::RandomNumberGenerator::getEngine(), edm::Event::streamID(), theHitWriter, theInputTag, theNEvents, useRandFlat, and writeHits().

73  {
74  CLHEP::HepRandomEngine* engine = nullptr;
75  if (useRandFlat) {
77  engine = &rng->getEngine(e.streamID());
78  }
80  ++theNEvents;
82  e.getByLabel(theInputTag, hits);
83 
84  // sort hits by chamber
85  map<int, edm::PSimHitContainer> hitsByChamber;
86  for (edm::PSimHitContainer::const_iterator hitItr = hits->begin(); hitItr != hits->end(); ++hitItr) {
87  int chamberIndex = chamberId(hitItr->detUnitId());
88  hitsByChamber[chamberIndex].push_back(*hitItr);
89  }
90 
91  // now write out each chamber's contents
92  for (map<int, edm::PSimHitContainer>::iterator hitsByChamberItr = hitsByChamber.begin();
93  hitsByChamberItr != hitsByChamber.end();
94  ++hitsByChamberItr) {
95  int chambertype = chamberType(hitsByChamberItr->first);
96  writeHits(chambertype, hitsByChamberItr->second, engine);
97  }
99 }
const edm::EventSetup & c
virtual int chamberType(int globalDetId) const =0
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
virtual void beginEvent(edm::Event &e, const edm::EventSetup &es)
Definition: NeutronWriter.h:20
virtual int chamberId(int globalDetId) const =0
virtual void endEvent()
Definition: NeutronWriter.h:21
virtual void writeHits(int chamberType, edm::PSimHitContainer &chamberHits, CLHEP::HepRandomEngine *)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
StreamID streamID() const
Definition: Event.h:98
void SubsystemNeutronWriter::updateCount ( int  chamberType)
protected

updates the counter

Definition at line 172 of file SubsystemNeutronWriter.cc.

References mps_splice::entry, and theCountPerChamberType.

Referenced by writeCluster().

172  {
173  map<int, int>::iterator entry = theCountPerChamberType.find(chamberType);
174  if (entry == theCountPerChamberType.end()) {
175  theCountPerChamberType.insert(pair<int, int>(chamberType, 1));
176  } else {
177  ++(entry->second);
178  }
179 }
virtual int chamberType(int globalDetId) const =0
std::map< int, int > theCountPerChamberType
list entry
Definition: mps_splice.py:68
void SubsystemNeutronWriter::writeCluster ( int  chamberType,
const edm::PSimHitContainer cluster 
)
protected

Definition at line 146 of file SubsystemNeutronWriter.cc.

References accept(), theHitWriter, updateCount(), and NeutronWriter::writeCluster().

Referenced by writeHits().

146  {
147  if (accept(cluster)) {
150  }
151 }
virtual int chamberType(int globalDetId) const =0
virtual bool accept(const edm::PSimHitContainer &cluster) const =0
decides whether this cluster is good enough to be included
virtual void writeCluster(int detType, const edm::PSimHitContainer &simHits)=0
writes out a list of SimHits.
void updateCount(int chamberType)
updates the counter
void SubsystemNeutronWriter::writeHits ( int  chamberType,
edm::PSimHitContainer chamberHits,
CLHEP::HepRandomEngine *  engine 
)
protectedvirtual

Definition at line 107 of file SubsystemNeutronWriter.cc.

References adjust(), chamberType(), PSimHit::detUnitId(), mps_fire::i, LogDebug, PSimHit::pabs(), PSimHit::particleType(), startTime, theNeutronTimeCut, theT0, theTimeWindow, PSimHit::tof(), useRandFlat, and writeCluster().

Referenced by produce().

109  {
110  sort(chamberHits.begin(), chamberHits.end(), SortByTime());
111  edm::PSimHitContainer cluster;
112  float startTime = -1000.;
113  float smearing = 0.;
114  for (size_t i = 0; i < chamberHits.size(); ++i) {
115  PSimHit hit = chamberHits[i];
116  float tof = hit.tof();
117  LogDebug("SubsystemNeutronWriter") << "found hit from part type " << hit.particleType() << " at tof " << tof
118  << " p " << hit.pabs() << " on det " << hit.detUnitId() << " chamber type "
119  << chamberType;
120  if (tof > theNeutronTimeCut) {
121  if (tof > (startTime + theTimeWindow)) { // 1st in cluster
122  startTime = tof;
123  // set the time to be [t0, t0+25] at start of event
124  smearing = theT0;
125  if (useRandFlat) {
126  smearing += CLHEP::RandFlat::shoot(engine, 25.);
127  }
128  if (!cluster.empty()) {
129  LogDebug("SubsystemNeutronWriter") << "filling old cluster";
130  writeCluster(chamberType, cluster);
131  cluster.clear();
132  }
133  LogDebug("SubsystemNeutronWriter")
134  << "starting neutron cluster at time " << startTime << " on detType " << chamberType;
135  }
136  adjust(hit, -1. * startTime, smearing);
137  cluster.push_back(hit);
138  }
139  }
140  // get any leftover cluster
141  if (!cluster.empty()) {
142  writeCluster(chamberType, cluster);
143  }
144 }
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:76
virtual int chamberType(int globalDetId) const =0
void adjust(PSimHit &h, float timeOffset, float smearing)
helper to add time offsets and local det ID
string startTime
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
void writeCluster(int chamberType, const edm::PSimHitContainer &cluster)
int particleType() const
Definition: PSimHit.h:89
std::vector< PSimHit > PSimHitContainer
unsigned int detUnitId() const
Definition: PSimHit.h:97
#define LogDebug(id)

Member Data Documentation

bool SubsystemNeutronWriter::initialized
private

Definition at line 70 of file SubsystemNeutronWriter.h.

std::map<int, int> SubsystemNeutronWriter::theCountPerChamberType
private

Definition at line 73 of file SubsystemNeutronWriter.h.

Referenced by printStats(), and updateCount().

NeutronWriter* SubsystemNeutronWriter::theHitWriter
private
edm::InputTag SubsystemNeutronWriter::theInputTag
private

Definition at line 65 of file SubsystemNeutronWriter.h.

Referenced by produce().

double SubsystemNeutronWriter::theNeutronTimeCut
private

Definition at line 66 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

int SubsystemNeutronWriter::theNEvents
private

Definition at line 69 of file SubsystemNeutronWriter.h.

Referenced by printStats(), and produce().

double SubsystemNeutronWriter::theT0
private

Definition at line 68 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

double SubsystemNeutronWriter::theTimeWindow
private

Definition at line 67 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

bool SubsystemNeutronWriter::useLocalDetId_
private

Definition at line 72 of file SubsystemNeutronWriter.h.

Referenced by adjust(), and SubsystemNeutronWriter().

bool SubsystemNeutronWriter::useRandFlat
private

Definition at line 64 of file SubsystemNeutronWriter.h.

Referenced by produce(), SubsystemNeutronWriter(), and writeHits().