CMS 3D CMS Logo

TkAlV0sAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: TkAlV0sAnalyzer
5 //
6 /*
7  *\class TkAlV0sAnalyzer TkAlV0sAnalyzer.cc Alignment/TkAlV0sAnalyzer/plugins/TkAlV0sAnalyzer.cc
8 
9  Description: [one line class summary]
10 
11  Implementation:
12  [Notes on implementation]
13 */
14 //
15 // Original Author: Marco Musich
16 // Created: Thu, 14 Dec 2023 15:10:34 GMT
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
38 
39 #include "TLorentzVector.h"
40 
41 //
42 // class declaration
43 //
44 
45 // If the analyzer does not use TFileService, please remove
46 // the template argument to the base class so the class inherits
47 // from edm::one::EDAnalyzer<>
48 // This will improve performance in multithreaded jobs.
49 
51 
52 class TkAlV0sAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
53 public:
54  explicit TkAlV0sAnalyzer(const edm::ParameterSet&);
55  ~TkAlV0sAnalyzer() override = default;
56 
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58 
59 private:
60  template <typename T, typename... Args>
61  T* book(const Args&... args) const;
62  void beginJob() override;
63  void analyze(const edm::Event&, const edm::EventSetup&) override;
64 
65  // ----------member data ---------------------------
66  const edm::EDGetTokenT<TrackCollection> tracksToken_; //used to select what tracks to read from configuration file
69 
71  TH1F* h_V0Mass;
72 };
73 
74 static constexpr double piMass2 = 0.13957018 * 0.13957018;
75 
76 //
77 // constructors and destructor
78 //
80  : tracksToken_(consumes<TrackCollection>(iConfig.getUntrackedParameter<edm::InputTag>("tracks"))),
81  vccToken_(consumes<reco::VertexCompositeCandidateCollection>(
82  iConfig.getParameter<edm::InputTag>("vertexCompositeCandidates"))) {
83  usesResource(TFileService::kSharedResource);
84 }
85 
86 template <typename T, typename... Args>
87 T* TkAlV0sAnalyzer::book(const Args&... args) const {
88  T* t = fs_->make<T>(args...);
89  return t;
90 }
91 
92 //
93 // member functions
94 //
95 // ------------ method called for each event ------------
97  using namespace edm;
98 
99  std::vector<const reco::Track*> myTracks;
100 
102  iEvent.getByToken(vccToken_, vccHandle);
103 
104  if (vccHandle->empty())
105  return;
106 
108 
109  for (const auto& track : iEvent.get(tracksToken_)) {
110  myTracks.emplace_back(&track);
111  }
112 
113  // exclude multiple candidates
114  if (myTracks.size() != 2)
115  return;
116 
117  for (const auto& v0 : v0s) {
118  float mass = v0.mass();
119  h_V0Mass->Fill(mass);
120 
121  for (size_t i = 0; i < v0.numberOfDaughters(); ++i) {
122  //LogPrint("AlignmentTrackFromVertexCompositeCandidateSelector") << "daughter: " << i << std::endl;
123  const reco::Candidate* daughter = v0.daughter(i);
124  const reco::RecoChargedCandidate* chargedDaughter = dynamic_cast<const reco::RecoChargedCandidate*>(daughter);
125  if (chargedDaughter) {
126  //LogPrint("AlignmentTrackFromVertexCompositeCandidateSelector") << "charged daughter: " << i << std::endl;
127  const reco::TrackRef trackRef = chargedDaughter->track();
128  if (trackRef.isNonnull()) {
129  // LogPrint("AlignmentTrackFromVertexCompositeCandidateSelector")
130  // << "charged daughter has non-null trackref: " << i << std::endl;
131  }
132  }
133  }
134  }
135 
136  const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
137  const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
138 
139  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + piMass2));
140  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + piMass2));
141 
142  const auto& V0p4 = p4_tplus + p4_tminus;
143  float track_invMass = V0p4.M();
144  h_diTrackMass->Fill(track_invMass);
145 }
146 
148  h_diTrackMass = book<TH1F>("diTrackMass", "V0 mass from tracks in Event", 100, 0.400, 0.600);
149  h_V0Mass = book<TH1F>("V0kMass", "Reconstructed V0 mass in Event", 100, 0.400, 0.600);
150 }
151 
152 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
155  desc.add<edm::InputTag>("vertexCompositeCandidates", edm::InputTag("generalV0Candidates:Kshort"));
156  desc.addUntracked<edm::InputTag>("tracks", edm::InputTag("ALCARECOTkAlKShortTracks"));
157  descriptions.addWithDefaultLabel(desc);
158 }
159 
160 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
T const * product() const
Definition: Handle.h:70
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
const edm::EDGetTokenT< TrackCollection > tracksToken_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::Service< TFileService > fs_
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > vccToken_
T * book(const Args &... args) const
T sqrt(T t)
Definition: SSEVec.h:19
reco::TrackRef track() const override
reference to a track
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~TkAlV0sAnalyzer() override=default
static constexpr double piMass2
fixed size matrix
HLT enums.
TkAlV0sAnalyzer(const edm::ParameterSet &)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void beginJob() override
long double T
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)