CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/SimCalorimetry/HcalTestBeam/src/HcalTBDigiProducer.cc

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