ATLAS Offline Software
Loading...
Searching...
No Matches
AddFlowByShifting.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
5// File: Generators/FlowAfterburner/AddFlowByShifting.h
6// Description:
7// This code is used to introduce particle flow
8// to particles from generated events
9// It works by modifying phi angles of particles
10// according to requested flow type and magnitude
11//
12// It takes from SG a container of tracks as input
13// and registers in SG a new container with modified tracks on output
14//
15// It currently uses Hijing generator specific class HijingEventParams
16// with truth event parameters information
17//
18// AuthorList:
19// Andrzej Olszewski: Initial Code February 2006
20
21#ifndef ADDFLOWBYSHIFTING_H
22#define ADDFLOWBYSHIFTING_H
23
25#include "GaudiKernel/ServiceHandle.h"
26//
28#include "AtlasHepMC/GenParticle.h" //typedef GenParticlePtr
29#include "AtlasHepMC/GenVertex.h" //GenParticlePtr
30
31
32
34namespace CLHEP{
35 class HepRandomEngine;
36}
37class TGraph;
38
40public:
41 AddFlowByShifting(const std::string& name, ISvcLocator* pSvcLocator);
42 StatusCode initialize();
43 StatusCode execute();
44 StatusCode finalize();
45
46 //functions used for root finding when using the "exact"(and not the "approximate") method
47 static double vn_func (double x, void *params);
48
49
50private:
51 CLHEP::HepRandomEngine* getRandomEngine(const std::string& streamName, const EventContext& ctx) const;
52
53 double SetParentToRanPhi(HepMC::GenParticlePtr parent, CLHEP::HepRandomEngine *rndmEngine);
55 const HijingEventParams *hijing_pars);
56 void MoveDescendantsToParent(HepMC::GenParticlePtr parent, double phishift);
57
58
59 // flow functions to set the vn values
60 void (AddFlowByShifting::*m_flow_function) (double b, double eta, double pt);//function pointer which is set to one of the functions below
61 void jjia_minbias_new (double b, double eta, double pt);
62 void jjia_minbias_new_v2only(double b, double eta, double pt);
63 void fixed_vn (double b, double eta, double pt);
64 void fixed_v2 (double b, double eta, double pt);
65 void jjia_minbias_old (double b, double eta, double pt);
66 void ao_test (double b, double eta, double pt);
67 void custom_vn (double b, double eta, double pt);
68 void p_Pb_cent_eta_indep (double b, double eta, double pt); //for p_Pb
69 void OO_eta_indep (double b, double eta, double pt); //for OO
70
71 TGraph *m_graph_fluc{};//TGraph storing the v2_RP/delta Vs b_imp
72 void Set_EbE_Fluctuation_Multipliers(HepMC::GenVertexPtr mainvtx, float b, CLHEP::HepRandomEngine *rndmEngine);
73
74 // Random number service
75 ServiceHandle<IAthRNGSvc> m_rndmSvc{this, "RndmSvc", "AthRNGSvc"};
76
77 // Setable Properties:-
78 StringProperty m_inkey {this, "McTruthKey" , "GEN_EVENT" }; //FIXME use Read/WriteHandles
79 StringProperty m_outkey {this, "McFlowKey" , "FLOW_EVENT" }; //FIXME use Read/WriteHandles
80 IntegerProperty m_ranphi_sw {this, "RandomizePhi" , 0 };
81 StringProperty m_flow_function_name {this, "FlowFunctionName" , "jjia_minbias_new"};
82 StringProperty m_flow_implementation{this, "FlowImplementation", "exact" };
84 BooleanProperty m_flow_fluctuations {this, "FlowFluctuations" , false };
85 IntegerProperty m_floweta_sw {this, "FlowEtaSwitch" , 0 };
86 FloatProperty m_flow_maxeta {this, "FlowMaxEtaCut" , 10.0 };
87 FloatProperty m_flow_mineta {this, "FlowMinEtaCut" , 0.f };
88 IntegerProperty m_flowpt_sw {this, "FlowPtSwitch" , 0 };
89 FloatProperty m_flow_maxpt {this, "FlowMaxPtCut" , 1000000.f };
90 FloatProperty m_flow_minpt {this, "FlowMinPtCut" , 0.f };
91 IntegerProperty m_flowb_sw {this, "FlowBSwitch" , 0 };//currently not used
92 FloatProperty m_custom_v1 {this, "custom_v1" , 0.f };
93 FloatProperty m_custom_v2 {this, "custom_v2" , 0.f };
94 FloatProperty m_custom_v3 {this, "custom_v3" , 0.f };
95 FloatProperty m_custom_v4 {this, "custom_v4" , 0.f };
96 FloatProperty m_custom_v5 {this, "custom_v5" , 0.f };
97 FloatProperty m_custom_v6 {this, "custom_v6" , 0.f };
98
100
103 v1=0,
104 v2=1,
105 v3=2,
106 v4=3,
107 v5=4,
108 v6=5,
109 };
113};
114
115#endif
Scalar eta() const
pseudorapidity method
#define x
StringProperty m_flow_function_name
StringProperty m_outkey
IntegerProperty m_flowpt_sw
AddFlowByShifting(const std::string &name, ISvcLocator *pSvcLocator)
FloatProperty m_flow_mineta
void Set_EbE_Fluctuation_Multipliers(HepMC::GenVertexPtr mainvtx, float b, CLHEP::HepRandomEngine *rndmEngine)
static double vn_func(double x, void *params)
double AddFlowToParent(HepMC::GenParticlePtr parent, const HijingEventParams *hijing_pars)
FloatProperty m_custom_v3
void jjia_minbias_new(double b, double eta, double pt)
void ao_test(double b, double eta, double pt)
double SetParentToRanPhi(HepMC::GenParticlePtr parent, CLHEP::HepRandomEngine *rndmEngine)
ServiceHandle< IAthRNGSvc > m_rndmSvc
FloatProperty m_custom_v6
StringProperty m_flow_implementation
FloatProperty m_custom_v1
FloatProperty m_custom_v2
IntegerProperty m_floweta_sw
void p_Pb_cent_eta_indep(double b, double eta, double pt)
FloatProperty m_flow_minpt
FloatProperty m_flow_maxeta
IntegerProperty m_flowb_sw
void jjia_minbias_new_v2only(double b, double eta, double pt)
void MoveDescendantsToParent(HepMC::GenParticlePtr parent, double phishift)
CLHEP::HepRandomEngine * getRandomEngine(const std::string &streamName, const EventContext &ctx) const
void jjia_minbias_old(double b, double eta, double pt)
void fixed_v2(double b, double eta, double pt)
StringProperty m_inkey
IntegerProperty m_ranphi_sw
FloatProperty m_flow_maxpt
BooleanProperty m_flow_fluctuations
void OO_eta_indep(double b, double eta, double pt)
void custom_vn(double b, double eta, double pt)
float m_psi_n[Harmonic::NumHar]
void fixed_vn(double b, double eta, double pt)
void(AddFlowByShifting::* m_flow_function)(double b, double eta, double pt)
float m_EbE_Multiplier_vn[Harmonic::NumHar]
FloatProperty m_custom_v4
FloatProperty m_custom_v5
float m_v_n[Harmonic::NumHar]
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
GenParticle * GenParticlePtr
Definition GenParticle.h:37