69 template <
typename TEle,
typename TCand>
72 if (iConfig.
existsAs<std::vector<edm::InputTag>>(
"electronTags")) {
73 electronTags_ = iConfig.
getParameter<std::vector<edm::InputTag>>(
"electronTags");
74 if (electronTags_.empty())
76 <<
"[SelectedElectronFEDListProducer] empty electron collection is given --> at least one \n";
79 <<
"[SelectedElectronFEDListProducer] no electron collection are given --> need at least one \n";
82 LogDebug(
"SelectedElectronFEDListProducer") <<
" Electron Collections" << std::endl;
83 for (std::vector<edm::InputTag>::const_iterator itEleTag = electronTags_.begin(); itEleTag != electronTags_.end();
85 electronToken_.push_back(consumes<TEleColl>(*itEleTag));
86 LogDebug(
"SelectedElectronFEDListProducer") <<
" Ele collection: " << *(itEleTag) << std::endl;
90 if (iConfig.
existsAs<std::vector<edm::InputTag>>(
"recoEcalCandidateTags")) {
91 recoEcalCandidateTags_ = iConfig.
getParameter<std::vector<edm::InputTag>>(
"recoEcalCandidateTags");
92 if (recoEcalCandidateTags_.empty())
93 throw cms::Exception(
"Configuration") <<
"[SelectedElectronFEDListProducer] empty ecal candidate collections "
94 "collection is given --> at least one \n";
96 throw cms::Exception(
"Configuration") <<
"[SelectedElectronFEDListProducer] no electron reco ecal candidate "
97 "collection are given --> need at least one \n";
100 for (std::vector<edm::InputTag>::const_iterator itEcalCandTag = recoEcalCandidateTags_.begin();
101 itEcalCandTag != recoEcalCandidateTags_.end();
103 recoEcalCandidateToken_.push_back(consumes<trigger::TriggerFilterObjectWithRefs>(*itEcalCandTag));
104 LogDebug(
"SelectedElectronFEDListProducer") <<
" Reco ecal candidate collection: " << *(itEcalCandTag) << std::endl;
108 if (iConfig.
existsAs<std::vector<int>>(
"isGsfElectronCollection")) {
109 isGsfElectronCollection_ = iConfig.
getParameter<std::vector<int>>(
"isGsfElectronCollection");
110 if (isGsfElectronCollection_.empty())
112 <<
"[SelectedElectronFEDListProducer] empty electron flag collection --> at least one \n";
115 <<
"[SelectedElectronFEDListProducer] no electron flag are given --> need at least one \n";
117 if (isGsfElectronCollection_.size() != electronTags_.size()
or
118 isGsfElectronCollection_.size() != recoEcalCandidateTags_.size())
119 throw cms::Exception(
"Configuration") <<
"[SelectedElectronFEDListProducer] electron flag , electron collection "
120 "and reco ecal cand collection must have the same size ! \n";
129 beamSpotToken_ = consumes<reco::BeamSpot>(beamSpotTag_);
131 LogDebug(
"SelectedElectronFEDListProducer") <<
" Beam Spot Tag " << beamSpotTag_ << std::endl;
140 hbheRecHitToken_ = consumes<HBHERecHitCollection>(HBHERecHitTag_);
149 rawDataToken_ = consumes<FEDRawDataCollection>(rawDataTag_);
151 LogDebug(
"SelectedElectronFEDListProducer") <<
" RawDataInput " << rawDataTag_ << std::endl;
154 if (iConfig.
existsAs<std::vector<int>>(
"addThisSelectedFEDs")) {
155 addThisSelectedFEDs_ = iConfig.
getParameter<std::vector<int>>(
"addThisSelectedFEDs");
156 if (addThisSelectedFEDs_.empty())
157 addThisSelectedFEDs_.push_back(-1);
159 addThisSelectedFEDs_.push_back(-1);
161 std::vector<int>::const_iterator AddFed = addThisSelectedFEDs_.begin();
162 for (; AddFed != addThisSelectedFEDs_.end(); ++AddFed)
163 LogDebug(
"SelectedElectronFEDListProducer") <<
" Additional FED: " << *(AddFed) << std::endl;
169 ESLookupTable_ =
edm::FileInPath(
"EventFilter/ESDigiToRaw/data/ES_lookup_table.dat");
175 outputLabelModule_ =
"streamElectronRawData";
177 LogDebug(
"SelectedElectronFEDListProducer") <<
" Output Label " << outputLabelModule_ << std::endl;
180 if (iConfig.
existsAs<
double>(
"dRStripRegion"))
181 dRStripRegion_ = iConfig.
getParameter<
double>(
"dRStripRegion");
183 dRStripRegion_ = 0.5;
185 LogDebug(
"SelectedElectronFEDListProducer") <<
" dRStripRegion " << dRStripRegion_ << std::endl;
188 if (iConfig.
existsAs<
double>(
"dRHcalRegion"))
189 dRHcalRegion_ = iConfig.
getParameter<
double>(
"dRHcalRegion");
194 if (iConfig.
existsAs<
double>(
"dPhiPixelRegion"))
195 dPhiPixelRegion_ = iConfig.
getParameter<
double>(
"dPhiPixelRegion");
197 dPhiPixelRegion_ = 0.5;
199 if (iConfig.
existsAs<
double>(
"dEtaPixelRegion"))
200 dEtaPixelRegion_ = iConfig.
getParameter<
double>(
"dEtaPixelRegion");
202 dEtaPixelRegion_ = 0.5;
204 if (iConfig.
existsAs<
double>(
"maxZPixelRegion"))
205 maxZPixelRegion_ = iConfig.
getParameter<
double>(
"maxZPixelRegion");
207 maxZPixelRegion_ = 24.;
209 LogDebug(
"SelectedElectronFEDListProducer")
210 <<
" dPhiPixelRegion " << dPhiPixelRegion_ <<
" dEtaPixelRegion " << dEtaPixelRegion_ <<
" MaxZPixelRegion "
211 << maxZPixelRegion_ << std::endl;
214 if (iConfig.
existsAs<
bool>(
"dumpSelectedEcalFed"))
215 dumpSelectedEcalFed_ = iConfig.
getParameter<
bool>(
"dumpSelectedEcalFed");
217 dumpSelectedEcalFed_ =
true;
219 if (iConfig.
existsAs<
bool>(
"dumpSelectedSiStripFed"))
220 dumpSelectedSiStripFed_ = iConfig.
getParameter<
bool>(
"dumpSelectedSiStripFed");
222 dumpSelectedSiStripFed_ =
true;
224 if (iConfig.
existsAs<
bool>(
"dumpSelectedSiPixelFed"))
225 dumpSelectedSiPixelFed_ = iConfig.
getParameter<
bool>(
"dumpSelectedSiPixelFed");
227 dumpSelectedSiPixelFed_ =
true;
229 if (iConfig.
existsAs<
bool>(
"dumpSelectedHCALFed"))
230 dumpSelectedHCALFed_ = iConfig.
getParameter<
bool>(
"dumpSelectedHCALFed");
232 dumpSelectedHCALFed_ =
true;
234 LogDebug(
"SelectedElectronFEDListProducer")
235 <<
" DumpEcalFedList set to " << dumpSelectedEcalFed_ <<
" DumpSelectedSiStripFed " << dumpSelectedSiStripFed_
236 <<
" DumpSelectedSiPixelFed " << dumpSelectedSiPixelFed_ << std::endl;
238 if (iConfig.
existsAs<
bool>(
"dumpAllEcalFed"))
239 dumpAllEcalFed_ = iConfig.
getParameter<
bool>(
"dumpAllEcalFed");
241 dumpAllEcalFed_ =
false;
243 if (iConfig.
existsAs<
bool>(
"dumpAllTrackerFed"))
244 dumpAllTrackerFed_ = iConfig.
getParameter<
bool>(
"dumpAllTrackerFed");
246 dumpAllTrackerFed_ =
false;
248 if (iConfig.
existsAs<
bool>(
"dumpAllHCALFed"))
249 dumpAllHCALFed_ = iConfig.
getParameter<
bool>(
"dumpAllHCALFed");
251 dumpAllHCALFed_ =
false;
253 LogDebug(
"SelectedElectronFEDListProducer")
254 <<
" DumpAllEcalFed " << dumpAllEcalFed_ <<
" DumpAllTrackerFed " << dumpAllTrackerFed_ <<
" Dump all HCAL fed "
255 << dumpAllHCALFed_ << std::endl;
258 for (
int i = 0;
i < 2; ++
i)
259 for (
int j = 0;
j < 2; ++
j)
260 for (
int k = 0;
k < 40; ++
k)
261 for (
int m = 0;
m < 40;
m++)
262 ES_fedId_[
i][
j][
k][
m] = -1;
265 int nLines, iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
266 std::ifstream ES_file;
267 ES_file.open(ESLookupTable_.fullPath().c_str());
268 LogDebug(
"SelectedElectronFEDListProducer")
269 <<
" Look Up table for ES " << ESLookupTable_.fullPath().c_str() << std::endl;
270 if (ES_file.is_open()) {
273 ES_file >> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
274 ES_fedId_[(3 - iz) / 2 - 1][ip - 1][ix - 1][iy - 1] = fed;
277 LogDebug(
"SelectedElectronFEDListProducer")
278 <<
" Look up table file can not be found in " << ESLookupTable_.fullPath().c_str() << std::endl;
282 produces<FEDRawDataCollection>(outputLabelModule_);
285 template <
typename TEle,
typename TCand>
287 if (!electronTags_.empty())
288 electronTags_.clear();
289 if (!recoEcalCandidateTags_.empty())
290 recoEcalCandidateTags_.clear();
291 if (!recoEcalCandidateToken_.empty())
292 recoEcalCandidateToken_.clear();
293 if (!electronToken_.empty())
294 electronToken_.clear();
295 if (!fedList_.empty())
297 if (!pixelModuleVector_.empty())
298 pixelModuleVector_.clear();
301 template <
typename TEle,
typename TCand>
303 LogDebug(
"SelectedElectronFEDListProducer") <<
" Begin of the Job " << std::endl;
306 template <
typename TEle,
typename TCand>
316 EcalMapping_ = ecalmapping.
product();
321 GeometryCalo_ = caloGeometry.
product();
329 PixelCabling_.reset();
335 if (pixelModuleVector_.empty()) {
337 std::vector<const GeomDet*>::const_iterator itTracker = trackerGeometry->
dets().begin();
338 for (; itTracker != trackerGeometry->
dets().end(); ++itTracker) {
339 int subdet = (*itTracker)->geographicalId().subdetId();
343 module.
x = (*itTracker)->position().x();
344 module.
y = (*itTracker)->position().y();
345 module.
z = (*itTracker)->position().z();
346 module.
Phi = (*itTracker)->position().phi();
347 module.
Eta = (*itTracker)->position().eta();
348 module.
DetId = (*itTracker)->geographicalId().rawId();
349 const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = PixelCabling_->pathToDetUnit(module.
DetId);
350 module.
Fed = path2det[0].fed;
352 pixelModuleVector_.push_back(module);
354 std::sort(pixelModuleVector_.begin(), pixelModuleVector_.end());
359 StripRegionCabling_ = SiStripCablingHandle.
product();
363 regionDimension_ = StripRegionCabling_->regionDimensions();
369 iEvent.getByToken(rawDataToken_, rawdata);
376 beamSpotPosition_ =
beamSpot->position();
378 beamSpotPosition_.SetXYZ(0, 0, 0);
383 iEvent.getByToken(hbheRecHitToken_, hbheRecHitHandle);
386 hcalRecHitCollection = hbheRecHitHandle.
product();
388 double radTodeg = 180. /
Geom::pi();
390 if (dumpAllEcalFed_) {
392 fedList_.push_back(iEcalFed);
394 fedList_.push_back(iESFed);
397 if (dumpAllTrackerFed_) {
399 fedList_.push_back(iPixelFed);
401 fedList_.push_back(iStripFed);
404 if (dumpAllHCALFed_) {
406 fedList_.push_back(iHcalFed);
414 std::vector<edm::Ref<TCandColl>> recoEcalCandColl;
417 typename std::vector<edm::EDGetTokenT<TEleColl>>::const_iterator itElectronColl = electronToken_.begin();
418 std::vector<int>::const_iterator itElectronCollFlag = isGsfElectronCollection_.begin();
419 std::vector<edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>>::const_iterator itRecoEcalCandColl =
420 recoEcalCandidateToken_.begin();
423 if (!dumpAllTrackerFed_ || !dumpAllEcalFed_) {
425 for (; itRecoEcalCandColl != recoEcalCandidateToken_.end() and itElectronColl != electronToken_.end() and
426 itElectronCollFlag != isGsfElectronCollection_.end();
427 ++itElectronColl, ++itElectronCollFlag, ++itRecoEcalCandColl) {
429 iEvent.getByToken(*itRecoEcalCandColl, triggerRecoEcalCandidateCollection);
430 if (triggerRecoEcalCandidateCollection.
failedToGet())
439 if (recoEcalCandColl.empty())
441 if (recoEcalCandColl.empty())
444 typename std::vector<edm::Ref<TCandColl>>::const_iterator itRecoEcalCand =
445 recoEcalCandColl.begin();
448 for (; itRecoEcalCand != recoEcalCandColl.end(); ++itRecoEcalCand) {
449 recoEcalCand = (*itRecoEcalCand);
451 recoEcalCand->superCluster();
453 typename TEleColl::const_iterator itEle =
electrons->begin();
454 for (; itEle !=
electrons->end(); ++itEle) {
458 if (scRefRecoEcalCand != scRef)
461 const std::vector<std::pair<DetId, float>>&
hits = scRef->hitsAndFractions();
463 std::vector<std::pair<DetId, float>>::const_iterator itSChits =
hits.begin();
464 if (!dumpAllEcalFed_) {
465 for (; itSChits !=
hits.end(); ++itSChits) {
466 if ((*itSChits).first.subdetId() ==
EcalBarrel) {
467 EBDetId idEBRaw((*itSChits).first);
470 EcalMapping_->GetFED(
double(
point.eta()),
double(
point.phi()) * radTodeg);
474 LogDebug(
"SelectedElectronFEDListProducer")
475 <<
" electron hit detID Barrel " << (*itSChits).first.rawId() <<
" eta " << double(
point.eta())
476 <<
" phi " <<
double(
point.phi()) * radTodeg <<
" FED " << hitFED << std::endl;
478 if (dumpSelectedEcalFed_) {
479 if (!fedList_.empty()) {
480 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
481 fedList_.push_back(hitFED);
483 fedList_.push_back(hitFED);
485 }
else if ((*itSChits).first.subdetId() ==
EcalEndcap) {
486 EEDetId idEERaw((*itSChits).first);
489 EcalMapping_->GetFED(
double(
point.eta()),
double(
point.phi()) * radTodeg);
493 LogDebug(
"SelectedElectronFEDListProducer")
494 <<
" electron hit detID Endcap " << (*itSChits).first.rawId() <<
" eta " << double(
point.eta())
495 <<
" phi " <<
double(
point.phi()) * radTodeg <<
" FED " << hitFED << std::endl;
496 if (dumpSelectedEcalFed_) {
497 if (!fedList_.empty()) {
498 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
499 fedList_.push_back(hitFED);
501 fedList_.push_back(hitFED);
505 (dynamic_cast<const EcalPreshowerGeometry*>(GeometryES_))->getClosestCellInPlane(
point, 1);
508 ES_fedId_[(3 - stripX.
zside()) / 2 - 1][stripX.
plane() - 1][stripX.
six() - 1][stripX.
siy() - 1];
509 LogDebug(
"SelectedElectronFEDListProducer")
510 <<
" ES hit plane X (deiID) " << stripX.
rawId() <<
" six " << stripX.
six() <<
" siy "
511 << stripX.
siy() <<
" plane " << stripX.
plane() <<
" FED ID " << hitFED << std::endl;
516 if (!fedList_.empty()) {
517 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
518 fedList_.push_back(hitFED);
520 fedList_.push_back(hitFED);
523 (dynamic_cast<const EcalPreshowerGeometry*>(GeometryES_))->getClosestCellInPlane(
point, 2);
526 ES_fedId_[(3 - stripY.
zside()) / 2 - 1][stripY.
plane() - 1][stripY.
six() - 1][stripY.
siy() - 1];
529 LogDebug(
"SelectedElectronFEDListProducer")
530 <<
" ES hit plane Y (deiID) " << stripY.
rawId() <<
" six " << stripY.
six() <<
" siy "
531 << stripY.
siy() <<
" plane " << stripY.
plane() <<
" FED ID " << hitFED << std::endl;
534 if (!fedList_.empty()) {
535 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
536 fedList_.push_back(hitFED);
538 fedList_.push_back(hitFED);
544 if (dumpSelectedHCALFed_) {
546 for (; itHcalRecHit != hcalRecHitCollection->
end(); ++itHcalRecHit) {
549 static_cast<const HcalGeometry*>(GeometryCalo_->getSubdetectorGeometry(recHitId));
552 cellGeometry->getPosition(recHitId).eta(),
553 cellGeometry->getPosition(recHitId).phi());
554 if (
dR <= dRHcalRegion_) {
557 LogDebug(
"SelectedElectronFEDListProducer")
558 <<
" matched hcal recHit : HcalDetId " << recHitId <<
" HcalElectronicsId " << electronicId
559 <<
" dcc id " << electronicId.
dccid() <<
" spigot " << electronicId.
spigot() <<
" fiber channel "
565 if (!fedList_.empty()) {
566 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
567 fedList_.push_back(hitFED);
569 fedList_.push_back(hitFED);
576 if (!dumpAllTrackerFed_) {
578 if (dumpSelectedSiStripFed_) {
581 if (*itElectronCollFlag) {
588 for (uint32_t iCabling = 0; iCabling < SiStripCabling.size(); iCabling++) {
590 double dphi = fabs(
pos.second - phi);
592 dphi = 2 * acos(-1) - dphi;
594 if (
R -
sqrt(
pow(regionDimension_.first / 2, 2) +
pow(regionDimension_.second / 2, 2)) > dRStripRegion_)
605 regSubdetLayers[ilayer];
606 SiStripRegionCabling::ElementCabling::const_iterator itFedMap = fedVectorMap.begin();
607 for (; itFedMap != fedVectorMap.end(); itFedMap++) {
608 for (uint32_t op = 0; op < (itFedMap->second).
size(); op++) {
609 int hitFED = (itFedMap->second)[op].
fedId();
612 LogDebug(
"SelectedElectronFEDListProducer") <<
" SiStrip (FedID) " << hitFED << std::endl;
613 if (!fedList_.empty()) {
614 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
615 fedList_.push_back(hitFED);
617 fedList_.push_back(hitFED);
624 if (dumpSelectedSiPixelFed_) {
626 if (*itElectronCollFlag)
627 momentum =
electron.gsfTrack()->momentum();
629 momentum =
electron.track()->momentum();
630 PixelRegion region(momentum, dPhiPixelRegion_, dEtaPixelRegion_, maxZPixelRegion_);
634 std::vector<PixelModule>::const_iterator itUp, itDn;
635 if (lowerBound.Phi >= -
M_PI && upperBound.Phi <=
M_PI) {
636 itDn =
std::lower_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), lowerBound);
637 itUp =
std::upper_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), upperBound);
638 pixelFedDump(itDn, itUp,
region);
640 if (lowerBound.Phi < -
M_PI)
641 lowerBound.Phi = lowerBound.Phi + 2 *
M_PI;
643 itDn =
std::lower_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), lowerBound);
644 itUp =
std::upper_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), phi_p);
645 pixelFedDump(itDn, itUp,
region);
647 if (upperBound.Phi < -
M_PI)
648 upperBound.Phi = upperBound.Phi - 2 *
M_PI;
650 itDn =
std::lower_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), phi_m);
651 itUp =
std::upper_bound(pixelModuleVector_.begin(), pixelModuleVector_.end(), upperBound);
652 pixelFedDump(itDn, itUp,
region);
661 for (
unsigned int iFed = 0; iFed < addThisSelectedFEDs_.size(); iFed++) {
662 if (addThisSelectedFEDs_.at(iFed) == -1)
664 fedList_.push_back(addThisSelectedFEDs_.at(iFed));
668 auto streamFEDRawProduct = std::make_unique<FEDRawDataCollection>();
669 std::sort(fedList_.begin(), fedList_.end());
670 std::vector<uint32_t>::const_iterator itfedList = fedList_.begin();
671 for (; itfedList != fedList_.end(); ++itfedList) {
672 LogDebug(
"SelectedElectronFEDListProducer") <<
" fed point " << *itfedList <<
" ";
674 if (
data.size() > 0) {
675 FEDRawData& fedData = streamFEDRawProduct->FEDData(*itfedList);
683 if (!fedList_.empty())
687 template <
typename TEle,
typename TCand>
689 LogDebug(
"SelectedElectronFEDListProducer") <<
" End of the Job " << std::endl;
692 template <
typename TEle,
typename TCand>
694 std::vector<PixelModule>::const_iterator& itUp,
696 for (; itDn != itUp; ++itDn) {
697 float zmodule = itDn->z - ((itDn->x - beamSpotPosition_.x()) *
region.cosphi +
698 (itDn->y - beamSpotPosition_.y()) *
region.sinphi) *
702 int hitFED = itDn->Fed;
705 LogDebug(
"SelectedElectronFEDListProducer")
706 <<
" electron pixel hit " << itDn->DetId <<
" hitFED " << hitFED << std::endl;
707 if (!fedList_.empty()) {
708 if (
std::find(fedList_.begin(), fedList_.end(), hitFED) == fedList_.end())
709 fedList_.push_back(hitFED);
711 fedList_.push_back(hitFED);
717 template <
typename TEle,
typename TCand>
720 desc.add<vector<edm::InputTag>>(
"electronTags", {
edm::InputTag(
"hltEgammaGsfElectrons")});
721 desc.add<vector<edm::InputTag>>(
"recoEcalCandidateTags", {
edm::InputTag(
"hltL1EG25Ele27WP85GsfTrackIsoFilter")});
726 desc.add<vector<int>>(
"addThisSelectedFEDs", {812, 813});
727 desc.add<vector<int>>(
"isGsfElectronCollection", {
true});
729 desc.add<
bool>(
"dumpSelectedSiPixelFed",
true);
730 desc.add<
bool>(
"dumpSelectedSiStripFed",
true);
731 desc.add<
bool>(
"dumpSelectedEcalFed",
true);
732 desc.add<
bool>(
"dumpSelectedHCALFed",
true);
733 desc.add<
double>(
"dPhiPixelRegion", 0.3);
734 desc.add<
double>(
"dEtaPixelRegion", 0.3);
735 desc.add<
double>(
"dRStripRegion", 0.3);
736 desc.add<
double>(
"dRHcalRegion", 0.3);
737 desc.add<
double>(
"maxZPixelRegion", 24);
738 desc.add<
bool>(
"dumpAllTrackerFed",
false);
739 desc.add<
bool>(
"dumpAllEcalFed",
false);
740 desc.add<
bool>(
"dumpAllHcalFed",
false);