CMS 3D CMS Logo

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

const edm::EDGetTokenT< edm::PSimHitContainerhitToken_
 
bool initialized
 
std::map< int, int > theCountPerChamberType
 
NeutronWritertheHitWriter
 
const edm::InputTag theInputTag
 
const double theNeutronTimeCut
 
int theNEvents
 
const double theT0
 
const 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::SubsystemNeutronWriter ( edm::ParameterSet const &  pset)
explicit

Definition at line 24 of file SubsystemNeutronWriter.cc.

References Exception, edm::Service< T >::isAvailable(), convertSQLitetoXML_cfg::output, muonDTDigis_cfi::pset, theHitWriter, useLocalDetId_, useRandFlat, and cscNeutronWriter_cfi::writer.

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  hitToken_(consumes<edm::PSimHitContainer>(theInputTag)),
32  theNEvents(0),
33  initialized(false),
34  useLocalDetId_(true) {
35  string writer = pset.getParameter<string>("writer");
36  string output = pset.getParameter<string>("output");
37  if (writer == "ASCII") {
39  } else if (writer == "ROOT") {
41  } else if (writer == "EDM") {
42  produces<edm::PSimHitContainer>();
44  // write out the real DetId, not the local one
45  useLocalDetId_ = false;
46  // smear the times
48  if (!rng.isAvailable()) {
49  throw cms::Exception("Configuration")
50  << "SubsystemNeutronWriter requires the RandomNumberGeneratorService\n"
51  "which is not present in the configuration file. You must add the service\n"
52  "in the configuration file or remove the modules that require it.";
53  }
54  useRandFlat = true;
55  } else {
56  throw cms::Exception("NeutronWriter") << "Bad writer: " << writer;
57  }
58 }
const edm::InputTag theInputTag
const edm::EDGetTokenT< edm::PSimHitContainer > hitToken_
bool isAvailable() const
Definition: Service.h:40

◆ ~SubsystemNeutronWriter()

SubsystemNeutronWriter::~SubsystemNeutronWriter ( )
override

destructor prints statistics on number of events written

Definition at line 60 of file SubsystemNeutronWriter.cc.

References printStats(), and theHitWriter.

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

Member Function Documentation

◆ accept()

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().

◆ adjust()

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 h, localDetId(), hgcalLayerClusters_cfi::timeOffset, 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 }
virtual int localDetId(int globalDetId) const =0
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4

◆ chamberId()

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

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by produce().

◆ chamberType()

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

◆ initialize()

void SubsystemNeutronWriter::initialize ( int  chamberType)

good practice to do once for each chamber type

Definition at line 101 of file SubsystemNeutronWriter.cc.

References chamberType(), 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

◆ localDetId()

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

Implemented in CSCNeutronWriter, RPCNeutronWriter, and DTNeutronWriter.

Referenced by adjust().

◆ printStats()

void SubsystemNeutronWriter::printStats ( )

Definition at line 65 of file SubsystemNeutronWriter.cc.

References theCountPerChamberType, and theNEvents.

Referenced by ~SubsystemNeutronWriter().

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

◆ produce()

void SubsystemNeutronWriter::produce ( edm::Event e,
edm::EventSetup const &  c 
)
override

Definition at line 74 of file SubsystemNeutronWriter.cc.

References NeutronWriter::beginEvent(), c, chamberId(), chamberType(), MillePedeFileConverter_cfg::e, NeutronWriter::endEvent(), edm::RandomNumberGenerator::getEngine(), hfClusterShapes_cfi::hits, hitToken_, theHitWriter, theNEvents, useRandFlat, and writeHits().

74  {
75  CLHEP::HepRandomEngine* engine = nullptr;
76  if (useRandFlat) {
78  engine = &rng->getEngine(e.streamID());
79  }
81  ++theNEvents;
83 
84  // sort hits by chamber
85  std::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 (std::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 }
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 *)
const edm::EDGetTokenT< edm::PSimHitContainer > hitToken_

◆ updateCount()

void SubsystemNeutronWriter::updateCount ( int  chamberType)
protected

updates the counter

Definition at line 172 of file SubsystemNeutronWriter.cc.

References chamberType(), 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

◆ writeCluster()

void SubsystemNeutronWriter::writeCluster ( int  chamberType,
const edm::PSimHitContainer cluster 
)
protected

Definition at line 146 of file SubsystemNeutronWriter.cc.

References accept(), chamberType(), 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

◆ writeHits()

void SubsystemNeutronWriter::writeHits ( int  chamberType,
edm::PSimHitContainer chamberHits,
CLHEP::HepRandomEngine *  engine 
)
protectedvirtual

Definition at line 107 of file SubsystemNeutronWriter.cc.

References adjust(), chamberType(), mps_fire::i, LogDebug, jetUpdater_cfi::sort, startTime, theNeutronTimeCut, theT0, theTimeWindow, 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 }
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
void writeCluster(int chamberType, const edm::PSimHitContainer &cluster)
std::vector< PSimHit > PSimHitContainer
#define LogDebug(id)

Member Data Documentation

◆ hitToken_

const edm::EDGetTokenT<edm::PSimHitContainer> SubsystemNeutronWriter::hitToken_
private

Definition at line 69 of file SubsystemNeutronWriter.h.

Referenced by produce().

◆ initialized

bool SubsystemNeutronWriter::initialized
private

Definition at line 71 of file SubsystemNeutronWriter.h.

◆ theCountPerChamberType

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

Definition at line 74 of file SubsystemNeutronWriter.h.

Referenced by printStats(), and updateCount().

◆ theHitWriter

NeutronWriter* SubsystemNeutronWriter::theHitWriter
private

◆ theInputTag

const edm::InputTag SubsystemNeutronWriter::theInputTag
private

Definition at line 65 of file SubsystemNeutronWriter.h.

◆ theNeutronTimeCut

const double SubsystemNeutronWriter::theNeutronTimeCut
private

Definition at line 66 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

◆ theNEvents

int SubsystemNeutronWriter::theNEvents
private

Definition at line 70 of file SubsystemNeutronWriter.h.

Referenced by printStats(), and produce().

◆ theT0

const double SubsystemNeutronWriter::theT0
private

Definition at line 68 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

◆ theTimeWindow

const double SubsystemNeutronWriter::theTimeWindow
private

Definition at line 67 of file SubsystemNeutronWriter.h.

Referenced by writeHits().

◆ useLocalDetId_

bool SubsystemNeutronWriter::useLocalDetId_
private

Definition at line 73 of file SubsystemNeutronWriter.h.

Referenced by adjust(), and SubsystemNeutronWriter().

◆ useRandFlat

bool SubsystemNeutronWriter::useRandFlat
private

Definition at line 64 of file SubsystemNeutronWriter.h.

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