1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
|
#include "Alignment/CommonAlignment/interface/Alignable.h"
#include "Alignment/CommonAlignment/interface/Utilities.h"
#include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
#include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
#include "DataFormats/Alignment/interface/AliClusterValueMap.h"
#include "DataFormats/Alignment/interface/AlignmentClusterFlag.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventPrincipal.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/one/EDProducer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "Geometry/Records/interface/IdealGeometryRecord.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
#include "RecoTracker/TransientTrackingRecHit/interface/TSiPixelRecHit.h"
#include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
#include "Utilities/General/interface/ClassName.h"
#include "TFile.h"
#include "TTree.h"
#include "TRandom3.h"
#include "TH1F.h"
#include <string>
class TrackerTopology;
class AlignmentPrescaler : public edm::one::EDProducer<> {
public:
AlignmentPrescaler(const edm::ParameterSet& iConfig);
~AlignmentPrescaler() override;
void beginJob() override;
void endJob() override;
void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
static void fillDescriptions(edm::ConfigurationDescriptions&);
private:
// tokens
edm::EDGetTokenT<reco::TrackCollection> tracksToken_; //tracks in input
edm::EDGetTokenT<AliClusterValueMap> aliClusterToken_; //Hit-quality association map
std::string prescfilename_; //name of the file containing the TTree with the prescaling factors
std::string presctreename_; //name of the TTree with the prescaling factors
TFile* fpresc_;
TTree* tpresc_;
TRandom3* myrand_;
int layerFromId(const DetId& id, const TrackerTopology* tTopo) const;
unsigned int detid_;
float hitPrescFactor_, overlapPrescFactor_;
int totnhitspxl_;
};
using namespace std;
AlignmentPrescaler::AlignmentPrescaler(const edm::ParameterSet& iConfig)
: tracksToken_(consumes(iConfig.getParameter<edm::InputTag>("src"))),
aliClusterToken_(consumes(iConfig.getParameter<edm::InputTag>("assomap"))),
prescfilename_(iConfig.getParameter<std::string>("PrescFileName")),
presctreename_(iConfig.getParameter<std::string>("PrescTreeName")) {
// issue the produce<>
produces<AliClusterValueMap>();
produces<AliTrackTakenClusterValueMap>();
}
AlignmentPrescaler::~AlignmentPrescaler() = default;
void AlignmentPrescaler::beginJob() {
//
edm::LogPrint("AlignmentPrescaler") << "in AlignmentPrescaler::beginJob" << std::flush;
fpresc_ = new TFile(prescfilename_.c_str(), "READ");
tpresc_ = (TTree*)fpresc_->Get(presctreename_.c_str());
tpresc_->BuildIndex("DetId");
tpresc_->SetBranchStatus("*", false);
tpresc_->SetBranchStatus("DetId", true);
tpresc_->SetBranchStatus("PrescaleFactor", true);
tpresc_->SetBranchStatus("PrescaleFactorOverlap", true);
edm::LogPrint("AlignmentPrescaler") << " Branches activated " << std::flush;
detid_ = 0;
hitPrescFactor_ = 99.0;
overlapPrescFactor_ = 88.0;
tpresc_->SetBranchAddress("DetId", &detid_);
tpresc_->SetBranchAddress("PrescaleFactor", &hitPrescFactor_);
tpresc_->SetBranchAddress("PrescaleFactorOverlap", &overlapPrescFactor_);
edm::LogPrint("AlignmentPrescaler") << " addressed " << std::flush;
myrand_ = new TRandom3();
// myrand_->SetSeed();
edm::LogPrint("AlignmentPrescaler") << " ok ";
}
void AlignmentPrescaler::endJob() {
delete tpresc_;
fpresc_->Close();
delete fpresc_;
delete myrand_;
}
void AlignmentPrescaler::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
LogDebug("AlignmentPrescaler") << "\n\n#################\n### Starting the AlignmentPrescaler::produce ; Event: "
<< iEvent.id().run() << ", " << iEvent.id().event() << std::endl;
edm::Handle<reco::TrackCollection> Tracks;
iEvent.getByToken(tracksToken_, Tracks);
//take HitAssomap
AliClusterValueMap InValMap = iEvent.get(aliClusterToken_);
//prepare the output of the ValueMap flagging tracks
std::vector<int> trackflags(Tracks->size(), 0);
//loop on tracks
for (std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk;
++ittrk) {
//loop on tracking rechits
LogDebug("AlignmentPrescaler") << "Loop on hits of track #" << (ittrk - Tracks->begin()) << std::endl;
int ntakenhits = 0;
bool firstTakenHit = false;
for (auto const& hit : ittrk->recHits()) {
if (!hit->isValid()) {
continue;
}
uint32_t tmpdetid = hit->geographicalId().rawId();
tpresc_->GetEntryWithIndex(tmpdetid);
//-------------
//decide whether to take this hit or not
bool takeit = false;
int subdetId = hit->geographicalId().subdetId();
//check first if the cluster is also in the overlap asso map
bool isOverlapHit = false;
// bool first=true;
//ugly...
const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
AlignmentClusterFlag tmpflag(hit->geographicalId());
int stripType = 0;
if (subdetId > 2) { // SST case
const std::type_info& type = typeid(*hit);
if (type == typeid(SiStripRecHit1D))
stripType = 1;
else if (type == typeid(SiStripRecHit2D))
stripType = 2;
else
stripType = 3;
if (stripType == 1) {
// const SiStripRecHit1D* stripHit1D = dynamic_cast<const SiStripRecHit1D*>(hit);
if (stripHit1D != nullptr) {
SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
tmpflag = InValMap[stripclust];
tmpflag.SetDetId(hit->geographicalId());
if (tmpflag.isOverlap())
isOverlapHit = true;
LogDebug("AlignmentPrescaler")
<< "~*~*~* Prescale (1D) for module " << tmpflag.detId().rawId() << "("
<< InValMap[stripclust].detId().rawId() << ") is " << hitPrescFactor_ << std::flush;
if (tmpflag.isOverlap())
LogDebug("AlignmentPrescaler") << " (it is Overlap)";
} //end if striphit1D!=0
} else if (stripType == 2) {
//const SiStripRecHit2D* stripHit2D = dynamic_cast<const SiStripRecHit2D*>(hit);
if (stripHit2D != nullptr) {
SiStripRecHit2D::ClusterRef stripclust(stripHit2D->cluster());
tmpflag = InValMap[stripclust];
tmpflag.SetDetId(hit->geographicalId());
if (tmpflag.isOverlap())
isOverlapHit = true;
LogDebug("AlignmentPrescaler")
<< "~*~*~* Prescale (2D) for module " << tmpflag.detId().rawId() << "("
<< InValMap[stripclust].detId().rawId() << ") is " << hitPrescFactor_ << std::flush;
if (tmpflag.isOverlap())
LogDebug("AlignmentPrescaler") << " (it is Overlap)" << endl;
} //end if striphit2D!=0
}
} //end if is a strip hit
else {
// const SiPixelRecHit* pixelhit= dynamic_cast<const SiPixelRecHit*>(hit);
if (pixelhit != nullptr) {
SiPixelClusterRefNew pixclust(pixelhit->cluster());
tmpflag = InValMap[pixclust];
tmpflag.SetDetId(hit->geographicalId());
if (tmpflag.isOverlap())
isOverlapHit = true;
}
} //end else is a pixel hit
// tmpflag.SetDetId(hit->geographicalId());
if (isOverlapHit) {
LogDebug("AlignmentPrescaler") << " DetId=" << tmpdetid << " is Overlap! " << flush;
takeit = (float(myrand_->Rndm()) <= overlapPrescFactor_);
}
if (!takeit) {
float rr = float(myrand_->Rndm());
takeit = (rr <= hitPrescFactor_);
}
if (takeit) { //HIT TAKEN !
LogDebug("AlignmentPrescaler") << " DetId=" << tmpdetid << " taken!" << flush;
tmpflag.SetTakenFlag();
if (subdetId > 2) {
if (stripType == 1) {
SiStripRecHit1D::ClusterRef stripclust(stripHit1D->cluster());
InValMap[stripclust] = tmpflag; //.SetTakenFlag();
} else if (stripType == 2) {
SiStripRecHit1D::ClusterRef stripclust(stripHit2D->cluster());
InValMap[stripclust] = tmpflag; //.SetTakenFlag();
} else
std::cout << "Unknown type of strip hit" << std::endl;
} else {
SiPixelClusterRefNew pixclust(pixelhit->cluster());
InValMap[pixclust] = tmpflag; //.SetTakenFlag();
}
if (!firstTakenHit) {
firstTakenHit = true;
LogDebug("AlignmentPrescaler") << "Index of the track iterator is " << ittrk - Tracks->begin();
}
ntakenhits++;
} //end if take this hit
} //end loop on RecHits
trackflags[ittrk - Tracks->begin()] = ntakenhits;
} //end loop on tracks
//save the asso map, tracks...
// prepare output
auto OutVM = std::make_unique<AliClusterValueMap>();
*OutVM = InValMap;
iEvent.put(std::move(OutVM));
auto trkVM = std::make_unique<AliTrackTakenClusterValueMap>();
AliTrackTakenClusterValueMap::Filler trkmapfiller(*trkVM);
trkmapfiller.insert(Tracks, trackflags.begin(), trackflags.end());
trkmapfiller.fill();
iEvent.put(std::move(trkVM));
} //end produce
int AlignmentPrescaler::layerFromId(const DetId& id, const TrackerTopology* tTopo) const {
if (uint32_t(id.subdetId()) == PixelSubdetector::PixelBarrel) {
return tTopo->pxbLayer(id);
} else if (uint32_t(id.subdetId()) == PixelSubdetector::PixelEndcap) {
return tTopo->pxfDisk(id) + (3 * (tTopo->pxfSide(id) - 1));
} else if (id.subdetId() == StripSubdetector::TIB) {
return tTopo->tibLayer(id);
} else if (id.subdetId() == StripSubdetector::TOB) {
return tTopo->tobLayer(id);
} else if (id.subdetId() == StripSubdetector::TEC) {
return tTopo->tecWheel(id) + (9 * (tTopo->pxfSide(id) - 1));
} else if (id.subdetId() == StripSubdetector::TID) {
return tTopo->tidWheel(id) + (3 * (tTopo->tidSide(id) - 1));
}
return -1;
} //end layerfromId
void AlignmentPrescaler::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.setComment("Prescale Tracker Alignment hits for alignment procedures");
desc.add<edm::InputTag>("src", edm::InputTag("generalTracks"));
desc.add<edm::InputTag>("assomap", edm::InputTag("OverlapAssoMap"));
desc.add<std::string>("PrescFileName", "PrescaleFactors.root");
desc.add<std::string>("PrescTreeName", "AlignmentHitMap");
descriptions.addWithDefaultLabel(desc);
}
// ========= MODULE DEF ==============
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(AlignmentPrescaler);
|