Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:31:05

0001 #include "TauAnalysis/MCEmbeddingTools/plugins/CaloCleaner.h"
0002 
0003 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
0004 #include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
0005 #include "DataFormats/HcalRecHit/interface/HFRecHit.h"
0006 #include "DataFormats/HcalRecHit/interface/HORecHit.h"
0007 #include "DataFormats/HcalRecHit/interface/ZDCRecHit.h"
0008 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
0009 
0010 typedef CaloCleaner<EcalRecHit> EcalRecHitColCleaner;
0011 typedef CaloCleaner<HBHERecHit> HBHERecHitColCleaner;
0012 typedef CaloCleaner<HFRecHit> HFRecHitColCleaner;
0013 typedef CaloCleaner<HORecHit> HORecHitColCleaner;
0014 typedef CaloCleaner<CastorRecHit> CastorRecHitColCleaner;
0015 typedef CaloCleaner<ZDCRecHit> ZDCRecHitColCleaner;
0016 
0017 //-------------------------------------------------------------------------------
0018 // define 'buildRecHit' functions used for different types of recHits
0019 //-------------------------------------------------------------------------------
0020 
0021 template <typename T>
0022 void CaloCleaner<T>::fill_correction_map(TrackDetMatchInfo *, std::map<uint32_t, float> *) {
0023   assert(0);  // CV: make sure general function never gets called;
0024               //     always use template specializations
0025 }
0026 
0027 template <>
0028 void CaloCleaner<EcalRecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0029   if (is_preshower_) {
0030     for (std::vector<DetId>::const_iterator detId = info->crossedPreshowerIds.begin();
0031          detId != info->crossedPreshowerIds.end();
0032          ++detId) {
0033       (*cor_map)[detId->rawId()] = 9999999;  // just remove all energy (Below 0 is not possible)
0034     }
0035   } else {
0036     for (std::vector<const EcalRecHit *>::const_iterator hit = info->crossedEcalRecHits.begin();
0037          hit != info->crossedEcalRecHits.end();
0038          hit++) {
0039       //    (*cor_map) [(*hit)->detid().rawId()] +=(*hit)->energy();
0040       (*cor_map)[(*hit)->detid().rawId()] = (*hit)->energy();
0041     }
0042   }
0043 }
0044 
0045 template <>
0046 void CaloCleaner<HBHERecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0047   for (std::vector<const HBHERecHit *>::const_iterator hit = info->crossedHcalRecHits.begin();
0048        hit != info->crossedHcalRecHits.end();
0049        hit++) {
0050     (*cor_map)[(*hit)->detid().rawId()] = (*hit)->energy();
0051   }
0052 }
0053 
0054 template <>
0055 void CaloCleaner<HORecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0056   for (std::vector<const HORecHit *>::const_iterator hit = info->crossedHORecHits.begin();
0057        hit != info->crossedHORecHits.end();
0058        hit++) {
0059     (*cor_map)[(*hit)->detid().rawId()] = (*hit)->energy();
0060   }
0061 }
0062 
0063 template <>
0064 void CaloCleaner<HFRecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0065   return;  // No corrections for HF
0066 }
0067 
0068 template <>
0069 void CaloCleaner<CastorRecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0070   return;  // No corrections for Castor
0071 }
0072 
0073 template <>
0074 void CaloCleaner<ZDCRecHit>::fill_correction_map(TrackDetMatchInfo *info, std::map<uint32_t, float> *cor_map) {
0075   return;  // No corrections for Castor
0076 }
0077 
0078 DEFINE_FWK_MODULE(EcalRecHitColCleaner);
0079 DEFINE_FWK_MODULE(HBHERecHitColCleaner);
0080 DEFINE_FWK_MODULE(HORecHitColCleaner);
0081 // no  need for cleaning outside of tracker, so just a copy of the old collection
0082 DEFINE_FWK_MODULE(HFRecHitColCleaner);
0083 DEFINE_FWK_MODULE(CastorRecHitColCleaner);
0084 DEFINE_FWK_MODULE(ZDCRecHitColCleaner);