CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoLocalMuon/DTSegment/src/DTClusterer.cc

Go to the documentation of this file.
00001 /******* \class DTClusterer *******
00002  *
00003  * Description:
00004  *  
00005  *  detailed description
00006  *
00007  * \author : Stefano Lacaprara - INFN LNL <stefano.lacaprara@pd.infn.it>
00008  * $date   : 17/04/2008 14:55:53 CEST $
00009  *
00010  * Modification:
00011  *
00012  *********************************/
00013 
00014 /* This Class Header */
00015 #include "RecoLocalMuon/DTSegment/src/DTClusterer.h"
00016 
00017 /* Collaborating Class Header */
00018 #include "FWCore/Framework/interface/Event.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00021 using namespace edm;
00022 
00023 #include "Geometry/DTGeometry/interface/DTLayer.h"
00024 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00025 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00026 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00027 
00028 #include "DataFormats/DTRecHit/interface/DTRecHit1DPair.h"
00029 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00030 #include "DataFormats/DTRecHit/interface/DTRecClusterCollection.h"
00031 #include "DataFormats/DTRecHit/interface/DTRangeMapAccessor.h"
00032 
00033 /* C++ Headers */
00034 #include <iostream>
00035 using namespace std;
00036 
00037 /* ====================================================================== */
00038 
00039 /* Constructor */ 
00040 DTClusterer::DTClusterer(const edm::ParameterSet& pset) {
00041   // Set verbose output
00042   debug = pset.getUntrackedParameter<bool>("debug"); 
00043 
00044   // the name of the 1D rec hits collection
00045   theRecHits1DLabel = pset.getParameter<InputTag>("recHits1DLabel");
00046   // min number of hits to build a cluster
00047   theMinHits = pset.getParameter<unsigned int>("minHits");
00048   // min number of hits to build a cluster
00049   theMinLayers = pset.getParameter<unsigned int>("minLayers");
00050 
00051   if(debug)
00052     cout << "[DTClusterer] Constructor called" << endl;
00053 
00054   produces<DTRecClusterCollection>();
00055 }
00056 
00057 /* Destructor */ 
00058 DTClusterer::~DTClusterer() {
00059 }
00060 
00061 /* Operations */ 
00062 void DTClusterer::produce(edm::Event& event, const edm::EventSetup& setup) {
00063   if(debug)
00064     cout << "[DTClusterer] produce called" << endl;
00065   // Get the DT Geometry
00066   ESHandle<DTGeometry> dtGeom;
00067   setup.get<MuonGeometryRecord>().get(dtGeom);
00068 
00069   // Get the 1D rechits from the event
00070   Handle<DTRecHitCollection> allHits; 
00071   event.getByLabel(theRecHits1DLabel, allHits);
00072 
00073   // Create the pointer to the collection which will store the rechits
00074   auto_ptr<DTRecClusterCollection> clusters(new DTRecClusterCollection());
00075 
00076   // Iterate through all hit collections ordered by LayerId
00077   DTRecHitCollection::id_iterator dtLayerIt;
00078   DTSuperLayerId oldSlId;
00079   for (dtLayerIt = allHits->id_begin(); dtLayerIt != allHits->id_end(); ++dtLayerIt){
00080     // The layerId
00081     DTLayerId layerId = (*dtLayerIt);
00082     const DTSuperLayerId SLId = layerId.superlayerId();
00083     if (SLId==oldSlId) continue; // I'm on the same SL as before
00084     oldSlId = SLId;
00085 
00086     if(debug) cout <<"Reconstructing the clusters in "<< SLId << endl;
00087 
00088     const DTSuperLayer* sl = dtGeom->superLayer(SLId);
00089 
00090     // Get all the rec hit in the same superLayer in which layerId relies 
00091     DTRecHitCollection::range range =
00092       allHits->get(DTRangeMapAccessor::layersBySuperLayer(SLId));
00093 
00094     // Fill the vector with the 1D RecHit
00095     vector<DTRecHit1DPair> pairs(range.first,range.second);
00096     if(debug) cout << "Number of 1D-RecHit pairs " << pairs.size() << endl;
00097     vector<DTSLRecCluster> clus=buildClusters(sl, pairs);
00098     if(debug) cout << "Number of clusters build " << clus.size() << endl;
00099     if (clus.size() > 0 )
00100       clusters->put(sl->id(), clus.begin(), clus.end());
00101   }
00102   event.put(clusters);
00103 }
00104 
00105 vector<DTSLRecCluster> DTClusterer::buildClusters(const DTSuperLayer* sl,
00106                                                   vector<DTRecHit1DPair>& pairs) {
00107   // create a vector of hits with wire position in SL frame
00108   vector<pair<float, DTRecHit1DPair> > hits = initHits(sl, pairs);
00109 
00110   vector<DTSLRecCluster> result;
00111   // loop over pairs
00112   vector<DTRecHit1DPair> adiacentPairs;
00113   float lastPos=hits.front().first;
00114   const float cellWidth=4.2; // cm
00115   float sum=0.;
00116   float sum2=0.;
00117 
00118   for (vector<pair<float, DTRecHit1DPair> >::const_iterator hit=hits.begin(); hit!=hits.end();
00119        ++hit) {
00120     if(debug) cout << "Hit: " << (*hit).first << " lastPos: " << lastPos << endl;
00121     // start from first hits
00122     // two cells are adiacente if their position is closer than cell width
00123     if(abs((*hit).first-lastPos)>cellWidth ) {
00124       if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
00125         // if not, build the cluster with so far collection hits and go on
00126         float mean=sum/adiacentPairs.size();
00127         float err2=sum2/adiacentPairs.size()-mean*mean;
00128         DTSLRecCluster cluster(sl->id(),
00129                                LocalPoint(mean,0.,0.),
00130                                LocalError(err2,0.,0.),
00131                                adiacentPairs);
00132         if(debug) cout << "Cluster " << cluster << endl;
00133         result.push_back(cluster);
00134       }
00135       // clean the vector
00136       adiacentPairs.clear();
00137       sum=0.;
00138       sum2=0.;
00139     }
00140     // if adiacente, add them to a vector
00141     adiacentPairs.push_back((*hit).second);
00142     if(debug) cout << "adiacentPairs " << adiacentPairs.size() << endl;
00143     sum+=(*hit).first;
00144     sum2+=(*hit).first*(*hit).first;
00145 
00146     lastPos=(*hit).first;
00147   }
00148   // build the last cluster
00149   if (adiacentPairs.size() >= theMinHits && differentLayers(adiacentPairs)>=theMinLayers) {
00150     float mean=sum/adiacentPairs.size();
00151     float err2=sum2/adiacentPairs.size()-mean*mean;
00152     DTSLRecCluster cluster(sl->id(),
00153                            LocalPoint(mean,0.,0.),
00154                            LocalError(err2,0.,0.),
00155                            adiacentPairs);
00156     if(debug) cout << "Cluster " << cluster << endl;
00157     result.push_back(cluster);
00158   }
00159 
00160   return result;
00161 }
00162 
00163 vector<pair<float, DTRecHit1DPair> > DTClusterer::initHits(const DTSuperLayer* sl,
00164                                                            vector<DTRecHit1DPair>& pairs) {
00165   vector<pair<float, DTRecHit1DPair> > result;
00166   for (vector<DTRecHit1DPair>::const_iterator pair=pairs.begin();
00167        pair!=pairs.end(); ++pair ) {
00168 
00169     // get wire
00170     DTWireId wid = (*pair).wireId();
00171     // get Layer
00172     const DTLayer* lay= sl->layer(wid.layer());
00173     // get wire position in SL (only x)
00174     LocalPoint
00175       posInLayer(lay->specificTopology().wirePosition(wid.wire()),0.,0.);
00176     LocalPoint posInSL = sl->toLocal(lay->toGlobal(posInLayer));
00177     // put the pair into result
00178     result.push_back(make_pair(posInSL.x(), *pair));
00179   }
00180   // sorted by x
00181   sort(result.begin(), result.end(), sortClusterByX());
00182 
00183   return result;
00184 }
00185 
00186 unsigned int DTClusterer::differentLayers(vector<DTRecHit1DPair>& hits) {
00187 // Count the number of different layers
00188   int layers=0;
00189   unsigned int result=0;
00190   for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin(); hit!=hits.end();
00191        ++hit) {
00192     int pos = (1 << ((*hit).wireId().layer() -1 ));
00193     if (! (pos & layers)) {
00194       result++;
00195       layers |= pos;
00196     }
00197   }
00198   return result;
00199 }