CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentPrescaler.cc
Go to the documentation of this file.
2 
12 
15 
23 
24 // #include "Riostream.h"
25 
26 using namespace std;
27 
29  src_(iConfig.getParameter<edm::InputTag>("src")),
30  srcQualityMap_(iConfig.getParameter<edm::InputTag>("assomap")),
31  prescfilename_(iConfig.getParameter<std::string>("PrescFileName")),
32  presctreename_(iConfig.getParameter<std::string>("PrescTreeName"))
33 {
34  // issue the produce<>
35  produces<AliClusterValueMap>();
36  produces<AliTrackTakenClusterValueMap>();
37 
38 }
39 
41  //
42 }
43 
45  //
46  std::cout<<"in AlignmentPrescaler::beginJob"<<std::flush;
47  fpresc_=new TFile(prescfilename_.c_str(),"READ");
48  tpresc_=(TTree*)fpresc_->Get(presctreename_.c_str());
49  tpresc_->BuildIndex("DetId");
50  tpresc_->SetBranchStatus("*",0);
51  tpresc_->SetBranchStatus("DetId",1);
52  tpresc_->SetBranchStatus("PrescaleFactor",1);
53  tpresc_->SetBranchStatus("PrescaleFactorOverlap",1);
54  cout<<" Branches activated "<<std::flush;
55  detid_=0;
56  hitPrescFactor_=99.0;
58 
59  tpresc_->SetBranchAddress("DetId",&detid_);
60  tpresc_->SetBranchAddress("PrescaleFactor",&hitPrescFactor_);
61  tpresc_->SetBranchAddress("PrescaleFactorOverlap",&overlapPrescFactor_);
62  cout<<" addressed "<<std::flush;
63  myrand_=new TRandom3();
64  // myrand_->SetSeed();
65  cout<<" ok "<<std::endl;
66 
67 }
68 
70 
71  delete tpresc_;
72  fpresc_->Close();
73  delete fpresc_;
74  delete myrand_;
75 }
76 
78  // std::cout<<"\n\n#################\n### Starting the AlignmentPrescaler::produce ; Event: "<<iEvent.id().run() <<", "<<iEvent.id().event()<<std::endl;
80  iEvent.getByLabel(src_, Tracks);
81 
82  //take HitAssomap
84  iEvent.getByLabel(srcQualityMap_, hMap);
85  AliClusterValueMap InValMap=*hMap;
86 
87  //prepare the output of the ValueMap flagging tracks
88  std::vector<int> trackflags(Tracks->size(),0);
89 
90 
91  //int npxlhits=0;
92 
93  //loop on tracks
94  for(std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk; ++ittrk){
95  //loop on tracking rechits
96  // std::cout << "Loop on hits of track #" << (ittrk - Tracks->begin()) << std::endl;
97  int nhit=0;
98  int ntakenhits=0;
99  bool firstTakenHit=false;
100 
101  for (trackingRecHit_iterator ith = ittrk->recHitsBegin(), edh = ittrk->recHitsEnd(); ith != edh; ++ith) {
102  const TrackingRecHit *hit = *ith; // ith is an iterator on edm::Ref to rechit
103  if(! hit->isValid()){
104  nhit++;
105  continue;
106  }
107  uint32_t tmpdetid = hit->geographicalId().rawId();
108  tpresc_->GetEntryWithIndex(tmpdetid);
109 
110 
111  //-------------
112  //decide whether to take this hit or not
113  bool takeit=false;
114  int subdetId=hit->geographicalId().subdetId();
115 
116 
117  //check first if the cluster is also in the overlap asso map
118  bool isOverlapHit=false;
119  // bool first=true;
120  //ugly...
121  const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
122  const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
123  const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
124 
125  AlignmentClusterFlag tmpflag(hit->geographicalId());
126  int stripType=0;
127  if(subdetId>2){// SST case
128  const std::type_info &type = typeid(*hit);
129  if (type == typeid(SiStripRecHit1D)) stripType=1;
130  else if (type == typeid(SiStripRecHit2D)) stripType=2;
131  else stripType=3;
132 
133  if(stripType==1) {
134  // const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
135 
136  if(stripHit1D!=0){
137  SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
138  tmpflag=InValMap[stripclust];
139  tmpflag.SetDetId(hit->geographicalId());
140  if(tmpflag.isOverlap())isOverlapHit=true;
141  // std::cout<<"~*~*~* Prescale (1D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
142  // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
143  // else cout<<endl;
144 
145  }//end if striphit1D!=0
146  }
147  else if (stripType==2) {
148  //const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
149  if(stripHit2D!=0){
150  SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster());
151  tmpflag=InValMap[stripclust];
152  tmpflag.SetDetId(hit->geographicalId());
153  if(tmpflag.isOverlap())isOverlapHit=true;
154  // std::cout<<"~*~*~* Prescale (2D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
155  // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
156  // else cout<<endl;
157 
158  }//end if striphit2D!=0
159  }
160  }//end if is a strip hit
161  else{
162  // const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
163  if(pixelhit!=0){
164  //npxlhits++;
165  SiPixelClusterRefNew pixclust(pixelhit->cluster());
166  tmpflag=InValMap[pixclust];
167  tmpflag.SetDetId(hit->geographicalId());
168  if(tmpflag.isOverlap())isOverlapHit=true;
169  }
170  }//end else is a pixel hit
171  // tmpflag.SetDetId(hit->geographicalId());
172 
173  if( isOverlapHit ){
174  //cout<<" DetId="<<tmpdetid<<" is Overlap! "<<flush;
175  takeit=(float(myrand_->Rndm())<=overlapPrescFactor_);
176  }
177  if( !takeit ){
178  float rr=float(myrand_->Rndm());
179  takeit=(rr<=hitPrescFactor_);
180  }
181  if(takeit){//HIT TAKEN !
182  //cout<<" DetId="<<tmpdetid<<" taken!"<<flush;
183  tmpflag.SetTakenFlag();
184 
185  if(subdetId>2){
186  if(stripType==1){
187  SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
188  InValMap[stripclust]=tmpflag;//.SetTakenFlag();
189  }
190  else if(stripType==2){
191  SiStripRecHit1D::ClusterRef stripclust(stripHit2D->cluster());
192  InValMap[stripclust]=tmpflag;//.SetTakenFlag();
193  }
194  else std::cout<<"Unknown type of strip hit"<<std::endl;
195  }
196  else{
197  SiPixelClusterRefNew pixclust(pixelhit->cluster());
198  InValMap[pixclust]=tmpflag;//.SetTakenFlag();
199  }
200 
201  if(!firstTakenHit){
202  firstTakenHit=true;
203  //std::cout<<"Index of the track iterator is "<< ittrk-Tracks->begin() <<endl;
204 
205  }
206  ntakenhits++;
207  }//end if take this hit
208  //cout<<endl;
209 
210  nhit++;
211  //cout<<endl;
212  }//end loop on RecHits
213  trackflags[ittrk-Tracks->begin()]=ntakenhits;
214 
215  }//end loop on tracks
216 
217 
218 
219  // totnhitspxl_+=ntakenhits;
220  //cout<<"AlignmentPrescaler::produce says that in this event "<<ntakenhits<<" pixel clusters were taken (out of "<<npxlhits<<" total pixel hits."<<endl;
221 
222 
223 
224  //save the asso map, tracks...
225  // prepare output
226  std::auto_ptr<AliClusterValueMap> OutVM( new AliClusterValueMap);
227  *OutVM=InValMap;
228 
229  iEvent.put(OutVM);
230 
231 
232  std::auto_ptr<AliTrackTakenClusterValueMap> trkVM( new AliTrackTakenClusterValueMap);
233  AliTrackTakenClusterValueMap::Filler trkmapfiller(*trkVM);
234  trkmapfiller.insert(Tracks,trackflags.begin(),trackflags.end() );
235  trkmapfiller.fill();
236  iEvent.put(trkVM);
237 
238 
239 }//end produce
240 
241 
242 int AlignmentPrescaler::layerFromId (const DetId& id, const TrackerTopology* tTopo) const
243 {
244  if ( uint32_t(id.subdetId())==PixelSubdetector::PixelBarrel ) {
245 
246  return tTopo->pxbLayer(id);
247  }
248  else if ( uint32_t(id.subdetId())==PixelSubdetector::PixelEndcap ) {
249 
250  return tTopo->pxfDisk(id) + (3*(tTopo->pxfSide(id)-1));
251  }
252  else if ( id.subdetId()==StripSubdetector::TIB ) {
253 
254  return tTopo->tibLayer(id);
255  }
256  else if ( id.subdetId()==StripSubdetector::TOB ) {
257 
258  return tTopo->tobLayer(id);
259  }
260  else if ( id.subdetId()==StripSubdetector::TEC ) {
261 
262  return tTopo->tecWheel(id) + (9*(tTopo->pxfSide(id)-1));
263  }
264  else if ( id.subdetId()==StripSubdetector::TID ) {
265 
266  return tTopo->tidWheel(id) + (3*(tTopo->tidSide(id)-1));
267  }
268  return -1;
269 
270 }//end layerfromId
271 
272 // ========= MODULE DEF ==============
ClusterRef cluster() const
type
Definition: HCALResponse.h:21
unsigned int tibLayer(const DetId &id) const
AlignmentPrescaler(const edm::ParameterSet &iConfig)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
unsigned int pxfDisk(const DetId &id) const
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
unsigned int tidWheel(const DetId &id) const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
edm::InputTag srcQualityMap_
int iEvent
Definition: GenABIO.cc:230
unsigned int tidSide(const DetId &id) const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
ClusterRef cluster() const
std::string prescfilename_
std::string presctreename_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:418
unsigned int pxbLayer(const DetId &id) const
Definition: DetId.h:18
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
bool isValid() const
unsigned int pxfSide(const DetId &id) const
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
tuple cout
Definition: gather_cfg.py:145
DetId geographicalId() const
unsigned int tecWheel(const DetId &id) const
int layerFromId(const DetId &id, const TrackerTopology *tTopo) const
unsigned int tobLayer(const DetId &id) const
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Our base class.
Definition: SiPixelRecHit.h:23