CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DTRecSegment2DExtendedProducer Class Reference

#include <DTRecSegment2DExtendedProducer.h>

Inheritance diagram for DTRecSegment2DExtendedProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 DTRecSegment2DExtendedProducer (const edm::ParameterSet &)
 Constructor. More...
 
virtual void produce (edm::Event &event, const edm::EventSetup &setup)
 The method which produces the 2D-segments. More...
 
virtual ~DTRecSegment2DExtendedProducer ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

bool debug
 
edm::EDGetTokenT
< DTRecClusterCollection
recClusToken_
 
edm::EDGetTokenT
< DTRecHitCollection
recHits1DToken_
 
DTCombinatorialExtendedPatternRecotheAlgo
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- 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

Producer for DT segment in one projection.

Author
Stefano Lacaprara - INFN Legnaro stefa.nosp@m.no.l.nosp@m.acapr.nosp@m.ara@.nosp@m.pd.in.nosp@m.fn.i.nosp@m.t
Riccardo Bellan - INFN TO ricca.nosp@m.rdo..nosp@m.bella.nosp@m.n@ce.nosp@m.rn.ch

Definition at line 32 of file DTRecSegment2DExtendedProducer.h.

Constructor & Destructor Documentation

DTRecSegment2DExtendedProducer::DTRecSegment2DExtendedProducer ( const edm::ParameterSet pset)

Constructor.

Definition at line 35 of file DTRecSegment2DExtendedProducer.cc.

References gather_cfg::cout, debug, edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

35  {
36  // Set verbose output
37  debug = pset.getUntrackedParameter<bool>("debug");
38 
39  // the name of the 1D rec hits collection
40  recHits1DToken_ = consumes<DTRecHitCollection>(pset.getParameter<InputTag>("recHits1DLabel"));
41  recClusToken_ = consumes<DTRecClusterCollection>(pset.getParameter<InputTag>("recClusLabel"));
42 
43  if(debug)
44  cout << "[DTRecSegment2DExtendedProducer] Constructor called" << endl;
45 
46  produces<DTRecSegment2DCollection>();
47 
48  // Get the concrete reconstruction algo from the factory
50 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTCombinatorialExtendedPatternReco * theAlgo
edm::EDGetTokenT< DTRecHitCollection > recHits1DToken_
tuple cout
Definition: gather_cfg.py:121
edm::EDGetTokenT< DTRecClusterCollection > recClusToken_
DTRecSegment2DExtendedProducer::~DTRecSegment2DExtendedProducer ( )
virtual

Destructor.

Definition at line 53 of file DTRecSegment2DExtendedProducer.cc.

References gather_cfg::cout, and debug.

53  {
54  if(debug)
55  cout << "[DTRecSegment2DExtendedProducer] Destructor called" << endl;
56  delete theAlgo;
57 }
DTCombinatorialExtendedPatternReco * theAlgo
tuple cout
Definition: gather_cfg.py:121

Member Function Documentation

void DTRecSegment2DExtendedProducer::produce ( edm::Event event,
const edm::EventSetup setup 
)
virtual

The method which produces the 2D-segments.

Implements edm::EDProducer.

Definition at line 60 of file DTRecSegment2DExtendedProducer.cc.

References edm::OwnVector< T, P >::begin(), filterCSVwithJSON::copy, gather_cfg::cout, debug, edm::OwnVector< T, P >::end(), edm::EventSetup::get(), DTRangeMapAccessor::layersBySuperLayer(), filesave_online::pairs, edm::OwnVector< T, P >::size(), and DTLayerId::superlayerId().

61  {
62  if(debug)
63  cout << "[DTRecSegment2DExtendedProducer] produce called" << endl;
64  // Get the DT Geometry
65  ESHandle<DTGeometry> dtGeom;
66  setup.get<MuonGeometryRecord>().get(dtGeom);
67 
68  theAlgo->setES(setup);
69 
70  // Get the 1D rechits from the event
72  event.getByToken(recHits1DToken_, allHits);
73 
74  // Get the 1D clusters from the event
76  event.getByToken(recClusToken_, dtClusters);
77  theAlgo->setClusters(vector<DTSLRecCluster>(dtClusters->begin(),
78  dtClusters->end()));
79 
80  // Create the pointer to the collection which will store the rechits
81  auto_ptr<DTRecSegment2DCollection> segments(new DTRecSegment2DCollection());
82 
83  // Iterate through all hit collections ordered by LayerId
84  DTRecHitCollection::id_iterator dtLayerIt;
85  DTSuperLayerId oldSlId;
86  for (dtLayerIt = allHits->id_begin(); dtLayerIt != allHits->id_end(); ++dtLayerIt){
87  // The layerId
88  DTLayerId layerId = (*dtLayerIt);
89  const DTSuperLayerId SLId = layerId.superlayerId();
90  if (SLId==oldSlId) continue; // I'm on the same SL as before
91  oldSlId = SLId;
92 
93  if(debug) cout <<"Reconstructing the 2D segments in "<< SLId << endl;
94 
95  const DTSuperLayer* sl = dtGeom->superLayer(SLId);
96 
97  // Get all the rec hit in the same superLayer in which layerId relies
99  allHits->get(DTRangeMapAccessor::layersBySuperLayer(SLId));
100 
101  // Fill the vector with the 1D RecHit
102  vector<DTRecHit1DPair> pairs(range.first,range.second);
103 
104  if(debug) cout << "Number of 1D-RecHit pairs " << pairs.size() << endl;
105 
106  //if(debug) cout << "Start the 2D-segments Reco "<< endl;
108  if(debug) {
109  cout << "Number of Reconstructed segments: " << segs.size() << endl;
110  copy(segs.begin(), segs.end(),
111  ostream_iterator<DTSLRecSegment2D>(cout, "\n"));
112  }
113 
114  if (segs.size() > 0 )
115  segments->put(SLId, segs.begin(),segs.end());
116  }
117  event.put(segments);
118 }
void setClusters(const std::vector< DTSLRecCluster > &clusters)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
static std::pair< DTLayerId, DTSuperLayerIdComparator > layersBySuperLayer(DTSuperLayerId slId)
Access by SL objects written into a RangeMap by layer.
DTCombinatorialExtendedPatternReco * theAlgo
virtual edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
this function is called in the producer
size_type size() const
Definition: OwnVector.h:254
edm::EDGetTokenT< DTRecHitCollection > recHits1DToken_
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
virtual void setES(const edm::EventSetup &setup)
edm::RangeMap< DTSuperLayerId, edm::OwnVector< DTSLRecSegment2D > > DTRecSegment2DCollection
iterator begin()
Definition: OwnVector.h:234
iterator end()
Definition: OwnVector.h:239
const T & get() const
Definition: EventSetup.h:55
tuple cout
Definition: gather_cfg.py:121
list pairs
sort tag files by run number
edm::EDGetTokenT< DTRecClusterCollection > recClusToken_

Member Data Documentation

bool DTRecSegment2DExtendedProducer::debug
private

Definition at line 51 of file DTRecSegment2DExtendedProducer.h.

edm::EDGetTokenT<DTRecClusterCollection> DTRecSegment2DExtendedProducer::recClusToken_
private

Definition at line 58 of file DTRecSegment2DExtendedProducer.h.

edm::EDGetTokenT<DTRecHitCollection> DTRecSegment2DExtendedProducer::recHits1DToken_
private

Definition at line 57 of file DTRecSegment2DExtendedProducer.h.

DTCombinatorialExtendedPatternReco* DTRecSegment2DExtendedProducer::theAlgo
private

Definition at line 54 of file DTRecSegment2DExtendedProducer.h.