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