Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:39

0001 /****************************************************************************
0002  *
0003  * This is a part of TOTEM offline software.
0004  * Authors:
0005  *   Jan Kašpar (jan.kaspar@gmail.com)
0006  *
0007  ****************************************************************************/
0008 
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/global/EDProducer.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 
0015 #include "DataFormats/CTPPSReco/interface/CTPPSLocalTrackLite.h"
0016 #include "DataFormats/CTPPSReco/interface/CTPPSLocalTrackLiteFwd.h"
0017 
0018 #include "CondFormats/AlignmentRecord/interface/RPRealAlignmentRecord.h"
0019 #include "CondFormats/PPSObjects/interface/CTPPSRPAlignmentCorrectionsData.h"
0020 
0021 //----------------------------------------------------------------------------------------------------
0022 
0023 class PPSLocalTrackLiteReAligner : public edm::global::EDProducer<> {
0024 public:
0025   explicit PPSLocalTrackLiteReAligner(const edm::ParameterSet &);
0026 
0027   void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
0028 
0029   static void fillDescriptions(edm::ConfigurationDescriptions &);
0030 
0031 private:
0032   const edm::EDGetTokenT<CTPPSLocalTrackLiteCollection> inputTrackToken_;
0033 
0034   const edm::ESGetToken<CTPPSRPAlignmentCorrectionsData, RPRealAlignmentRecord> alignmentToken_;
0035 
0036   const std::string outputTrackTag_;
0037 };
0038 
0039 //----------------------------------------------------------------------------------------------------
0040 
0041 PPSLocalTrackLiteReAligner::PPSLocalTrackLiteReAligner(const edm::ParameterSet &iConfig)
0042     : inputTrackToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("inputTrackTag"))),
0043       alignmentToken_(esConsumes<CTPPSRPAlignmentCorrectionsData, RPRealAlignmentRecord>(
0044           iConfig.getParameter<edm::ESInputTag>("alignmentTag"))),
0045       outputTrackTag_(iConfig.getParameter<std::string>("outputTrackTag")) {
0046   produces<CTPPSLocalTrackLiteCollection>(outputTrackTag_);
0047 }
0048 
0049 //----------------------------------------------------------------------------------------------------
0050 
0051 void PPSLocalTrackLiteReAligner::fillDescriptions(edm::ConfigurationDescriptions &descr) {
0052   edm::ParameterSetDescription desc;
0053 
0054   desc.add<edm::InputTag>("inputTrackTag", edm::InputTag("ctppsLocalTrackLiteProducer"))
0055       ->setComment("tag of the input CTPPSLocalTrackLiteCollection");
0056 
0057   desc.add<edm::ESInputTag>("alignmentTag", edm::ESInputTag(""))->setComment("tag of the alignment data");
0058 
0059   desc.add<std::string>("outputTrackTag", "")->setComment("tag of the output CTPPSLocalTrackLiteCollection");
0060 
0061   descr.add("ppsLocalTrackLiteReAligner", desc);
0062 }
0063 
0064 //----------------------------------------------------------------------------------------------------
0065 
0066 void PPSLocalTrackLiteReAligner::produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const {
0067   // get alignment corrections
0068   auto const &alignment = iSetup.getData(alignmentToken_);
0069 
0070   // prepare output
0071   auto output = std::make_unique<CTPPSLocalTrackLiteCollection>();
0072 
0073   // apply alignment corrections
0074   auto const &inputTracks = iEvent.get(inputTrackToken_);
0075   for (const auto &tr : inputTracks) {
0076     auto it = alignment.getRPMap().find(tr.rpId());
0077     if (it == alignment.getRPMap().end()) {
0078       edm::LogError("PPS") << "Cannot find alignment correction for RP " << tr.rpId() << ". The track will be skipped.";
0079     } else {
0080       output->emplace_back(tr.rpId(),
0081                            tr.x() + it->second.getShX(),
0082                            tr.xUnc(),
0083                            tr.y() + it->second.getShY(),
0084                            tr.yUnc(),
0085                            tr.tx(),
0086                            tr.txUnc(),
0087                            tr.ty(),
0088                            tr.tyUnc(),
0089                            tr.chiSquaredOverNDF(),
0090                            tr.pixelTrackRecoInfo(),
0091                            tr.numberOfPointsUsedForFit(),
0092                            tr.time(),
0093                            tr.timeUnc());
0094     }
0095   }
0096 
0097   // save output to event
0098   iEvent.put(std::move(output), outputTrackTag_);
0099 }
0100 
0101 //----------------------------------------------------------------------------------------------------
0102 
0103 DEFINE_FWK_MODULE(PPSLocalTrackLiteReAligner);