ATLAS Offline Software
mergingDJRs.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
5 #include "UserHooksUtils.h"
6 #include "UserSetting.h"
7 
8 // ROOT, for saving file.
9 #include "TFile.h"
10 
11 // ROOT, for histogramming.
12 #include "TF1.h"
13 #include "TH1.h"
14 
15 #include <iostream>
16 
17 namespace Pythia8 {
18 
19 // Class to compute the DJRs
20 class mergingDJRs : public UserHooks {
21 
22 private:
23  // Settings for SlowJet jet finder, with kT clustering
24  int m_power = 1;
25  double m_etaMax = 10.;
26  double m_radius = 1.0;
27  double m_pTjetMin = 10.;
28  Pythia8::SlowJetHook *m_sjHookPtrIn = 0;
29  bool m_useFJcoreIn = false; // SHOULD BE FALSE TO ALLOW THE SUB PROCESS BEFORE
30  // CLUSTERING THE PAARTICLES
31  bool m_useStandardRin = true;
32  int m_nSel = 2; // Exclude neutrinos (and other invisible) from study.
33  int m_massSetIn = 1;
34  std::unique_ptr<Pythia8::SlowJet> m_slowJet;
35 
36  // ROOT output files
37  std::unique_ptr<TH1F> m_HistDjr, m_HistDjr2;
38  std::unique_ptr<TFile> m_outFile;
39 
40  vector<double> m_result, m_DJR;
42 
43 public:
45  // Slowjet pointer
46  m_slowJet = std::make_unique<Pythia8::SlowJet>(
49 
50  // ROOT histograms and the output file where to save them
51  m_HistDjr = std::make_unique<TH1F>("HistDjr", "The first DJR", 100, 0.0, 3.0);
52  m_HistDjr2 =
53  std::make_unique<TH1F>("HistDjr2", "The second DJR", 100, 0.0, 3.0);
54  m_outFile = std::make_unique<TFile>("hist-DJR.root", "RECREATE");
55 
56  std::cout << "**********************************************************"
57  << std::endl;
58  std::cout << "* *"
59  << std::endl;
60  std::cout << "* the jet merging userhook CKKWL DJRS is working *"
61  << std::endl;
62  std::cout << "* *"
63  << std::endl;
64  std::cout << "**********************************************************"
65  << std::endl;
66  }
67 
68  // Destructor prints histogram
70  m_HistDjr->Write();
71  m_HistDjr2->Write();
72  }
73 
74  // Parton level vetoing
75  virtual bool canVetoPartonLevel() override { return true; }
76  virtual bool doVetoPartonLevel(const Event &event) override;
77 
78  // Function to compute the DJRs
79  virtual void getDJR(const Event &event);
80 
81  // To initialize the variables
82  virtual bool initAfterBeams() override {
83 
84  m_workEventJet.init("(workEventJet)", particleDataPtr);
85 
86  return true;
87  }
88 };
89 
90 // parton level veto (before beam remnants and resonance showers)
92 
93  // subEvent method extract a list of the current partons list and save the
94  // output in workEvent.
95  subEvent(event);
96  m_workEventJet = workEvent;
97 
98  // The selected particles to pass to slowjet
99  for (int j = 0; j < m_workEventJet.size(); ++j) {
100  if (!(m_workEventJet[j].isFinal()) || m_workEventJet[j].isLepton() ||
101  m_workEventJet[j].id() == 23 || std::abs(m_workEventJet[j].id()) == 24 ||
102  m_workEventJet[j].id() == 22) {
103  m_workEventJet[j].statusNeg();
104  continue;
105  }
106  }
107 
108  // slowjet analyze the events
109  if (not m_slowJet->setup(m_workEventJet)){
110  std::cout<< "setup failed in mergingDJRs::doVetoPartonLevel\n";
111  return false;
112  }
113 
114  // Call getDJR and store the DJRs vector
116 
117  // if we reached here then no veto
118  return false;
119 }
120 
121 // Compute DJR vector
122 inline void mergingDJRs::getDJR(const Event &event) {
123 
124  // setup slowjet pointer
125  if (not m_slowJet->setup(event)){
126  std::cout<< "setup failed in mergingDJRs::getDJR\n";
127  return;
128  };
129 
130  // Clear members.
131  m_DJR.clear();
132  m_result.clear();
133 
134  while (m_slowJet->sizeAll() - m_slowJet->sizeJet() > 0) {
135  // Save the next clustering scale.
136  m_result.push_back(sqrt(m_slowJet->dNext()));
137  // Perform step.
138  m_slowJet->doStep();
139  }
140 
141  // Save clustering scales in reverse order
142  for (int i = m_result.size() - 1; i >= 0; --i) {
143  m_DJR.push_back(m_result[i]);
144  }
145 
146  // Fill the histogram and Normalize them
147  double eventWeight = infoPtr->mergingWeight() * infoPtr->weight();
148 
149  if (m_DJR.size() > 0) {
150  m_HistDjr->Fill(log10(m_DJR[0]), eventWeight);
151  m_HistDjr2->Fill(log10(m_DJR[1]), eventWeight);
152  }
153 }
154 
156  Ckkwl("mergingDJRs");
157 
158 } // namespace Pythia8
Pythia8::mergingDJRs::m_useStandardRin
bool m_useStandardRin
Definition: mergingDJRs.cxx:31
Pythia8::mergingDJRs::m_HistDjr2
std::unique_ptr< TH1F > m_HistDjr2
Definition: mergingDJRs.cxx:37
Pythia8::mergingDJRs::doVetoPartonLevel
virtual bool doVetoPartonLevel(const Event &event) override
Definition: mergingDJRs.cxx:91
Event
Definition: trigbs_orderedMerge.cxx:42
Pythia8::mergingDJRs::m_DJR
vector< double > m_DJR
Definition: mergingDJRs.cxx:40
Pythia8::mergingDJRs::m_nSel
int m_nSel
Definition: mergingDJRs.cxx:32
Pythia8::mergingDJRs::m_etaMax
double m_etaMax
Definition: mergingDJRs.cxx:25
Pythia8::mergingDJRs::m_sjHookPtrIn
Pythia8::SlowJetHook * m_sjHookPtrIn
Definition: mergingDJRs.cxx:28
Pythia8::mergingDJRs::m_slowJet
std::unique_ptr< Pythia8::SlowJet > m_slowJet
Definition: mergingDJRs.cxx:34
Pythia8::mergingDJRs::m_radius
double m_radius
Definition: mergingDJRs.cxx:26
Pythia8::mergingDJRs::m_massSetIn
int m_massSetIn
Definition: mergingDJRs.cxx:33
Pythia8_UserHooks::UserHooksFactory::Creator
Definition: UserHooksFactory.h:54
UserHooksFactory.h
Pythia8
Author: James Monk (jmonk@cern.ch)
Definition: IPythia8Custom.h:13
Pythia8::mergingDJRs::m_HistDjr
std::unique_ptr< TH1F > m_HistDjr
Definition: mergingDJRs.cxx:37
Pythia8::mergingDJRs::m_outFile
std::unique_ptr< TFile > m_outFile
Definition: mergingDJRs.cxx:38
event
POOL::TEvent event(POOL::TEvent::kClassAccess)
lumiFormat.i
int i
Definition: lumiFormat.py:92
Pythia8::mergingDJRs::m_pTjetMin
double m_pTjetMin
Definition: mergingDJRs.cxx:27
Pythia8::mergingDJRs
Definition: mergingDJRs.cxx:20
Pythia8::mergingDJRs::m_workEventJet
Event m_workEventJet
Definition: mergingDJRs.cxx:41
Pythia8::mergingDJRs::getDJR
virtual void getDJR(const Event &event)
Definition: mergingDJRs.cxx:122
UserSetting.h
Pythia8::mergingDJRs::m_useFJcoreIn
bool m_useFJcoreIn
Definition: mergingDJRs.cxx:29
Pythia8::mergingDJRs::m_power
int m_power
Definition: mergingDJRs.cxx:24
Pythia8::mergingDJRs::m_result
vector< double > m_result
Definition: mergingDJRs.cxx:40
Pythia8::Ckkwl
Pythia8_UserHooks::UserHooksFactory::Creator< Pythia8::mergingDJRs > Ckkwl("mergingDJRs")
UserHooksUtils.h
Pythia8::mergingDJRs::canVetoPartonLevel
virtual bool canVetoPartonLevel() override
Definition: mergingDJRs.cxx:75
isLepton
bool isLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:130
Pythia8::mergingDJRs::initAfterBeams
virtual bool initAfterBeams() override
Definition: mergingDJRs.cxx:82
Pythia8::mergingDJRs::mergingDJRs
mergingDJRs()
Definition: mergingDJRs.cxx:44
Pythia8::mergingDJRs::~mergingDJRs
~mergingDJRs()
Definition: mergingDJRs.cxx:69