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