CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ShallowGainCalibration.cc
Go to the documentation of this file.
4 
5 using namespace edm;
6 using namespace reco;
7 using namespace std;
8 
9 
10 
11 
13  : theTracksLabel( iConfig.getParameter<edm::InputTag>("Tracks") ),
14  Suffix ( iConfig.getParameter<std::string>("Suffix") ),
15  Prefix ( iConfig.getParameter<std::string>("Prefix") )
16 {
17  produces <std::vector<int> > ( Prefix + "trackindex" + Suffix );
18  produces <std::vector<unsigned int> > ( Prefix + "rawid" + Suffix );
19  produces <std::vector<float> > ( Prefix + "localdirx" + Suffix );
20  produces <std::vector<float> > ( Prefix + "localdiry" + Suffix );
21  produces <std::vector<float> > ( Prefix + "localdirz" + Suffix );
22  produces <std::vector<unsigned short> > ( Prefix + "firststrip" + Suffix );
23  produces <std::vector<unsigned short> > ( Prefix + "nstrips" + Suffix );
24  produces <std::vector<bool> > ( Prefix + "saturation" + Suffix );
25  produces <std::vector<bool> > ( Prefix + "overlapping" + Suffix );
26  produces <std::vector<bool> > ( Prefix + "farfromedge" + Suffix );
27  produces <std::vector<unsigned int> > ( Prefix + "charge" + Suffix );
28  produces <std::vector<float> > ( Prefix + "path" + Suffix );
29  produces <std::vector<float> > ( Prefix + "chargeoverpath" + Suffix );
30  produces <std::vector<unsigned char> > ( Prefix + "amplitude" + Suffix );
31  produces <std::vector<double> > ( Prefix + "gainused" + Suffix );
32  produces <std::vector<double> > ( Prefix + "gainusedTick" + Suffix );
33 }
34 
37  std::auto_ptr<std::vector<int> > trackindex ( new std::vector<int> );
38  std::auto_ptr<std::vector<unsigned int> > rawid ( new std::vector<unsigned int> );
39  std::auto_ptr<std::vector<float> > localdirx ( new std::vector<float> );
40  std::auto_ptr<std::vector<float> > localdiry ( new std::vector<float> );
41  std::auto_ptr<std::vector<float> > localdirz ( new std::vector<float> );
42  std::auto_ptr<std::vector<unsigned short> > firststrip ( new std::vector<unsigned short> );
43  std::auto_ptr<std::vector<unsigned short> > nstrips ( new std::vector<unsigned short> );
44  std::auto_ptr<std::vector<bool> > saturation ( new std::vector<bool> );
45  std::auto_ptr<std::vector<bool> > overlapping ( new std::vector<bool> );
46  std::auto_ptr<std::vector<bool> > farfromedge ( new std::vector<bool> );
47  std::auto_ptr<std::vector<unsigned int> > charge ( new std::vector<unsigned int> );
48  std::auto_ptr<std::vector<float> > path ( new std::vector<float> );
49  std::auto_ptr<std::vector<float> > chargeoverpath( new std::vector<float> );
50  std::auto_ptr<std::vector<unsigned char> > amplitude ( new std::vector<unsigned char> );
51  std::auto_ptr<std::vector<double> > gainused ( new std::vector<double> );
52  std::auto_ptr<std::vector<double> > gainusedTick ( new std::vector<double> );
53 
54  edm::ESHandle<TrackerGeometry> theTrackerGeometry; iSetup.get<TrackerDigiGeometryRecord>().get( theTrackerGeometry );
55  m_tracker=&(* theTrackerGeometry );
56  edm::ESHandle<SiStripGain> gainHandle; iSetup.get<SiStripGainRcd>().get(gainHandle);
58  edm::Handle<TrajTrackAssociationCollection> associations; iEvent.getByLabel(theTracksLabel, associations);
59 
60  for( TrajTrackAssociationCollection::const_iterator association = associations->begin(); association != associations->end(); association++) {
61  const Trajectory* traj = association->key.get();
62  const reco::Track* track = association->val.get();
63 
64  vector<TrajectoryMeasurement> measurements = traj->measurements();
65  for(vector<TrajectoryMeasurement>::const_iterator measurement_it = measurements.begin(); measurement_it!=measurements.end(); measurement_it++){
66  TrajectoryStateOnSurface trajState = measurement_it->updatedState();
67  if( !trajState.isValid() ) continue;
68 
69  const TrackingRecHit* hit = (*measurement_it->recHit()).hit();
70  const SiStripRecHit1D* sistripsimple1dhit = dynamic_cast<const SiStripRecHit1D*>(hit);
71  const SiStripRecHit2D* sistripsimplehit = dynamic_cast<const SiStripRecHit2D*>(hit);
72  const SiStripMatchedRecHit2D* sistripmatchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
73  const SiPixelRecHit* sipixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
74 
75  const SiPixelCluster* PixelCluster = NULL;
76  const SiStripCluster* StripCluster = NULL;
77  uint32_t DetId = 0;
78 
79  for(unsigned int h=0;h<2;h++){
80  if(!sistripmatchedhit && h==1){
81  continue;
82  }else if(sistripmatchedhit && h==0){
83  StripCluster = &sistripmatchedhit->monoCluster();
84  DetId = sistripmatchedhit->monoId();
85  }else if(sistripmatchedhit && h==1){
86  StripCluster = &sistripmatchedhit->stereoCluster();;
87  DetId = sistripmatchedhit->stereoId();
88  }else if(sistripsimplehit){
89  StripCluster = (sistripsimplehit->cluster()).get();
90  DetId = sistripsimplehit->geographicalId().rawId();
91  }else if(sistripsimple1dhit){
92  StripCluster = (sistripsimple1dhit->cluster()).get();
93  DetId = sistripsimple1dhit->geographicalId().rawId();
94  }else if(sipixelhit){
95  PixelCluster = (sipixelhit->cluster()).get();
96  DetId = sipixelhit->geographicalId().rawId();
97  }else{
98  continue;
99  }
100 
101  LocalVector trackDirection = trajState.localDirection();
102  double cosine = trackDirection.z()/trackDirection.mag();
103  bool Saturation = false;
104  bool Overlapping = false;
105  unsigned int Charge = 0;
106  double Path = (10.0*thickness(DetId))/fabs(cosine);
107  double PrevGain = -1;
108  double PrevGainTick = -1;
109  int FirstStrip = 0;
110  int NStrips = 0;
111 
112  if(StripCluster){
113  const auto & Ampls = StripCluster->amplitudes();
114  FirstStrip = StripCluster->firstStrip();
115  NStrips = Ampls.size();
116  int APVId = FirstStrip/128;
117 
118 
119  if(gainHandle.isValid()){
120  PrevGain = gainHandle->getApvGain(APVId,gainHandle->getRange(DetId, 1),1);
121  PrevGainTick = gainHandle->getApvGain(APVId,gainHandle->getRange(DetId, 0),1);
122  }
123 
124  for(unsigned int a=0;a<Ampls.size();a++){
125  Charge+=Ampls[a];
126  if(Ampls[a] >=254)Saturation =true;
127  amplitude->push_back( Ampls[a] );
128  }
129 
130  if(FirstStrip==0 )Overlapping=true;
131  if(FirstStrip==128 )Overlapping=true;
132  if(FirstStrip==256 )Overlapping=true;
133  if(FirstStrip==384 )Overlapping=true;
134  if(FirstStrip==512 )Overlapping=true;
135  if(FirstStrip==640 )Overlapping=true;
136 
137  if(FirstStrip<=127 && FirstStrip+Ampls.size()>127)Overlapping=true;
138  if(FirstStrip<=255 && FirstStrip+Ampls.size()>255)Overlapping=true;
139  if(FirstStrip<=383 && FirstStrip+Ampls.size()>383)Overlapping=true;
140  if(FirstStrip<=511 && FirstStrip+Ampls.size()>511)Overlapping=true;
141  if(FirstStrip<=639 && FirstStrip+Ampls.size()>639)Overlapping=true;
142 
143  if(FirstStrip+Ampls.size()==127 )Overlapping=true;
144  if(FirstStrip+Ampls.size()==255 )Overlapping=true;
145  if(FirstStrip+Ampls.size()==383 )Overlapping=true;
146  if(FirstStrip+Ampls.size()==511 )Overlapping=true;
147  if(FirstStrip+Ampls.size()==639 )Overlapping=true;
148  if(FirstStrip+Ampls.size()==767 )Overlapping=true;
149  }else if(PixelCluster){
150  const auto& Ampls = PixelCluster->pixelADC();
151  int FirstRow = PixelCluster->minPixelRow();
152  int FirstCol = PixelCluster->minPixelCol();
153  FirstStrip = ((FirstRow/80)<<3 | (FirstCol/52)) * 128; //Hack to save the APVId
154  NStrips = 0;
155  Saturation = false;
156  Overlapping = false;
157 
158  for(unsigned int a=0;a<Ampls.size();a++){
159  Charge+=Ampls[a];
160  if(Ampls[a] >=254)Saturation =true;
161  }
162  }
163  double ChargeOverPath = (double)Charge / Path ;
164 
165  trackindex ->push_back( shallow::findTrackIndex(tracks, track) );
166  rawid ->push_back( DetId );
167  localdirx ->push_back( trackDirection.x() );
168  localdiry ->push_back( trackDirection.y() );
169  localdirz ->push_back( trackDirection.z() );
170  firststrip ->push_back( FirstStrip );
171  nstrips ->push_back( NStrips );
172  saturation ->push_back( Saturation );
173  overlapping ->push_back( Overlapping );
174  farfromedge ->push_back( StripCluster ? IsFarFromBorder(&trajState,DetId, &iSetup) : true );
175  charge ->push_back( Charge );
176  path ->push_back( Path );
177  chargeoverpath->push_back( ChargeOverPath );
178  gainused ->push_back( PrevGain );
179  gainusedTick ->push_back( PrevGainTick );
180  }
181  }
182  }
183 
184  iEvent.put(trackindex, Prefix + "trackindex" + Suffix );
185  iEvent.put(rawid , Prefix + "rawid" + Suffix );
186  iEvent.put(localdirx , Prefix + "localdirx" + Suffix );
187  iEvent.put(localdiry , Prefix + "localdiry" + Suffix );
188  iEvent.put(localdirz , Prefix + "localdirz" + Suffix );
189  iEvent.put(firststrip, Prefix + "firststrip" + Suffix );
190  iEvent.put(nstrips, Prefix + "nstrips" + Suffix );
191  iEvent.put(saturation, Prefix + "saturation" + Suffix );
192  iEvent.put(overlapping, Prefix + "overlapping" + Suffix );
193  iEvent.put(farfromedge, Prefix + "farfromedge" + Suffix );
194  iEvent.put(charge, Prefix + "charge" + Suffix );
195  iEvent.put(path, Prefix + "path" + Suffix );
196  iEvent.put(chargeoverpath,Prefix + "chargeoverpath"+ Suffix );
197  iEvent.put(amplitude, Prefix + "amplitude" + Suffix );
198  iEvent.put(gainused, Prefix + "gainused" + Suffix );
199  iEvent.put(gainusedTick, Prefix + "gainusedTick" + Suffix );
200 }
201 
202 /*
203 void ShallowGainCalibration::beginJob(const edm::EventSetup& iSetup)
204 {
205  printf("Befin JOB\n");
206 
207  edm::ESHandle<TrackerGeometry> tkGeom;
208  iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
209  vector<GeomDet*> Det = tkGeom->dets();
210 
211  edm::ESHandle<SiStripGain> gainHandle;
212  iSetup.get<SiStripGainRcd>().get(gainHandle);
213  if(!gainHandle.isValid()){printf("\n#####################\n\nERROR --> gainHandle is not valid\n\n#####################\n\n");exit(0);}
214 
215  for(unsigned int i=0;i<Det.size();i++){
216  DetId Detid = Det[i]->geographicalId();
217  int SubDet = Detid.subdetId();
218 
219  if( SubDet == StripSubdetector::TIB || SubDet == StripSubdetector::TID ||
220  SubDet == StripSubdetector::TOB || SubDet == StripSubdetector::TEC ){
221 
222  StripGeomDetUnit* DetUnit = dynamic_cast<StripGeomDetUnit*> (Det[i]);
223  if(!DetUnit)continue;
224 
225  const StripTopology& Topo = DetUnit->specificTopology();
226  unsigned int NAPV = Topo.nstrips()/128;
227 
228  for(unsigned int j=0;j<NAPV;j++){
229  stAPVGain* APV = new stAPVGain;
230  APV->DetId = Detid.rawId();
231  APV->APVId = j;
232  APV->PreviousGain = 1;
233 
234  APVsCollOrdered.push_back(APV);
235  APVsColl[(APV->DetId<<3) | APV->APVId] = APV;
236  }
237  }
238  }
239 }
240 
241 
242 void ShallowGainCalibration::beginRun(edm::Run &, const edm::EventSetup &iSetup){
243  printf("BEFIN RUN\n");
244 
245  edm::ESHandle<SiStripGain> gainHandle;
246  iSetup.get<SiStripGainRcd>().get(gainHandle);
247  if(!gainHandle.isValid()){printf("\n#####################\n\nERROR --> gainHandle is not valid\n\n#####################\n\n");exit(0);}
248 
249  for(std::vector<stAPVGain*>::iterator it = APVsCollOrdered.begin();it!=APVsCollOrdered.end();it++){
250  stAPVGain* APV = *it;
251  SiStripApvGain::Range detGainRange = gainHandle->getRange(APV->DetId);
252  APV->PreviousGain = *(detGainRange.first + APV->APVId);
253  }
254 }
255 */
256 
258 {
259  edm::ESHandle<TrackerGeometry> tkGeom; iSetup->get<TrackerDigiGeometryRecord>().get( tkGeom );
260 
261  LocalPoint HitLocalPos = trajState->localPosition();
262  LocalError HitLocalError = trajState->localError().positionError() ;
263 
264  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
265  if (dynamic_cast<const StripGeomDetUnit*>(it)==0 && dynamic_cast<const PixelGeomDetUnit*>(it)==0) {
266  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
267  return false;
268  }
269 
270  const BoundPlane plane = it->surface();
271  const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
272  const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
273 
274  double DistFromBorder = 1.0;
275  double HalfLength = it->surface().bounds().length() /2.0;
276 
277  if(trapezoidalBounds)
278  {
279  std::array<const float, 4> const & parameters = (*trapezoidalBounds).parameters();
280  HalfLength = parameters[3];
281  }else if(rectangularBounds){
282  HalfLength = it->surface().bounds().length() /2.0;
283  }else{return false;}
284 
285  if (fabs(HitLocalPos.y())+HitLocalError.yy() >= (HalfLength - DistFromBorder) ) return false;
286 
287  return true;
288 }
289 
290 
292 {
293  map<DetId,double>::iterator th=m_thicknessMap.find(id);
294  if(th!=m_thicknessMap.end())
295  return (*th).second;
296  else {
297  double detThickness=1.;
298  //compute thickness normalization
299  const GeomDetUnit* it = m_tracker->idToDetUnit(DetId(id));
300  bool isPixel = dynamic_cast<const PixelGeomDetUnit*>(it)!=0;
301  bool isStrip = dynamic_cast<const StripGeomDetUnit*>(it)!=0;
302  if (!isPixel && ! isStrip) {
303  //FIXME throw exception
304  edm::LogWarning("DeDxHitsProducer") << "\t\t this detID doesn't seem to belong to the Tracker";
305  detThickness = 1.;
306  }else{
307  detThickness = it->surface().bounds().thickness();
308  }
309 
310  m_thicknessMap[id]=detThickness;//computed value
311  return detThickness;
312  }
313 
314 }
315 
ClusterRef cluster() const
virtual const TrackerGeomDet * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
dictionary parameters
Definition: Parameters.py:2
int minPixelCol() const
unsigned int stereoId() const
virtual float length() const =0
int findTrackIndex(const edm::Handle< edm::View< reco::Track > > &h, const reco::Track *t)
Definition: ShallowTools.cc:29
SiStripCluster const & monoCluster() const
LocalVector localDirection() const
std::map< DetId, double > m_thicknessMap
T y() const
Definition: PV3DBase.h:63
#define NULL
Definition: scimark2.h:8
const Bounds & bounds() const
Definition: Surface.h:128
uint16_t firstStrip() const
ShallowGainCalibration(const edm::ParameterSet &)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
LocalError positionError() const
bool IsFarFromBorder(TrajectoryStateOnSurface *trajState, const uint32_t detid, const edm::EventSetup *iSetup)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
DataContainer const & measurements() const
Definition: Trajectory.h:203
int iEvent
Definition: GenABIO.cc:230
T mag() const
Definition: PV3DBase.h:67
int minPixelRow() const
float yy() const
Definition: LocalError.h:26
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
Definition: Path.h:40
T z() const
Definition: PV3DBase.h:64
ClusterRef cluster() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
const LocalTrajectoryError & localError() const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
const std::vector< uint16_t > & pixelADC() const
Definition: DetId.h:18
const TrackerGeometry * m_tracker
tuple tracks
Definition: testEve_cfg.py:39
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
const T & get() const
Definition: EventSetup.h:56
void produce(edm::Event &, const edm::EventSetup &)
SiStripCluster const & stereoCluster() const
Pixel cluster – collection of neighboring pixels above threshold.
double a
Definition: hdecay.h:121
bool isPixel(HitType hitType)
tuple cout
Definition: gather_cfg.py:121
unsigned int monoId() const
DetId geographicalId() const
bool isValid() const
Definition: ESHandle.h:47
T x() const
Definition: PV3DBase.h:62
const std::vector< uint8_t > & amplitudes() const
Our base class.
Definition: SiPixelRecHit.h:23