CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DeDxEstimatorProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DeDxEstimatorProducer
4 // Class: DeDxEstimatorProducer
5 //
13 //
14 // Original Author: andrea
15 // Created: Thu May 31 14:09:02 CEST 2007
16 // Code Updates: loic Quertenmont (querten)
17 // Created: Thu May 10 14:09:02 CEST 2008
18 // $Id: DeDxEstimatorProducer.cc,v 1.30 2010/12/16 17:45:03 innocent Exp $
19 //
20 //
21 
22 
23 // system include files
24 #include <memory>
26 
28 //#include "RecoTracker/DeDx/interface/DeDxTools.h"
34 
39 
43 
46 
47 
48 using namespace reco;
49 using namespace std;
50 using namespace edm;
51 
53 {
54 
55  produces<ValueMap<DeDxData> >();
56 
57 
58  string estimatorName = iConfig.getParameter<string>("estimator");
59  if(estimatorName == "median") m_estimator = new MedianDeDxEstimator(-2.);
60  if(estimatorName == "generic") m_estimator = new GenericAverageDeDxEstimator (iConfig.getParameter<double>("exponent"));
61  if(estimatorName == "truncated") m_estimator = new TruncatedAverageDeDxEstimator(iConfig.getParameter<double>("fraction"));
62  if(estimatorName == "unbinnedFit") m_estimator = new UnbinnedFitDeDxEstimator();
63 
64  MaxNrStrips = iConfig.getUntrackedParameter<unsigned>("maxNrStrips" , 255);
65  MinTrackHits = iConfig.getUntrackedParameter<unsigned>("MinTrackHits" , 4);
66 
67  m_tracksTag = iConfig.getParameter<edm::InputTag>("tracks");
68  m_trajTrackAssociationTag = iConfig.getParameter<edm::InputTag>("trajectoryTrackAssociation");
69 
70  usePixel = iConfig.getParameter<bool>("UsePixel");
71  useStrip = iConfig.getParameter<bool>("UseStrip");
72  MeVperADCPixel = iConfig.getParameter<double>("MeVperADCPixel");
73  MeVperADCStrip = iConfig.getParameter<double>("MeVperADCStrip");
74 
75  shapetest = iConfig.getParameter<bool>("ShapeTest");
76  useCalibration = iConfig.getParameter<bool>("UseCalibration");
77  m_calibrationPath = iConfig.getParameter<string>("calibrationPath");
78 
79  if(!usePixel && !useStrip)
80  edm::LogWarning("DeDxHitsProducer") << "Pixel Hits AND Strip Hits will not be used to estimate dEdx --> BUG, Please Update the config file";
81 }
82 
83 
85 {
86  delete m_estimator;
87 }
88 
89 
90 
91 // ------------ method called once each job just before starting event loop ------------
93 {
94  if(MODsColl.size()!=0)return;
95 
96 
98  iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
99 
100  vector<GeomDet*> Det = tkGeom->dets();
101  for(unsigned int i=0;i<Det.size();i++){
102  DetId Detid = Det[i]->geographicalId();
103 
104  StripGeomDetUnit* StripDetUnit = dynamic_cast<StripGeomDetUnit*> (Det[i]);
105  PixelGeomDetUnit* PixelDetUnit = dynamic_cast<PixelGeomDetUnit*> (Det[i]);
106 
107  double Thick=-1, Dist=-1, Norma=-1;
108  if(StripDetUnit){
109  Dist = StripDetUnit->position().mag();
110  Thick = StripDetUnit->surface().bounds().thickness();
111  Norma = MeVperADCStrip/Thick;
112  }else if(PixelDetUnit){
113  Dist = PixelDetUnit->position().mag();
114  Thick = PixelDetUnit->surface().bounds().thickness();
115  Norma = MeVperADCPixel/Thick;
116  }
117 
118  stModInfo* MOD = new stModInfo;
119  MOD->DetId = Detid.rawId();
120  MOD->Thickness = Thick;
121  MOD->Distance = Dist;
122  MOD->Normalization = Norma;
123  MOD->Gain = 1;
124  MODsColl[MOD->DetId] = MOD;
125  }
126 
127  MakeCalibrationMap();
128 }
129 
130 // ------------ method called once each job just after ending the event loop ------------
132  MODsColl.clear();
133 }
134 
135 
136 
137 
139 {
140  auto_ptr<ValueMap<DeDxData> > trackDeDxEstimateAssociation(new ValueMap<DeDxData> );
141  ValueMap<DeDxData>::Filler filler(*trackDeDxEstimateAssociation);
142 
143  Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
144  iEvent.getByLabel(m_trajTrackAssociationTag, trajTrackAssociationHandle);
145  const TrajTrackAssociationCollection & TrajToTrackMap = *trajTrackAssociationHandle.product();
146 
147  edm::Handle<reco::TrackCollection> trackCollectionHandle;
148  iEvent.getByLabel(m_tracksTag,trackCollectionHandle);
149 
150  size_t n = TrajToTrackMap.size();
151  std::vector<DeDxData> dedxEstimate(n);
152 
153  //assume trajectory collection size is equal to track collection size and that order is kept
154  int j=0;
155  for(TrajTrackAssociationCollection::const_iterator cit=TrajToTrackMap.begin(); cit!=TrajToTrackMap.end(); cit++,j++){
156 
157  const edm::Ref<std::vector<Trajectory> > traj = cit->key;
158  const reco::TrackRef track = cit->val;
159 
160  DeDxHitCollection dedxHits;
161  vector<DeDxTools::RawHits> hits;
162 // DeDxTools::trajectoryRawHits(traj, hits, usePixel, useStrip);
163 
164  const vector<TrajectoryMeasurement> & measurements = traj->measurements();
165  for(vector<TrajectoryMeasurement>::const_iterator it = measurements.begin(); it!=measurements.end(); it++){
166  TrajectoryStateOnSurface trajState=it->updatedState();
167  if( !trajState.isValid()) continue;
168 
169  const TrackingRecHit * recHit=(*it->recHit()).hit();
170  LocalVector trackDirection = trajState.localDirection();
171  double cosine = trackDirection.z()/trackDirection.mag();
172 
173  if(const SiStripMatchedRecHit2D* matchedHit=dynamic_cast<const SiStripMatchedRecHit2D*>(recHit)){
174  if(!useStrip) continue;
175  DeDxTools::RawHits mono,stereo;
176  mono.trajectoryMeasurement = &(*it);
177  stereo.trajectoryMeasurement = &(*it);
178  mono.angleCosine = cosine;
179  stereo.angleCosine = cosine;
180 
181  mono.charge = getCharge((matchedHit->monoHit()->cluster()).get(),mono.NSaturating);
182  stereo.charge = getCharge((matchedHit->stereoHit()->cluster()).get(),stereo.NSaturating);
183 
184  mono.detId= matchedHit->monoHit()->geographicalId();
185  stereo.detId= matchedHit->stereoHit()->geographicalId();
186 
187  if(shapetest && !(DeDxTools::shapeSelection(((matchedHit->stereoHit()->cluster()).get())->amplitudes()))) hits.push_back(stereo);
188  if(shapetest && !(DeDxTools::shapeSelection(((matchedHit->monoHit ()->cluster()).get())->amplitudes()))) hits.push_back(mono);
189  }else if(const ProjectedSiStripRecHit2D* projectedHit=dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit)) {
190  if(!useStrip) continue;
191  const SiStripRecHit2D* singleHit=&(projectedHit->originalHit());
192  DeDxTools::RawHits mono;
193 
194  mono.trajectoryMeasurement = &(*it);
195  mono.angleCosine = cosine;
196  mono.charge = getCharge((singleHit->cluster()).get(),mono.NSaturating);
197  mono.detId= singleHit->geographicalId();
198  if(shapetest && !(DeDxTools::shapeSelection(((singleHit->cluster()).get())->amplitudes()))) continue;
199  hits.push_back(mono);
200  }else if(const SiStripRecHit2D* singleHit=dynamic_cast<const SiStripRecHit2D*>(recHit)){
201  if(!useStrip) continue;
202  DeDxTools::RawHits mono;
203 
204  mono.trajectoryMeasurement = &(*it);
205  mono.angleCosine = cosine;
206  mono.charge = getCharge((singleHit->cluster()).get(),mono.NSaturating);
207  mono.detId= singleHit->geographicalId();
208  if(shapetest && !(DeDxTools::shapeSelection(((singleHit->cluster()).get())->amplitudes()))) continue;
209  hits.push_back(mono);
210  }else if(const SiStripRecHit1D* single1DHit=dynamic_cast<const SiStripRecHit1D*>(recHit)){
211  if(!useStrip) continue;
212  DeDxTools::RawHits mono;
213 
214  mono.trajectoryMeasurement = &(*it);
215  mono.angleCosine = cosine;
216  mono.charge = getCharge((single1DHit->cluster()).get(),mono.NSaturating);
217  mono.detId= single1DHit->geographicalId();
218  if(shapetest && !(DeDxTools::shapeSelection(((single1DHit->cluster()).get())->amplitudes()))) continue;
219  hits.push_back(mono);
220  }else if(const SiPixelRecHit* pixelHit=dynamic_cast<const SiPixelRecHit*>(recHit)){
221  if(!usePixel) continue;
222 
223  DeDxTools::RawHits pixel;
224 
225  pixel.trajectoryMeasurement = &(*it);
226  pixel.angleCosine = cosine;
227  pixel.charge = pixelHit->cluster()->charge();
228  pixel.NSaturating=-1;
229  pixel.detId= pixelHit->geographicalId();
230  hits.push_back(pixel);
231  }
232  }
234 
235  int NClusterSaturating = 0;
236  for(size_t i=0; i < hits.size(); i++)
237  {
238  stModInfo* MOD = MODsColl[hits[i].detId];
239  float pathLen = MOD->Thickness/std::abs(hits[i].angleCosine);
240  float charge = MOD->Normalization*hits[i].charge*std::abs(hits[i].angleCosine);
241  dedxHits.push_back( DeDxHit( charge, MOD->Distance, pathLen, hits[i].detId) );
242 
243  if(hits[i].NSaturating>0)NClusterSaturating++;
244  }
245 
246  sort(dedxHits.begin(),dedxHits.end(),less<DeDxHit>());
247  std::pair<float,float> val_and_error = m_estimator->dedx(dedxHits);
248 
249  //WARNING: Since the dEdX Error is not properly computed for the moment
250  //It was decided to store the number of saturating cluster in that dataformat
251  val_and_error.second = NClusterSaturating;
252 
253  dedxEstimate[j] = DeDxData(val_and_error.first, val_and_error.second, dedxHits.size() );
254  }
255  filler.insert(trackCollectionHandle, dedxEstimate.begin(), dedxEstimate.end());
256 
257  // really fill the association map
258  filler.fill();
259  // put into the event
260  iEvent.put(trackDeDxEstimateAssociation);
261 }
262 
263 
264 
266 {
267  if(!useCalibration)return;
268 
269  TChain* t1 = new TChain("SiStripCalib/APVGain");
270  t1->Add(m_calibrationPath.c_str());
271 
272  unsigned int tree_DetId;
273  unsigned char tree_APVId;
274  double tree_Gain;
275 
276  t1->SetBranchAddress("DetId" ,&tree_DetId );
277  t1->SetBranchAddress("APVId" ,&tree_APVId );
278  t1->SetBranchAddress("Gain" ,&tree_Gain );
279 
280 
281 
282  for (unsigned int ientry = 0; ientry < t1->GetEntries(); ientry++) {
283  t1->GetEntry(ientry);
284  stModInfo* MOD = MODsColl[tree_DetId];
285  MOD->Gain = tree_Gain;
286  }
287 
288  delete t1;
289 
290 }
291 
292 
293 int DeDxEstimatorProducer::getCharge(const SiStripCluster* Cluster, int& Saturating_Strips)
294 {
295  const vector<uint8_t>& Ampls = Cluster->amplitudes();
296  uint32_t DetId = Cluster->geographicalId();
297 
298 // float G=1.0f;
299 
300  int toReturn = 0;
301  Saturating_Strips = 0;
302  for(unsigned int i=0;i<Ampls.size();i++){
303  int CalibratedCharge = Ampls[i];
304 
305  if(useCalibration){
306  stModInfo* MOD = MODsColl[DetId];
307 // G = MOD->Gain;
308  CalibratedCharge = (int)(CalibratedCharge / MOD->Gain );
309  if(CalibratedCharge>=1024){
310  CalibratedCharge = 255;
311  }else if(CalibratedCharge>=255){
312  CalibratedCharge = 254;
313  }
314  }
315 
316  toReturn+=CalibratedCharge;
317  if(CalibratedCharge>=254)Saturating_Strips++;
318  }
319 
320 // printf("Charge = %i --> %i (Gain=%f)\n", accumulate(Ampls.begin(), Ampls.end(), 0), toReturn, G);
321 
322 
323  return toReturn;
324 }
325 
326 
327 
328 
329 //define this as a plug-in
T getParameter(std::string const &) const
DeDxEstimatorProducer(const edm::ParameterSet &)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const_iterator end() const
last iterator over the map (read only)
LocalVector localDirection() const
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::vector< DeDxHit > DeDxHitCollection
Definition: DeDxHit.h:49
#define abs(x)
Definition: mlp_lapack.h:159
double charge(const std::vector< uint8_t > &Ampls)
uint32_t geographicalId() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
virtual float thickness() const =0
int iEvent
Definition: GenABIO.cc:243
T mag() const
Definition: PV3DBase.h:61
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
double angleCosine
Definition: DeDxTools.h:14
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
int getCharge(const SiStripCluster *Cluster, int &Saturating_Strips)
T z() const
Definition: PV3DBase.h:58
int j
Definition: DBlmapReader.cc:9
virtual void beginRun(edm::Run &run, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
Definition: DetId.h:20
const Bounds & bounds() const
Definition: BoundSurface.h:89
bool shapeSelection(const std::vector< uint8_t > &ampls)
Definition: DeDxTools.cc:141
size_type size() const
map size
const TrajectoryMeasurement * trajectoryMeasurement
Definition: DeDxTools.h:16
const T & get() const
Definition: EventSetup.h:55
ClusterRef const & cluster() const
key_type key() const
Accessor for product key.
Definition: Ref.h:265
virtual void produce(edm::Event &, const edm::EventSetup &)
T const * product() const
Definition: Handle.h:74
const_iterator begin() const
first iterator over the map (read only)
DetId geographicalId() const
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const std::vector< uint8_t > & amplitudes() const
Definition: Run.h:32
Our base class.
Definition: SiPixelRecHit.h:27