CMS 3D CMS Logo

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

#include <DTClusterer.h>

Inheritance diagram for DTClusterer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Classes

struct  sortClusterByX
 

Public Member Functions

 DTClusterer (const edm::ParameterSet &)
 
virtual void produce (edm::Event &event, const edm::EventSetup &setup)
 
virtual ~DTClusterer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

std::vector< DTSLRecClusterbuildClusters (const DTSuperLayer *sl, std::vector< DTRecHit1DPair > &pairs)
 
unsigned int differentLayers (std::vector< DTRecHit1DPair > &hits)
 
std::vector< std::pair< float,
DTRecHit1DPair > > 
initHits (const DTSuperLayer *sl, std::vector< DTRecHit1DPair > &pairs)
 

Private Attributes

bool debug
 
unsigned int theMinHits
 
unsigned int theMinLayers
 
edm::InputTag theRecHits1DLabel
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- 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)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Description:

detailed description

Author
: Stefano Lacaprara - INFN LNL stefa.nosp@m.no.l.nosp@m.acapr.nosp@m.ara@.nosp@m.pd.in.nosp@m.fn.i.nosp@m.t $date : 17/04/2008 14:56:40 CEST $

Modification:

Definition at line 40 of file DTClusterer.h.

Constructor & Destructor Documentation

DTClusterer::DTClusterer ( const edm::ParameterSet pset)

Definition at line 40 of file DTClusterer.cc.

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

40  {
41  // Set verbose output
42  debug = pset.getUntrackedParameter<bool>("debug");
43 
44  // the name of the 1D rec hits collection
45  theRecHits1DLabel = pset.getParameter<InputTag>("recHits1DLabel");
46  // min number of hits to build a cluster
47  theMinHits = pset.getParameter<unsigned int>("minHits");
48  // min number of hits to build a cluster
49  theMinLayers = pset.getParameter<unsigned int>("minLayers");
50 
51  if(debug)
52  cout << "[DTClusterer] Constructor called" << endl;
53 
54  produces<DTRecClusterCollection>();
55 }
unsigned int theMinLayers
Definition: DTClusterer.h:77
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned int theMinHits
Definition: DTClusterer.h:76
edm::InputTag theRecHits1DLabel
Definition: DTClusterer.h:78
tuple cout
Definition: gather_cfg.py:41
DTClusterer::~DTClusterer ( )
virtual

Definition at line 58 of file DTClusterer.cc.

58  {
59 }

Member Function Documentation

vector< DTSLRecCluster > DTClusterer::buildClusters ( const DTSuperLayer sl,
std::vector< DTRecHit1DPair > &  pairs 
)
private

Definition at line 105 of file DTClusterer.cc.

References abs, gather_cfg::cout, debug, DTSuperLayer::id(), plotscripts::mean(), and query::result.

106  {
107  // create a vector of hits with wire position in SL frame
108  vector<pair<float, DTRecHit1DPair> > hits = initHits(sl, pairs);
109 
110  vector<DTSLRecCluster> result;
111  // loop over pairs
112  vector<DTRecHit1DPair> adiacentPairs;
113  float lastPos=hits.front().first;
114  const float cellWidth=4.2; // cm
115  float sum=0.;
116  float sum2=0.;
117 
118  for (vector<pair<float, DTRecHit1DPair> >::const_iterator hit=hits.begin(); hit!=hits.end();
119  ++hit) {
120  if(debug) cout << "Hit: " << (*hit).first << " lastPos: " << lastPos << endl;
121  // start from first hits
122  // two cells are adiacente if their position is closer than cell width
123  if(abs((*hit).first-lastPos)>cellWidth ) {
124  if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
125  // if not, build the cluster with so far collection hits and go on
126  float mean=sum/adiacentPairs.size();
127  float err2=sum2/adiacentPairs.size()-mean*mean;
128  DTSLRecCluster cluster(sl->id(),
129  LocalPoint(mean,0.,0.),
130  LocalError(err2,0.,0.),
131  adiacentPairs);
132  if(debug) cout << "Cluster " << cluster << endl;
133  result.push_back(cluster);
134  }
135  // clean the vector
136  adiacentPairs.clear();
137  sum=0.;
138  sum2=0.;
139  }
140  // if adiacente, add them to a vector
141  adiacentPairs.push_back((*hit).second);
142  if(debug) cout << "adiacentPairs " << adiacentPairs.size() << endl;
143  sum+=(*hit).first;
144  sum2+=(*hit).first*(*hit).first;
145 
146  lastPos=(*hit).first;
147  }
148  // build the last cluster
149  if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
150  float mean=sum/adiacentPairs.size();
151  float err2=sum2/adiacentPairs.size()-mean*mean;
152  DTSLRecCluster cluster(sl->id(),
153  LocalPoint(mean,0.,0.),
154  LocalError(err2,0.,0.),
155  adiacentPairs);
156  if(debug) cout << "Cluster " << cluster << endl;
157  result.push_back(cluster);
158  }
159 
160  return result;
161 }
unsigned int theMinLayers
Definition: DTClusterer.h:77
unsigned int differentLayers(std::vector< DTRecHit1DPair > &hits)
Definition: DTClusterer.cc:186
unsigned int theMinHits
Definition: DTClusterer.h:76
#define abs(x)
Definition: mlp_lapack.h:159
DTSuperLayerId id() const
Return the DetId of this SL.
Definition: DTSuperLayer.cc:38
tuple result
Definition: query.py:137
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
std::vector< std::pair< float, DTRecHit1DPair > > initHits(const DTSuperLayer *sl, std::vector< DTRecHit1DPair > &pairs)
Definition: DTClusterer.cc:163
tuple cout
Definition: gather_cfg.py:41
unsigned int DTClusterer::differentLayers ( std::vector< DTRecHit1DPair > &  hits)
private

Definition at line 186 of file DTClusterer.cc.

References pos, and query::result.

186  {
187 // Count the number of different layers
188  int layers=0;
189  unsigned int result=0;
190  for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin(); hit!=hits.end();
191  ++hit) {
192  int pos = (1 << ((*hit).wireId().layer() -1 ));
193  if (! (pos & layers)) {
194  result++;
195  layers |= pos;
196  }
197  }
198  return result;
199 }
tuple result
Definition: query.py:137
vector< pair< float, DTRecHit1DPair > > DTClusterer::initHits ( const DTSuperLayer sl,
std::vector< DTRecHit1DPair > &  pairs 
)
private

Definition at line 163 of file DTClusterer.cc.

References DTLayerId::layer(), DTSuperLayer::layer(), query::result, python.multivaluedict::sort(), DTLayer::specificTopology(), GeomDet::toGlobal(), GeomDet::toLocal(), DTWireId::wire(), DTTopology::wirePosition(), and PV3DBase< T, PVType, FrameType >::x().

164  {
165  vector<pair<float, DTRecHit1DPair> > result;
166  for (vector<DTRecHit1DPair>::const_iterator pair=pairs.begin();
167  pair!=pairs.end(); ++pair ) {
168 
169  // get wire
170  DTWireId wid = (*pair).wireId();
171  // get Layer
172  const DTLayer* lay= sl->layer(wid.layer());
173  // get wire position in SL (only x)
174  LocalPoint
175  posInLayer(lay->specificTopology().wirePosition(wid.wire()),0.,0.);
176  LocalPoint posInSL = sl->toLocal(lay->toGlobal(posInLayer));
177  // put the pair into result
178  result.push_back(make_pair(posInSL.x(), *pair));
179  }
180  // sorted by x
181  sort(result.begin(), result.end(), sortClusterByX());
182 
183  return result;
184 }
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:73
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:64
int layer() const
Return the layer number.
Definition: DTLayerId.h:55
const DTTopology & specificTopology() const
Definition: DTLayer.cc:44
tuple result
Definition: query.py:137
const DTLayer * layer(DTLayerId id) const
Return the layer corresponding to the given id.
Definition: DTSuperLayer.cc:70
int wire() const
Return the wire number.
Definition: DTWireId.h:58
T x() const
Definition: PV3DBase.h:56
void DTClusterer::produce ( edm::Event event,
const edm::EventSetup setup 
)
virtual

Implements edm::EDProducer.

Definition at line 62 of file DTClusterer.cc.

References gather_cfg::cout, debug, edm::EventSetup::get(), DTSuperLayer::id(), DTRangeMapAccessor::layersBySuperLayer(), and DTLayerId::superlayerId().

62  {
63  if(debug)
64  cout << "[DTClusterer] produce called" << endl;
65  // Get the DT Geometry
66  ESHandle<DTGeometry> dtGeom;
67  setup.get<MuonGeometryRecord>().get(dtGeom);
68 
69  // Get the 1D rechits from the event
71  event.getByLabel(theRecHits1DLabel, allHits);
72 
73  // Create the pointer to the collection which will store the rechits
74  auto_ptr<DTRecClusterCollection> clusters(new DTRecClusterCollection());
75 
76  // Iterate through all hit collections ordered by LayerId
78  DTSuperLayerId oldSlId;
79  for (dtLayerIt = allHits->id_begin(); dtLayerIt != allHits->id_end(); ++dtLayerIt){
80  // The layerId
81  DTLayerId layerId = (*dtLayerIt);
82  const DTSuperLayerId SLId = layerId.superlayerId();
83  if (SLId==oldSlId) continue; // I'm on the same SL as before
84  oldSlId = SLId;
85 
86  if(debug) cout <<"Reconstructing the clusters in "<< SLId << endl;
87 
88  const DTSuperLayer* sl = dtGeom->superLayer(SLId);
89 
90  // Get all the rec hit in the same superLayer in which layerId relies
92  allHits->get(DTRangeMapAccessor::layersBySuperLayer(SLId));
93 
94  // Fill the vector with the 1D RecHit
95  vector<DTRecHit1DPair> pairs(range.first,range.second);
96  if(debug) cout << "Number of 1D-RecHit pairs " << pairs.size() << endl;
97  vector<DTSLRecCluster> clus=buildClusters(sl, pairs);
98  if(debug) cout << "Number of clusters build " << clus.size() << endl;
99  if (clus.size() > 0 )
100  clusters->put(sl->id(), clus.begin(), clus.end());
101  }
102  event.put(clusters);
103 }
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
edm::RangeMap< DTSuperLayerId, edm::OwnVector< DTSLRecCluster > > DTRecClusterCollection
static std::pair< DTLayerId, DTSuperLayerIdComparator > layersBySuperLayer(DTSuperLayerId slId)
Access by SL objects written into a RangeMap by layer.
identifier iterator
Definition: RangeMap.h:138
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:61
DTSuperLayerId id() const
Return the DetId of this SL.
Definition: DTSuperLayer.cc:38
std::vector< DTSLRecCluster > buildClusters(const DTSuperLayer *sl, std::vector< DTRecHit1DPair > &pairs)
Definition: DTClusterer.cc:105
edm::InputTag theRecHits1DLabel
Definition: DTClusterer.h:78
const T & get() const
Definition: EventSetup.h:55
tuple cout
Definition: gather_cfg.py:41

Member Data Documentation

bool DTClusterer::debug
private

Definition at line 74 of file DTClusterer.h.

unsigned int DTClusterer::theMinHits
private

Definition at line 76 of file DTClusterer.h.

unsigned int DTClusterer::theMinLayers
private

Definition at line 77 of file DTClusterer.h.

edm::InputTag DTClusterer::theRecHits1DLabel
private

Definition at line 78 of file DTClusterer.h.