CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/SimCalorimetry/HcalTestBeam/src/HcalTBDigiProducer.cc

Go to the documentation of this file.
00001 #include "SimCalorimetry/HcalTestBeam/interface/HcalTBDigiProducer.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 #include "SimCalorimetry/CaloSimAlgos/interface/CaloTDigitizer.h"
00007 #include "SimCalorimetry/CaloSimAlgos/interface/CaloShapeIntegrator.h"
00008 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00009 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00010 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00011 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00012 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00013 #include "SimDataFormats/EcalTestBeam/interface/PEcalTBInfo.h"
00014 #include "SimGeneral/MixingModule/interface/PileUpEventPrincipal.h"
00015 
00016 HcalTBDigiProducer::HcalTBDigiProducer(const edm::ParameterSet& ps, edm::EDProducer& mixMod) :
00017   theParameterMap(new HcalTBSimParameterMap(ps)), 
00018   theHcalShape(new HcalShape()),
00019   theHcalIntegratedShape(new CaloShapeIntegrator(theHcalShape)),
00020   theHBHEResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
00021   theHOResponse(new CaloHitResponse(theParameterMap, theHcalIntegratedShape)),
00022   theAmplifier(0), theCoderFactory(0), theElectronicsSim(0), 
00023   theHitCorrection(0), theHBHEDigitizer(0), theHODigitizer(0), theHBHEHits(),
00024   theHOHits(), thisPhaseShift(0) {
00025   std::string const instance("simHcalDigis");
00026   mixMod.produces<HBHEDigiCollection>(instance);
00027   mixMod.produces<HODigiCollection>(instance);
00028 
00029   DetId detId(DetId::Hcal, 1);
00030   bool syncPhase = (theParameterMap->simParameters(detId)).syncPhase();
00031   doPhaseShift   = !syncPhase;
00032 
00033   theHBHEResponse->setHitFilter(&theHBHEHitFilter);
00034   theHOResponse->setHitFilter(&theHOHitFilter);
00035 
00036   bool doTimeSlew = ps.getParameter<bool>("doTimeSlew");
00037   if(doTimeSlew) {
00038     // no time slewing for HF
00039     theHitCorrection = new HcalHitCorrection(theParameterMap);
00040     theHBHEResponse->setHitCorrection(theHitCorrection);
00041     theHOResponse->setHitCorrection(theHitCorrection);
00042   }
00043 
00044   bool doNoise = ps.getParameter<bool>("doNoise");
00045   theAmplifier = new HcalAmplifier(theParameterMap, doNoise);
00046   theCoderFactory = new HcalCoderFactory(HcalCoderFactory::DB);
00047   theElectronicsSim = new HcalElectronicsSim(theAmplifier, theCoderFactory);
00048 
00049   theHBHEDigitizer = new HBHEDigitizer(theHBHEResponse, theElectronicsSim, doNoise);
00050   theHODigitizer = new HODigitizer(theHOResponse, theElectronicsSim, doNoise);
00051 
00052   tunePhaseShift =  ps.getUntrackedParameter<double>("tunePhaseShiftTB", 1.);
00053   ecalTBInfoLabel = ps.getUntrackedParameter<std::string>("EcalTBInfoLabel","SimEcalTBG4Object");
00054   edm::LogInfo("HcalSim") << "HcalTBDigiProducer initialized with doNoise = "
00055                           << doNoise << ", doTimeSlew = " << doTimeSlew
00056                           << " and doPhaseShift = " << doPhaseShift
00057                           << " tunePhasShift = " << tunePhaseShift;
00058 
00059 }
00060 
00061 HcalTBDigiProducer::~HcalTBDigiProducer() {
00062 
00063   if (theHBHEDigitizer)       delete theHBHEDigitizer;
00064   if (theHODigitizer)         delete theHODigitizer;
00065   if (theParameterMap)        delete theParameterMap;
00066   if (theHcalShape)           delete theHcalShape;
00067   if (theHcalIntegratedShape) delete theHcalIntegratedShape;
00068   if (theHBHEResponse)        delete theHBHEResponse;
00069   if (theHOResponse)          delete theHOResponse;
00070   if (theElectronicsSim)      delete theElectronicsSim;
00071   if (theAmplifier)           delete theAmplifier;
00072   if (theCoderFactory)        delete theCoderFactory;
00073   if (theHitCorrection)       delete theHitCorrection;
00074 }
00075 
00076 
00077 void HcalTBDigiProducer::initializeEvent(edm::Event const& e, edm::EventSetup const& eventSetup) {
00078   // get the appropriate gains, noises, & widths for this event
00079   edm::ESHandle<HcalDbService> conditions;
00080   eventSetup.get<HcalDbRecord>().get(conditions);
00081   theAmplifier->setDbService(conditions.product());
00082   theCoderFactory->setDbService(conditions.product());
00083 
00084   // get the correct geometry
00085   checkGeometry(eventSetup);
00086 
00087   theHBHEHits.clear();
00088   theHOHits.clear();
00089   if (doPhaseShift) {
00090     
00091     edm::Handle<PEcalTBInfo> theEcalTBInfo;
00092     e.getByLabel(ecalTBInfoLabel,theEcalTBInfo);
00093     thisPhaseShift = theEcalTBInfo->phaseShift();
00094 
00095     DetId detIdHB(DetId::Hcal, 1);
00096     setPhaseShift(detIdHB);
00097     DetId detIdHO(DetId::Hcal, 3);
00098     setPhaseShift(detIdHO);
00099   }
00100 
00101   theHBHEDigitizer->initializeHits();
00102   theHODigitizer->initializeHits();
00103 }
00104 
00105 void HcalTBDigiProducer::accumulateCaloHits(edm::Handle<std::vector<PCaloHit> > const& hcalHandle, int bunchCrossing) {
00106 
00107   LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate trying to get SimHit";
00108 
00109   if(hcalHandle.isValid()) {
00110     std::vector<PCaloHit> hits = *hcalHandle.product();
00111     if(theHitCorrection != 0) {
00112       theHitCorrection->fillChargeSums(hits);
00113     }
00114     LogDebug("HcalSim") << "HcalTBDigiProducer::accumulate Hits corrected";
00115     theHBHEDigitizer->add(hits, bunchCrossing);
00116     theHODigitizer->add(hits, bunchCrossing);
00117   }
00118 }
00119 
00120 void HcalTBDigiProducer::accumulate(edm::Event const& e, edm::EventSetup const&) {
00121   // Step A: Get Inputs, and accumulate digis
00122 
00123   edm::InputTag hcalTag("g4SimHits", "HcalHits");
00124   edm::Handle<std::vector<PCaloHit> > hcalHandle;
00125   e.getByLabel(hcalTag, hcalHandle);
00126 
00127   accumulateCaloHits(hcalHandle, 0);
00128 }
00129 
00130 void HcalTBDigiProducer::accumulate(PileUpEventPrincipal const& e, edm::EventSetup const&) {
00131   // Step A: Get Inputs, and accumulate digis
00132 
00133   edm::InputTag hcalTag("g4SimHits", "HcalHits");
00134   edm::Handle<std::vector<PCaloHit> > hcalHandle;
00135   e.getByLabel(hcalTag, hcalHandle);
00136 
00137   accumulateCaloHits(hcalHandle, e.bunchCrossing());
00138 }
00139 
00140 void HcalTBDigiProducer::finalizeEvent(edm::Event& e, const edm::EventSetup& eventSetup) {
00141   // Step B: Create empty output
00142   std::auto_ptr<HBHEDigiCollection> hbheResult(new HBHEDigiCollection());
00143   std::auto_ptr<HODigiCollection> hoResult(new HODigiCollection());
00144   LogDebug("HcalSim") << "HcalTBDigiProducer::produce Empty collection created";
00145   // Step C: Invoke the algorithm, getting back outputs.
00146   theHBHEDigitizer->run(*hbheResult);
00147   edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HBHE digis : " << hbheResult->size();
00148   theHODigitizer->run(*hoResult);
00149   edm::LogInfo("HcalSim") << "HcalTBDigiProducer: HO digis   : " << hoResult->size();
00150 
00151   // Step D: Put outputs into event
00152   std::string const instance("simHcalDigis");
00153   e.put(hbheResult, instance);
00154   e.put(hoResult, instance);
00155 
00156 }
00157 
00158 void HcalTBDigiProducer::sortHits(const edm::PCaloHitContainer & hits) {
00159 
00160   for (edm::PCaloHitContainer::const_iterator hitItr = hits.begin();
00161        hitItr != hits.end(); ++hitItr) {
00162     HcalSubdetector subdet = HcalDetId(hitItr->id()).subdet();
00163     if(subdet == HcalBarrel || subdet == HcalEndcap) {
00164       theHBHEHits.push_back(*hitItr);
00165     } else if(subdet == HcalOuter) {
00166       theHOHits.push_back(*hitItr);
00167     } else {
00168       edm::LogError("HcalSim") << "Bad HcalHit subdetector " << subdet;
00169     }
00170   }
00171 }
00172 
00173 void HcalTBDigiProducer::checkGeometry(const edm::EventSetup & eventSetup) {
00174 
00175   // TODO find a way to avoid doing this every event
00176   edm::ESHandle<CaloGeometry> geometry;
00177   eventSetup.get<CaloGeometryRecord>().get(geometry);
00178  
00179   const CaloGeometry * pGeometry = &*geometry;
00180 
00181   // see if we need to update
00182   if(pGeometry != theGeometry) {
00183     theGeometry = pGeometry;
00184     updateGeometry();
00185   }
00186 }
00187 
00188 void HcalTBDigiProducer::updateGeometry() {
00189 
00190   theHBHEResponse->setGeometry(theGeometry);
00191   theHOResponse->setGeometry(theGeometry);
00192 
00193   // Get cells for HB and HE
00194   hbheCells.clear();
00195   hbheCells = theGeometry->getValidDetIds(DetId::Hcal, HcalBarrel);
00196   std::vector<DetId> heCells = theGeometry->getValidDetIds(DetId::Hcal, HcalEndcap);
00197   // combine HB & HE
00198   hbheCells.insert(hbheCells.end(), heCells.begin(), heCells.end());
00199 
00200   // Get cells for HO
00201   hoCells.clear();
00202   hoCells = theGeometry->getValidDetIds(DetId::Hcal, HcalOuter);
00203 
00204   edm::LogInfo("HcalSim") << "HcalTBDigiProducer update Geometry with "
00205                           << hbheCells.size() << " cells in HB/HE and "
00206                           << hoCells.size() << " cells in HO";
00207 
00208   theHBHEDigitizer->setDetIds(hbheCells);
00209   LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HB/HE";
00210   theHODigitizer->setDetIds(hoCells);
00211   LogDebug("HcalSim") << "HcalTBDigiProducer: Set DetID's for HO";
00212 }
00213 
00214 void HcalTBDigiProducer::setPhaseShift(const DetId & detId) {
00215 
00216   const CaloSimParameters & parameters = theParameterMap->simParameters(detId);
00217   if ( !parameters.syncPhase() ) {
00218     int    myDet          = detId.subdetId();
00219     double passPhaseShift = thisPhaseShift + tunePhaseShift;
00220     if (myDet <= 2) {
00221       theHBHEResponse->setPhaseShift(passPhaseShift);
00222     } else {
00223       theHOResponse->setPhaseShift(passPhaseShift);
00224     }
00225   }
00226 }