CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EcalBarrelClusterFastTimer Class Reference
Inheritance diagram for EcalBarrelClusterFastTimer:
edm::global::EDProducer<> edm::global::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 EcalBarrelClusterFastTimer (const edm::ParameterSet &)
 
virtual void produce (edm::StreamID, edm::Event &, const edm::EventSetup &) const override
 
 ~EcalBarrelClusterFastTimer ()
 
- Public Member Functions inherited from edm::global::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, std::unordered_multimap< std::string, edm::ProductResolverIndex > const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

float correctTimeToVertex (const float intime, const DetId &timeDet, const reco::Vertex &vtx, const CaloSubdetectorGeometry *ecalGeom) const
 
std::pair< float, DetIdgetTimeForECALPFCluster (const EcalRecHitCollection &, const reco::PFCluster &) const
 

Private Attributes

std::vector< std::unique_ptr< const ResolutionModel > > _resolutions
 
edm::EDGetTokenT< std::vector< reco::PFCluster > > ebClustersToken_
 
edm::EDGetTokenT< EcalRecHitCollectionebTimeHitsToken_
 
const unsigned ecalDepth_
 
const float minEnergy_
 
const float minFraction_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::global::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 37 of file EcalBarrelClusterFastTimer.cc.

Constructor & Destructor Documentation

EcalBarrelClusterFastTimer::EcalBarrelClusterFastTimer ( const edm::ParameterSet conf)

Definition at line 78 of file EcalBarrelClusterFastTimer.cc.

References _resolutions, Exception, reco::get(), edm::ParameterSet::getParameterSetVector(), dataset::name, fftjetproducer_cfi::resolution, and AlCaHLTBitMon_QueryRunRegistry::string.

78  :
79  ebTimeHitsToken_(consumes<EcalRecHitCollection>( conf.getParameter<edm::InputTag>("ebTimeHits") ) ),
80  ebClustersToken_(consumes<std::vector<reco::PFCluster> >( conf.getParameter<edm::InputTag>("ebClusters") ) ),
81  minFraction_( conf.getParameter<double>("minFractionToConsider") ),
82  minEnergy_(conf.getParameter<double>("minEnergyToConsider") ),
83  ecalDepth_(conf.getParameter<double>("ecalDepth") )
84 {
85  // setup resolution models
86  const std::vector<edm::ParameterSet>& resos = conf.getParameterSetVector("resolutionModels");
87  for( const auto& reso : resos ) {
88  const std::string& name = reso.getParameter<std::string>("modelName");
89  ResolutionModel* resomod = ResolutionModelFactory::get()->create(name,reso);
90  _resolutions.emplace_back( resomod );
91 
92  // times and time resolutions for general tracks
93  produces<edm::ValueMap<float> >(name);
94  produces<edm::ValueMap<float> >(name+resolution);
95  }
96  // get RNG engine
98  if (!rng.isAvailable()){
99  throw cms::Exception("Configuration")
100  << "EcalBarrelClusterFastTimer::EcalBarrelClusterFastTimer() - RandomNumberGeneratorService is not present in configuration file.\n"
101  << "Add the service in the configuration file or remove the modules that require it.";
102  }
103 }
T getParameter(std::string const &) const
VParameterSet const & getParameterSetVector(std::string const &name) const
edm::EDGetTokenT< std::vector< reco::PFCluster > > ebClustersToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< std::unique_ptr< const ResolutionModel > > _resolutions
edm::EDGetTokenT< EcalRecHitCollection > ebTimeHitsToken_
T get(const Candidate &c)
Definition: component.h:55
EcalBarrelClusterFastTimer::~EcalBarrelClusterFastTimer ( )
inline

Definition at line 40 of file EcalBarrelClusterFastTimer.cc.

References produce().

40 { }

Member Function Documentation

float EcalBarrelClusterFastTimer::correctTimeToVertex ( const float  intime,
const DetId timeDet,
const reco::Vertex vtx,
const CaloSubdetectorGeometry ecalGeom 
) const
private

Definition at line 176 of file EcalBarrelClusterFastTimer.cc.

References ecalDepth_, Exception, CaloSubdetectorGeometry::getGeometry(), reco::Vertex::position(), DetId::rawId(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

177  {
178  if( timeDet.rawId() == 0 ) return -999.;
179  // correct the cluster time from 0,0,0 to the primary vertex given
180  const CaloCellGeometry* cellGeometry ( ecalGeom->getGeometry( timeDet ) ) ;
181  if( nullptr == cellGeometry ) {
182  throw cms::Exception("BadECALBarrelCell")
183  << timeDet << " is not a valid ECAL Barrel DetId!";
184  }
185  //depth in mm in the middle of the layer position;
186  GlobalPoint layerPos = (dynamic_cast<const TruncatedPyramid*>(cellGeometry))->getPosition( ecalDepth_+0.5 );
187  const math::XYZPoint layerPos_cm( layerPos.x(), layerPos.y(), layerPos.z() );
188  const math::XYZVector to_center = layerPos_cm - math::XYZPoint(0.,0.,0.);
189  const math::XYZVector to_vtx = layerPos_cm - vtx.position();
190 
191  /*
192  std::cout << intime << ' ' << to_center.r()/cm_per_ns << ' ' << to_vtx.r()/cm_per_ns
193  << ' ' << intime + to_center.r()/cm_per_ns - to_vtx.r()/cm_per_ns << std::endl;
194  */
195 
196  return intime + to_center.r()/cm_per_ns - to_vtx.r()/cm_per_ns;
197 }
T y() const
Definition: PV3DBase.h:63
const Point & position() const
position
Definition: Vertex.h:109
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
T z() const
Definition: PV3DBase.h:64
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
A base class to handle the particular shape of Ecal Xtals. Taken from ORCA Calorimetry Code...
T x() const
Definition: PV3DBase.h:62
std::pair< float, DetId > EcalBarrelClusterFastTimer::getTimeForECALPFCluster ( const EcalRecHitCollection timehits,
const reco::PFCluster cluster 
) const
private

Definition at line 147 of file EcalBarrelClusterFastTimer.cc.

References edm::SortedCollection< T, SORT >::find(), dedxEstimators_cff::fraction, hpstanc_transforms::max, minEnergy_, minFraction_, and reco::PFCluster::recHitFractions().

Referenced by produce().

147  {
148 
149  const auto& rhfs = cluster.recHitFractions();
150  unsigned best_hit = 0;
151  float best_energy = -1.f;
152  for( const auto& rhf : rhfs ) {
153  const auto& hitref = rhf.recHitRef();
154  const unsigned detid = hitref->detId();
155  const auto fraction = rhf.fraction();
156  const auto energy = hitref->energy();
157  if( fraction < minFraction_ || energy < minEnergy_ ) continue;
158  auto timehit = timehits.find(detid);
159  float e_hit = rhf.fraction() * hitref->energy();
160  if( e_hit > best_energy && timehit->isTimeValid() ) {
161  best_energy = e_hit;
162  best_hit = detid;
163  }
164  }
165 
166  float best_time_guess = std::numeric_limits<float>::max();
167  if( best_energy > 0. ) {
168  best_time_guess = timehits.find(best_hit)->time();
169  }
170 
171  //std::cout << "EcalBarrelFastTimer: " << best_time_guess << ' ' << best_energy << ' ' << best_hit << std::endl;
172 
173  return std::make_pair(best_time_guess,DetId(best_hit));
174 }
Definition: DetId.h:18
iterator find(key_type k)
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
void EcalBarrelClusterFastTimer::produce ( edm::StreamID  sid,
edm::Event evt,
const edm::EventSetup es 
) const
overridevirtual

Definition at line 105 of file EcalBarrelClusterFastTimer.cc.

References _resolutions, fastPrimaryVertexProducer_cfi::clusters, ebClustersToken_, ebTimeHitsToken_, plotBeamSpotDB::first, edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), getTimeForECALPFCluster(), mps_fire::i, dataset::name, fftjetproducer_cfi::resolution, electronProducer_cfi::resolutions, AlCaHLTBitMon_QueryRunRegistry::string, and create_public_lumi_plots::times.

Referenced by ~EcalBarrelClusterFastTimer().

105  {
106  // get RNG engine
108  auto rng_engine = &(rng->getEngine(sid));
109 
112 
113  evt.getByToken(ebClustersToken_,clustersH);
114  evt.getByToken(ebTimeHitsToken_,timehitsH);
115 
116  const auto& clusters = *clustersH;
117  const auto& timehits = *timehitsH;
118 
119  std::vector<std::pair<float,DetId> > times; // perfect times keyed to cluster index
120  times.reserve(clusters.size());
121 
122  for( const auto& cluster : clusters ) {
123  times.push_back( getTimeForECALPFCluster( timehits, cluster ) );
124  }
125 
126  for( const auto& reso : _resolutions ) {
127  const std::string& name = reso->name();
128  std::vector<float> resolutions;
129  std::vector<float> smeared_times;
130  resolutions.reserve(clusters.size());
131  smeared_times.reserve(clusters.size());
132 
133  // smear once then correct to multiple vertices
134  for( unsigned i = 0 ; i < clusters.size(); ++i ) {
135  const float theresolution = reso->getTimeResolution(clusters[i]);
136 
137  smeared_times.emplace_back( CLHEP::RandGauss::shoot(rng_engine, times[i].first, theresolution) );
138  resolutions.push_back( theresolution );
139  }
140 
141  writeValueMap(evt,clustersH,smeared_times,name);
142  writeValueMap(evt,clustersH,resolutions,name+resolution);
143  }
144 
145 }
edm::EDGetTokenT< std::vector< reco::PFCluster > > ebClustersToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
std::vector< std::unique_ptr< const ResolutionModel > > _resolutions
edm::EDGetTokenT< EcalRecHitCollection > ebTimeHitsToken_
std::pair< float, DetId > getTimeForECALPFCluster(const EcalRecHitCollection &, const reco::PFCluster &) const

Member Data Documentation

std::vector<std::unique_ptr<const ResolutionModel> > EcalBarrelClusterFastTimer::_resolutions
private

Definition at line 49 of file EcalBarrelClusterFastTimer.cc.

Referenced by EcalBarrelClusterFastTimer(), and produce().

edm::EDGetTokenT<std::vector<reco::PFCluster> > EcalBarrelClusterFastTimer::ebClustersToken_
private

Definition at line 47 of file EcalBarrelClusterFastTimer.cc.

Referenced by produce().

edm::EDGetTokenT<EcalRecHitCollection> EcalBarrelClusterFastTimer::ebTimeHitsToken_
private

Definition at line 46 of file EcalBarrelClusterFastTimer.cc.

Referenced by produce().

const unsigned EcalBarrelClusterFastTimer::ecalDepth_
private

Definition at line 51 of file EcalBarrelClusterFastTimer.cc.

Referenced by correctTimeToVertex().

const float EcalBarrelClusterFastTimer::minEnergy_
private

Definition at line 50 of file EcalBarrelClusterFastTimer.cc.

Referenced by getTimeForECALPFCluster().

const float EcalBarrelClusterFastTimer::minFraction_
private

Definition at line 50 of file EcalBarrelClusterFastTimer.cc.

Referenced by getTimeForECALPFCluster().