CMS 3D CMS Logo

PPSLocalTrackLiteReAligner.cc
Go to the documentation of this file.
1 /****************************************************************************
2  *
3  * This is a part of TOTEM offline software.
4  * Authors:
5  * Jan Kašpar (jan.kaspar@gmail.com)
6  *
7  ****************************************************************************/
8 
14 
17 
20 
21 //----------------------------------------------------------------------------------------------------
22 
24 public:
26 
27  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
28 
30 
31 private:
33 
35 
37 };
38 
39 //----------------------------------------------------------------------------------------------------
40 
42  : inputTrackToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("inputTrackTag"))),
44  iConfig.getParameter<edm::ESInputTag>("alignmentTag"))),
45  outputTrackTag_(iConfig.getParameter<std::string>("outputTrackTag")) {
46  produces<CTPPSLocalTrackLiteCollection>(outputTrackTag_);
47 }
48 
49 //----------------------------------------------------------------------------------------------------
50 
53 
54  desc.add<edm::InputTag>("inputTrackTag", edm::InputTag("ctppsLocalTrackLiteProducer"))
55  ->setComment("tag of the input CTPPSLocalTrackLiteCollection");
56 
57  desc.add<edm::ESInputTag>("alignmentTag", edm::ESInputTag(""))->setComment("tag of the alignment data");
58 
59  desc.add<std::string>("outputTrackTag", "")->setComment("tag of the output CTPPSLocalTrackLiteCollection");
60 
61  descr.add("ppsLocalTrackLiteReAligner", desc);
62 }
63 
64 //----------------------------------------------------------------------------------------------------
65 
67  // get alignment corrections
68  auto const &alignment = iSetup.getData(alignmentToken_);
69 
70  // prepare output
71  auto output = std::make_unique<CTPPSLocalTrackLiteCollection>();
72 
73  // apply alignment corrections
74  auto const &inputTracks = iEvent.get(inputTrackToken_);
75  for (const auto &tr : inputTracks) {
76  auto it = alignment.getRPMap().find(tr.rpId());
77  if (it == alignment.getRPMap().end()) {
78  edm::LogError("PPS") << "Cannot find alignment correction for RP " << tr.rpId() << ". The track will be skipped.";
79  } else {
80  output->emplace_back(tr.rpId(),
81  tr.x() + it->second.getShX(),
82  tr.xUnc(),
83  tr.y() + it->second.getShY(),
84  tr.yUnc(),
85  tr.tx(),
86  tr.txUnc(),
87  tr.ty(),
88  tr.tyUnc(),
89  tr.chiSquaredOverNDF(),
90  tr.pixelTrackRecoInfo(),
91  tr.numberOfPointsUsedForFit(),
92  tr.time(),
93  tr.timeUnc());
94  }
95  }
96 
97  // save output to event
99 }
100 
101 //----------------------------------------------------------------------------------------------------
102 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &)
Log< level::Error, false > LogError
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > inputTrackToken_
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
PPSLocalTrackLiteReAligner(const edm::ParameterSet &)
Container for CTPPS RP alignment corrections. The corrections are stored on two levels - RP and senso...
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
HLT enums.
const edm::ESGetToken< CTPPSRPAlignmentCorrectionsData, RPRealAlignmentRecord > alignmentToken_
def move(src, dest)
Definition: eostools.py:511