CMS 3D CMS Logo

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("*",false);
51  tpresc_->SetBranchStatus("DetId",true);
52  tpresc_->SetBranchStatus("PrescaleFactor",true);
53  tpresc_->SetBranchStatus("PrescaleFactorOverlap",true);
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(auto const& hit : ittrk->recHits()) {
102  if(! hit->isValid()){
103  nhit++;
104  continue;
105  }
106  uint32_t tmpdetid = hit->geographicalId().rawId();
107  tpresc_->GetEntryWithIndex(tmpdetid);
108 
109 
110  //-------------
111  //decide whether to take this hit or not
112  bool takeit=false;
113  int subdetId=hit->geographicalId().subdetId();
114 
115 
116  //check first if the cluster is also in the overlap asso map
117  bool isOverlapHit=false;
118  // bool first=true;
119  //ugly...
120  const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
121  const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
122  const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
123 
124  AlignmentClusterFlag tmpflag(hit->geographicalId());
125  int stripType=0;
126  if(subdetId>2){// SST case
127  const std::type_info &type = typeid(*hit);
128  if (type == typeid(SiStripRecHit1D)) stripType=1;
129  else if (type == typeid(SiStripRecHit2D)) stripType=2;
130  else stripType=3;
131 
132  if(stripType==1) {
133  // const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
134 
135  if(stripHit1D!=nullptr){
136  SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
137  tmpflag=InValMap[stripclust];
138  tmpflag.SetDetId(hit->geographicalId());
139  if(tmpflag.isOverlap())isOverlapHit=true;
140  // std::cout<<"~*~*~* Prescale (1D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
141  // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
142  // else cout<<endl;
143 
144  }//end if striphit1D!=0
145  }
146  else if (stripType==2) {
147  //const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
148  if(stripHit2D!=nullptr){
149  SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster());
150  tmpflag=InValMap[stripclust];
151  tmpflag.SetDetId(hit->geographicalId());
152  if(tmpflag.isOverlap())isOverlapHit=true;
153  // std::cout<<"~*~*~* Prescale (2D) for module "<<tmpflag.detId().rawId()<<"("<<InValMap[stripclust].detId().rawId() <<") is "<<hitPrescFactor_<<std::flush;
154  // if(tmpflag.isOverlap())cout<<" (it is Overlap)"<<endl;
155  // else cout<<endl;
156 
157  }//end if striphit2D!=0
158  }
159  }//end if is a strip hit
160  else{
161  // const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
162  if(pixelhit!=nullptr){
163  //npxlhits++;
164  SiPixelClusterRefNew pixclust(pixelhit->cluster());
165  tmpflag=InValMap[pixclust];
166  tmpflag.SetDetId(hit->geographicalId());
167  if(tmpflag.isOverlap())isOverlapHit=true;
168  }
169  }//end else is a pixel hit
170  // tmpflag.SetDetId(hit->geographicalId());
171 
172  if( isOverlapHit ){
173  //cout<<" DetId="<<tmpdetid<<" is Overlap! "<<flush;
174  takeit=(float(myrand_->Rndm())<=overlapPrescFactor_);
175  }
176  if( !takeit ){
177  float rr=float(myrand_->Rndm());
178  takeit=(rr<=hitPrescFactor_);
179  }
180  if(takeit){//HIT TAKEN !
181  //cout<<" DetId="<<tmpdetid<<" taken!"<<flush;
182  tmpflag.SetTakenFlag();
183 
184  if(subdetId>2){
185  if(stripType==1){
186  SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
187  InValMap[stripclust]=tmpflag;//.SetTakenFlag();
188  }
189  else if(stripType==2){
190  SiStripRecHit1D::ClusterRef stripclust(stripHit2D->cluster());
191  InValMap[stripclust]=tmpflag;//.SetTakenFlag();
192  }
193  else std::cout<<"Unknown type of strip hit"<<std::endl;
194  }
195  else{
196  SiPixelClusterRefNew pixclust(pixelhit->cluster());
197  InValMap[pixclust]=tmpflag;//.SetTakenFlag();
198  }
199 
200  if(!firstTakenHit){
201  firstTakenHit=true;
202  //std::cout<<"Index of the track iterator is "<< ittrk-Tracks->begin() <<endl;
203 
204  }
205  ntakenhits++;
206  }//end if take this hit
207  //cout<<endl;
208 
209  nhit++;
210  //cout<<endl;
211  }//end loop on RecHits
212  trackflags[ittrk-Tracks->begin()]=ntakenhits;
213 
214  }//end loop on tracks
215 
216 
217 
218  // totnhitspxl_+=ntakenhits;
219  //cout<<"AlignmentPrescaler::produce says that in this event "<<ntakenhits<<" pixel clusters were taken (out of "<<npxlhits<<" total pixel hits."<<endl;
220 
221 
222 
223  //save the asso map, tracks...
224  // prepare output
225  auto OutVM = std::make_unique<AliClusterValueMap>();
226  *OutVM=InValMap;
227 
228  iEvent.put(std::move(OutVM));
229 
230 
231  auto trkVM = std::make_unique<AliTrackTakenClusterValueMap>();
232  AliTrackTakenClusterValueMap::Filler trkmapfiller(*trkVM);
233  trkmapfiller.insert(Tracks,trackflags.begin(),trackflags.end() );
234  trkmapfiller.fill();
235  iEvent.put(std::move(trkVM));
236 
237 
238 }//end produce
239 
240 
241 int AlignmentPrescaler::layerFromId (const DetId& id, const TrackerTopology* tTopo) const
242 {
243  if ( uint32_t(id.subdetId())==PixelSubdetector::PixelBarrel ) {
244 
245  return tTopo->pxbLayer(id);
246  }
247  else if ( uint32_t(id.subdetId())==PixelSubdetector::PixelEndcap ) {
248 
249  return tTopo->pxfDisk(id) + (3*(tTopo->pxfSide(id)-1));
250  }
251  else if ( id.subdetId()==StripSubdetector::TIB ) {
252 
253  return tTopo->tibLayer(id);
254  }
255  else if ( id.subdetId()==StripSubdetector::TOB ) {
256 
257  return tTopo->tobLayer(id);
258  }
259  else if ( id.subdetId()==StripSubdetector::TEC ) {
260 
261  return tTopo->tecWheel(id) + (9*(tTopo->pxfSide(id)-1));
262  }
263  else if ( id.subdetId()==StripSubdetector::TID ) {
264 
265  return tTopo->tidWheel(id) + (3*(tTopo->tidSide(id)-1));
266  }
267  return -1;
268 
269 }//end layerfromId
270 
271 // ========= MODULE DEF ==============
ClusterRef cluster() const
type
Definition: HCALResponse.h:21
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
unsigned int tibLayer(const DetId &id) const
AlignmentPrescaler(const edm::ParameterSet &iConfig)
unsigned int pxfDisk(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
~AlignmentPrescaler() override
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void beginJob() override
edm::InputTag srcQualityMap_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int tidSide(const DetId &id) const
ClusterRef cluster() const
std::string prescfilename_
std::string presctreename_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
unsigned int pxbLayer(const DetId &id) const
Definition: DetId.h:18
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
void endJob() override
HLT enums.
unsigned int pxfSide(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
def move(src, dest)
Definition: eostools.py:511
int layerFromId(const DetId &id, const TrackerTopology *tTopo) const
unsigned int tobLayer(const DetId &id) const
Our base class.
Definition: SiPixelRecHit.h:23