CMS 3D CMS Logo

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

List of all members.

Classes

struct  sortClusterByX

Public Member Functions

 DTClusterer (const edm::ParameterSet &)
virtual void produce (edm::Event &event, const edm::EventSetup &setup)
virtual ~DTClusterer ()

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

Detailed Description

Description:

detailed description

Author:
: Stefano Lacaprara - INFN LNL <stefano.lacaprara@pd.infn.it> $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().

                                                    {
  // Set verbose output
  debug = pset.getUntrackedParameter<bool>("debug"); 

  // the name of the 1D rec hits collection
  theRecHits1DLabel = pset.getParameter<InputTag>("recHits1DLabel");
  // min number of hits to build a cluster
  theMinHits = pset.getParameter<unsigned int>("minHits");
  // min number of hits to build a cluster
  theMinLayers = pset.getParameter<unsigned int>("minLayers");

  if(debug)
    cout << "[DTClusterer] Constructor called" << endl;

  produces<DTRecClusterCollection>();
}
DTClusterer::~DTClusterer ( ) [virtual]

Definition at line 58 of file DTClusterer.cc.

                          {
}

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.

                                                                                 {
  // create a vector of hits with wire position in SL frame
  vector<pair<float, DTRecHit1DPair> > hits = initHits(sl, pairs);

  vector<DTSLRecCluster> result;
  // loop over pairs
  vector<DTRecHit1DPair> adiacentPairs;
  float lastPos=hits.front().first;
  const float cellWidth=4.2; // cm
  float sum=0.;
  float sum2=0.;

  for (vector<pair<float, DTRecHit1DPair> >::const_iterator hit=hits.begin(); hit!=hits.end();
       ++hit) {
    if(debug) cout << "Hit: " << (*hit).first << " lastPos: " << lastPos << endl;
    // start from first hits
    // two cells are adiacente if their position is closer than cell width
    if(abs((*hit).first-lastPos)>cellWidth ) {
      if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
        // if not, build the cluster with so far collection hits and go on
        float mean=sum/adiacentPairs.size();
        float err2=sum2/adiacentPairs.size()-mean*mean;
        DTSLRecCluster cluster(sl->id(),
                               LocalPoint(mean,0.,0.),
                               LocalError(err2,0.,0.),
                               adiacentPairs);
        if(debug) cout << "Cluster " << cluster << endl;
        result.push_back(cluster);
      }
      // clean the vector
      adiacentPairs.clear();
      sum=0.;
      sum2=0.;
    }
    // if adiacente, add them to a vector
    adiacentPairs.push_back((*hit).second);
    if(debug) cout << "adiacentPairs " << adiacentPairs.size() << endl;
    sum+=(*hit).first;
    sum2+=(*hit).first*(*hit).first;

    lastPos=(*hit).first;
  }
  // build the last cluster
  if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
    float mean=sum/adiacentPairs.size();
    float err2=sum2/adiacentPairs.size()-mean*mean;
    DTSLRecCluster cluster(sl->id(),
                           LocalPoint(mean,0.,0.),
                           LocalError(err2,0.,0.),
                           adiacentPairs);
    if(debug) cout << "Cluster " << cluster << endl;
    result.push_back(cluster);
  }

  return result;
}
unsigned int DTClusterer::differentLayers ( std::vector< DTRecHit1DPair > &  hits) [private]

Definition at line 186 of file DTClusterer.cc.

References pos, and query::result.

                                                                      {
// Count the number of different layers
  int layers=0;
  unsigned int result=0;
  for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin(); hit!=hits.end();
       ++hit) {
    int pos = (1 << ((*hit).wireId().layer() -1 ));
    if (! (pos & layers)) {
      result++;
      layers |= pos;
    }
  }
  return result;
}
std::vector<std::pair<float, DTRecHit1DPair> > DTClusterer::initHits ( const DTSuperLayer sl,
std::vector< DTRecHit1DPair > &  pairs 
) [private]
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().

                                                                     {
  if(debug)
    cout << "[DTClusterer] produce called" << endl;
  // Get the DT Geometry
  ESHandle<DTGeometry> dtGeom;
  setup.get<MuonGeometryRecord>().get(dtGeom);

  // Get the 1D rechits from the event
  Handle<DTRecHitCollection> allHits; 
  event.getByLabel(theRecHits1DLabel, allHits);

  // Create the pointer to the collection which will store the rechits
  auto_ptr<DTRecClusterCollection> clusters(new DTRecClusterCollection());

  // Iterate through all hit collections ordered by LayerId
  DTRecHitCollection::id_iterator dtLayerIt;
  DTSuperLayerId oldSlId;
  for (dtLayerIt = allHits->id_begin(); dtLayerIt != allHits->id_end(); ++dtLayerIt){
    // The layerId
    DTLayerId layerId = (*dtLayerIt);
    const DTSuperLayerId SLId = layerId.superlayerId();
    if (SLId==oldSlId) continue; // I'm on the same SL as before
    oldSlId = SLId;

    if(debug) cout <<"Reconstructing the clusters in "<< SLId << endl;

    const DTSuperLayer* sl = dtGeom->superLayer(SLId);

    // Get all the rec hit in the same superLayer in which layerId relies 
    DTRecHitCollection::range range =
      allHits->get(DTRangeMapAccessor::layersBySuperLayer(SLId));

    // Fill the vector with the 1D RecHit
    vector<DTRecHit1DPair> pairs(range.first,range.second);
    if(debug) cout << "Number of 1D-RecHit pairs " << pairs.size() << endl;
    vector<DTSLRecCluster> clus=buildClusters(sl, pairs);
    if(debug) cout << "Number of clusters build " << clus.size() << endl;
    if (clus.size() > 0 )
      clusters->put(sl->id(), clus.begin(), clus.end());
  }
  event.put(clusters);
}

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.

Definition at line 78 of file DTClusterer.h.