Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:26

0001 /**_________________________________________________________________
0002    class:   AlcaBeamSpotFromDB.cc
0003    package: RecoVertex/BeamSpotProducer
0004 
0005 
0006 
0007  author: Francisco Yumiceva, Fermilab (yumiceva@fnal.gov)
0008 
0009 
0010 ________________________________________________________________**/
0011 
0012 // C++ standard
0013 #include <string>
0014 // CMS
0015 #include "Calibration/TkAlCaRecoProducers/interface/AlcaBeamSpotFromDB.h"
0016 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 #include "FWCore/Framework/interface/MakerMacros.h"
0021 #include "FWCore/ServiceRegistry/interface/Service.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/IOVSyncValue.h"
0024 #include "FWCore/Framework/interface/LuminosityBlock.h"
0025 
0026 AlcaBeamSpotFromDB::AlcaBeamSpotFromDB(const edm::ParameterSet &iConfig)
0027     : beamSpotToken_(esConsumes<BeamSpotObjects, BeamSpotObjectsRcd, edm::Transition::EndLuminosityBlock>()) {
0028   produces<reco::BeamSpot, edm::Transition::EndLuminosityBlock>("alcaBeamSpot");
0029 }
0030 
0031 AlcaBeamSpotFromDB::~AlcaBeamSpotFromDB() {}
0032 
0033 //--------------------------------------------------------------------------------------------------
0034 void AlcaBeamSpotFromDB::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {}
0035 
0036 //--------------------------------------------------------------------------------------------------
0037 void AlcaBeamSpotFromDB::endLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) {
0038   // read DB object
0039   const BeamSpotObjects *spotDB = &iSetup.getData(beamSpotToken_);
0040 
0041   // translate from BeamSpotObjects to reco::BeamSpot
0042   reco::BeamSpot::Point apoint(spotDB->x(), spotDB->y(), spotDB->z());
0043 
0044   reco::BeamSpot::CovarianceMatrix matrix;
0045   for (int i = 0; i < 7; ++i) {
0046     for (int j = 0; j < 7; ++j) {
0047       matrix(i, j) = spotDB->covariance(i, j);
0048     }
0049   }
0050 
0051   reco::BeamSpot aSpot;
0052   // this assume beam width same in x and y
0053   aSpot = reco::BeamSpot(apoint, spotDB->sigmaZ(), spotDB->dxdz(), spotDB->dydz(), spotDB->beamWidthX(), matrix);
0054   aSpot.setBeamWidthY(spotDB->beamWidthY());
0055   aSpot.setEmittanceX(spotDB->emittanceX());
0056   aSpot.setEmittanceY(spotDB->emittanceY());
0057   aSpot.setbetaStar(spotDB->betaStar());
0058 
0059   if (spotDB->beamType() == 2) {
0060     aSpot.setType(reco::BeamSpot::Tracker);
0061   } else {
0062     aSpot.setType(reco::BeamSpot::Fake);
0063   }
0064 
0065   auto result = std::make_unique<reco::BeamSpot>();
0066   *result = aSpot;
0067   lumiSeg.put(std::move(result), std::string("alcaBeamSpot"));
0068 
0069   // std::cout << " for runs: " << iEvent.id().run() << " - " <<
0070   // iEvent.id().run() << std::endl;
0071   std::cout << aSpot << std::endl;
0072 }
0073 
0074 void AlcaBeamSpotFromDB::beginJob() {}
0075 
0076 void AlcaBeamSpotFromDB::endJob() {}
0077 
0078 // define this as a plug-in
0079 DEFINE_FWK_MODULE(AlcaBeamSpotFromDB);