Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79
/**_________________________________________________________________
   class:   AlcaBeamSpotFromDB.cc
   package: RecoVertex/BeamSpotProducer



 author: Francisco Yumiceva, Fermilab (yumiceva@fnal.gov)


________________________________________________________________**/

// C++ standard
#include <string>
// CMS
#include "Calibration/TkAlCaRecoProducers/interface/AlcaBeamSpotFromDB.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"

#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/IOVSyncValue.h"
#include "FWCore/Framework/interface/LuminosityBlock.h"

AlcaBeamSpotFromDB::AlcaBeamSpotFromDB(const edm::ParameterSet &iConfig)
    : beamSpotToken_(esConsumes<BeamSpotObjects, BeamSpotObjectsRcd, edm::Transition::EndLuminosityBlock>()) {
  produces<reco::BeamSpot, edm::Transition::EndLuminosityBlock>("alcaBeamSpot");
}

AlcaBeamSpotFromDB::~AlcaBeamSpotFromDB() {}

//--------------------------------------------------------------------------------------------------
void AlcaBeamSpotFromDB::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {}

//--------------------------------------------------------------------------------------------------
void AlcaBeamSpotFromDB::endLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) {
  // read DB object
  const BeamSpotObjects *spotDB = &iSetup.getData(beamSpotToken_);

  // translate from BeamSpotObjects to reco::BeamSpot
  reco::BeamSpot::Point apoint(spotDB->x(), spotDB->y(), spotDB->z());

  reco::BeamSpot::CovarianceMatrix matrix;
  for (int i = 0; i < 7; ++i) {
    for (int j = 0; j < 7; ++j) {
      matrix(i, j) = spotDB->covariance(i, j);
    }
  }

  reco::BeamSpot aSpot;
  // this assume beam width same in x and y
  aSpot = reco::BeamSpot(apoint, spotDB->sigmaZ(), spotDB->dxdz(), spotDB->dydz(), spotDB->beamWidthX(), matrix);
  aSpot.setBeamWidthY(spotDB->beamWidthY());
  aSpot.setEmittanceX(spotDB->emittanceX());
  aSpot.setEmittanceY(spotDB->emittanceY());
  aSpot.setbetaStar(spotDB->betaStar());

  if (spotDB->beamType() == 2) {
    aSpot.setType(reco::BeamSpot::Tracker);
  } else {
    aSpot.setType(reco::BeamSpot::Fake);
  }

  auto result = std::make_unique<reco::BeamSpot>();
  *result = aSpot;
  lumiSeg.put(std::move(result), std::string("alcaBeamSpot"));

  // std::cout << " for runs: " << iEvent.id().run() << " - " <<
  // iEvent.id().run() << std::endl;
  std::cout << aSpot << std::endl;
}

void AlcaBeamSpotFromDB::beginJob() {}

void AlcaBeamSpotFromDB::endJob() {}

// define this as a plug-in
DEFINE_FWK_MODULE(AlcaBeamSpotFromDB);