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