Project
Loading...
Searching...
No Matches
GeneratorPythia8.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
14#ifndef ALICEO2_EVENTGEN_GENERATORPYTHIA8_H_
15#define ALICEO2_EVENTGEN_GENERATORPYTHIA8_H_
16
18#include "Pythia8/Pythia.h"
19#include <functional>
21
22namespace o2
23{
24namespace eventgen
25{
26
27/*****************************************************************/
28/*****************************************************************/
87{
88
89 public:
95 GeneratorPythia8(const Char_t* name, const Char_t* title = "ALICEo2 Pythia8 Generator");
97 ~GeneratorPythia8() override = default;
98
102 Bool_t Init() override;
104 Bool_t generateEvent() override;
106 Bool_t importParticles() override { return importParticles(mPythia.event); };
112 void setConfig(std::string val) { mConfig = val; };
116 void setHooksFileName(std::string val) { mHooksFileName = val; };
119 void setHooksFuncName(std::string val) { mHooksFuncName = val; };
122 void setUserHooks(Pythia8::UserHooks* hooks);
128 bool readString(std::string val) { return mPythia.readString(val, true); };
130 bool readFile(std::string val) { return mPythia.readFile(val, true); };
137 void getNcoll(int& nColl)
138 {
139 getNcoll(mPythia.info, nColl);
140 };
143 void getNpart(int& nPart)
144 {
145 getNpart(mPythia.info, nPart);
146 };
150 void getNpart(int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg)
151 {
152 getNpart(mPythia.info, nProtonProj, nNeutronProj, nProtonTarg, nNeutronTarg);
153 };
155 void getNremn(int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg)
156 {
157 getNremn(mPythia.event, nProtonProj, nNeutronProj, nProtonTarg, nNeutronTarg);
158 };
162 void getNfreeSpec(int& nFreenProj, int& nFreepProj, int& nFreenTarg, int& nFreepTarg)
163 {
164 getNfreeSpec(mPythia.info, nFreenProj, nFreepProj, nFreenTarg, nFreepTarg);
165 };
166
167 typedef std::function<bool(const Pythia8::Particle&)> UserFilterFcn;
168
173 bool setInitialSeed(long seed);
174
175 protected:
180
187 typedef std::vector<int> (*GetRelatives)(const Pythia8::Particle&);
189 typedef void (*SetRelatives)(Pythia8::Particle&, int, int);
191 typedef std::pair<int, int> (*FirstLastRelative)(const Pythia8::Particle&);
201 void updateHeader(o2::dataformats::MCEventHeader* eventHeader) override;
207 Bool_t importParticles(Pythia8::Event& event);
213 void pruneEvent(Pythia8::Event& event, Select select);
227 void investigateRelatives(Pythia8::Event& event, // Event
228 const std::vector<int>& old2New, // Map from old to new idx
229 size_t index, // Current particle
230 std::vector<bool>& done, // cache flag
231 GetRelatives getter, // get relatives
232 SetRelatives setter, // set relatives
233 FirstLastRelative firstLast, // get first and last relative
234 const std::string& what, // what are we looking for
235 const std::string& ind = ""); // logging indent
244 void selectFromAncestor(int ancestor, Pythia8::Event& inputEvent, Pythia8::Event& outputEvent);
247 void getNcoll(const Pythia8::Info& info, int& nColl);
250 void getNpart(const Pythia8::Info& info, int& nPart);
254 void getNpart(const Pythia8::Info& info, int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg);
256 void getNremn(const Pythia8::Event& event, int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg);
260 void getNfreeSpec(const Pythia8::Info& info, int& nFreenProj, int& nFreepProj, int& nFreenTarg, int& nFreepTarg);
264 void seedGenerator();
265
267 Pythia8::Pythia mPythia;
268
272 std::string mConfig;
274 std::string mHooksFileName;
277 std::string mHooksFuncName;
280 UserFilterFcn mUserFilterFcn = [](Pythia8::Particle const&) -> bool { return true; };
282
283 bool mApplyPruning = false;
284 bool mIsInitialized = false; // if Init function has been called
285 long mInitialRNGSeed = -1; // initial seed for Pythia random number state;
286 // will be transported to Pythia in the Init function through the Pythia::readString("Random:seed") mechanism.
287 // Value of -1 means unitialized; 0 will be time-dependent and values >1 <= MAX_SEED concrete reproducible seeding
288 Pythia8GenConfig mGenConfig; // configuration object
289
290 static std::atomic<int> Pythia8InstanceCounter;
292
293 constexpr static long MAX_SEED = 900000000;
294
296
297};
299/*****************************************************************/
300/*****************************************************************/
301
302} // namespace eventgen
303} // namespace o2
304
305#endif /* ALICEO2_EVENTGEN_GENERATORPYTHIA8_H_ */
306//
307// EOF
308//
bool done
~GeneratorPythia8() override=default
GeneratorPythia8 & operator=(const GeneratorPythia8 &)
void(* SetRelatives)(Pythia8::Particle &, int, int)
std::pair< int, int >(* FirstLastRelative)(const Pythia8::Particle &)
std::vector< int >(* GetRelatives)(const Pythia8::Particle &)
GeneratorPythia8(const GeneratorPythia8 &)
void setHooksFuncName(std::string val)
void selectFromAncestor(int ancestor, Pythia8::Event &inputEvent, Pythia8::Event &outputEvent)
void pruneEvent(Pythia8::Event &event, Select select)
bool readString(std::string val)
static std::atomic< int > Pythia8InstanceCounter
void setHooksFileName(std::string val)
void investigateRelatives(Pythia8::Event &event, const std::vector< int > &old2New, size_t index, std::vector< bool > &done, GetRelatives getter, SetRelatives setter, FirstLastRelative firstLast, const std::string &what, const std::string &ind="")
ClassDefOverride(GeneratorPythia8, 1)
void getNremn(int &nProtonProj, int &nNeutronProj, int &nProtonTarg, int &nNeutronTarg)
void seedGenerator()
performs seeding of the random state of Pythia (called from Init)
void getNpart(int &nProtonProj, int &nNeutronProj, int &nProtonTarg, int &nNeutronTarg)
std::function< bool(const Pythia8::Particle &)> UserFilterFcn
void setUserHooks(Pythia8::UserHooks *hooks)
void updateHeader(o2::dataformats::MCEventHeader *eventHeader) override
void getNfreeSpec(int &nFreenProj, int &nFreepProj, int &nFreenTarg, int &nFreepTarg)
void setConfig(std::string val)
struct _cl_event * event
Definition glcorearb.h:2982
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint GLfloat * val
Definition glcorearb.h:1582
std::vector< InputSpec > select(char const *matcher="")
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...