Line data Source code
1 : // UserHooks.h is a part of the PYTHIA event generator.
2 : // Copyright (C) 2015 Torbjorn Sjostrand.
3 : // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 : // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 :
6 : // Header file to allow user access to program at different stages.
7 : // UserHooks: almost empty base class, with user to write the rela code.
8 : // MyUserHooks: derived class, only intended as an example.
9 :
10 : #ifndef Pythia8_UserHooks_H
11 : #define Pythia8_UserHooks_H
12 :
13 : #include "Pythia8/Event.h"
14 : #include "Pythia8/PartonSystems.h"
15 : #include "Pythia8/PythiaStdlib.h"
16 : #include "Pythia8/SigmaProcess.h"
17 :
18 : namespace Pythia8 {
19 :
20 : //==========================================================================
21 :
22 : // Forward reference to the PhaseSpace class.
23 : class PhaseSpace;
24 :
25 : //==========================================================================
26 :
27 : // UserHooks is base class for user access to program execution.
28 :
29 : class UserHooks {
30 :
31 : public:
32 :
33 : // Destructor.
34 0 : virtual ~UserHooks() {}
35 :
36 : // Initialize pointers and workEvent. Note: not virtual.
37 : void initPtr( Info* infoPtrIn, Settings* settingsPtrIn,
38 : ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
39 : BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
40 : BeamParticle* beamPomAPtrIn, BeamParticle* beamPomBPtrIn,
41 : CoupSM* coupSMPtrIn, PartonSystems* partonSystemsPtrIn,
42 0 : SigmaTotal* sigmaTotPtrIn) { infoPtr = infoPtrIn;
43 0 : settingsPtr = settingsPtrIn; particleDataPtr = particleDataPtrIn;
44 0 : rndmPtr = rndmPtrIn; beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;
45 0 : beamPomAPtr = beamPomAPtrIn; beamPomBPtr = beamPomBPtrIn;
46 0 : coupSMPtr = coupSMPtrIn; partonSystemsPtr = partonSystemsPtrIn;
47 0 : sigmaTotPtr = sigmaTotPtrIn;
48 0 : workEvent.init("(work event)", particleDataPtr);}
49 :
50 : // Initialisation after beams have been set by Pythia::init().
51 0 : virtual bool initAfterBeams() { return true; }
52 :
53 : // Possibility to modify cross section of process.
54 0 : virtual bool canModifySigma() {return false;}
55 :
56 : // Multiplicative factor modifying the cross section of a hard process.
57 : virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
58 : const PhaseSpace* phaseSpacePtr, bool inEvent);
59 :
60 : // Possibility to bias selection of events, compensated by a weight.
61 0 : virtual bool canBiasSelection() {return false;}
62 :
63 : // Multiplicative factor in the phase space selection of a hard process.
64 : virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
65 : const PhaseSpace* phaseSpacePtr, bool inEvent);
66 :
67 : // Event weight to compensate for selection weight above.
68 0 : virtual double biasedSelectionWeight() {return 1./selBias;}
69 :
70 : // Possibility to veto event after process-level selection.
71 0 : virtual bool canVetoProcessLevel() {return false;}
72 :
73 : // Decide whether to veto current process or not, based on process record.
74 : // Usage: doVetoProcessLevel( process).
75 0 : virtual bool doVetoProcessLevel(Event& ) {return false;}
76 :
77 : // Possibility to veto resonance decay chain.
78 0 : virtual bool canVetoResonanceDecays() {return false;}
79 :
80 : // Decide whether to veto current resonance decay chain or not, based on
81 : // process record. Usage: doVetoProcessLevel( process).
82 0 : virtual bool doVetoResonanceDecays(Event& ) {return false;}
83 :
84 : // Possibility to veto MPI + ISR + FSR evolution and kill event,
85 : // making decision at a fixed pT scale. Useful for MLM-style matching.
86 0 : virtual bool canVetoPT() {return false;}
87 :
88 : // Transverse-momentum scale for veto test.
89 0 : virtual double scaleVetoPT() {return 0.;}
90 :
91 : // Decide whether to veto current event or not, based on event record.
92 : // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
93 : // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
94 : // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
95 0 : virtual bool doVetoPT( int , const Event& ) {return false;}
96 :
97 : // Possibility to veto MPI + ISR + FSR evolution and kill event,
98 : // making decision after fixed number of ISR or FSR steps.
99 0 : virtual bool canVetoStep() {return false;}
100 :
101 : // Up to how many ISR + FSR steps of hardest interaction should be checked.
102 0 : virtual int numberVetoStep() {return 1;}
103 :
104 : // Decide whether to veto current event or not, based on event record.
105 : // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
106 : // nISR and nFSR number of emissions so far for hard interaction only.
107 0 : virtual bool doVetoStep( int , int , int , const Event& ) {return false;}
108 :
109 : // Possibility to veto MPI + ISR + FSR evolution and kill event,
110 : // making decision after fixed number of MPI steps.
111 0 : virtual bool canVetoMPIStep() {return false;}
112 :
113 : // Up to how many MPI steps should be checked.
114 0 : virtual int numberVetoMPIStep() {return 1;}
115 :
116 : // Decide whether to veto current event or not, based on event record.
117 : // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
118 0 : virtual bool doVetoMPIStep( int , const Event& ) {return false;}
119 :
120 : // Possibility to veto event after ISR + FSR + MPI in parton level,
121 : // but before beam remnants and resonance decays.
122 0 : virtual bool canVetoPartonLevelEarly() {return false;}
123 :
124 : // Decide whether to veto current partons or not, based on event record.
125 : // Usage: doVetoPartonLevelEarly( event).
126 0 : virtual bool doVetoPartonLevelEarly( const Event& ) {return false;}
127 :
128 : // Retry same ProcessLevel with a new PartonLevel after a veto in
129 : // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
130 : // if you overload this method to return true.
131 0 : virtual bool retryPartonLevel() {return false;}
132 :
133 : // Possibility to veto event after parton-level selection.
134 0 : virtual bool canVetoPartonLevel() {return false;}
135 :
136 : // Decide whether to veto current partons or not, based on event record.
137 : // Usage: doVetoPartonLevel( event).
138 0 : virtual bool doVetoPartonLevel( const Event& ) {return false;}
139 :
140 : // Possibility to set initial scale in TimeShower for resonance decay.
141 0 : virtual bool canSetResonanceScale() {return false;}
142 :
143 : // Initial scale for TimeShower evolution.
144 : // Usage: scaleResonance( iRes, event), where iRes is location
145 : // of decaying resonance in the event record.
146 0 : virtual double scaleResonance( int, const Event& ) {return 0.;}
147 :
148 : // Possibility to veto an emission in the ISR machinery.
149 0 : virtual bool canVetoISREmission() {return false;}
150 :
151 : // Decide whether to veto current emission or not, based on event record.
152 : // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
153 : // of event record before current emission-to-be-scrutinized was added,
154 : // and iSys is the system of the radiation (according to PartonSystems).
155 0 : virtual bool doVetoISREmission( int, const Event&, int ) {return false;}
156 :
157 : // Possibility to veto an emission in the FSR machinery.
158 0 : virtual bool canVetoFSREmission() {return false;}
159 :
160 : // Decide whether to veto current emission or not, based on event record.
161 : // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
162 : // sizeOld is size of event record before current emission-to-be-scrutinized
163 : // was added, iSys is the system of the radiation (according to
164 : // PartonSystems), and inResonance is true if the emission takes place in a
165 : // resonance decay.
166 : virtual bool doVetoFSREmission( int, const Event&, int, bool = false )
167 0 : {return false;}
168 :
169 : // Possibility to veto an MPI.
170 0 : virtual bool canVetoMPIEmission() { return false; }
171 :
172 : // Decide whether to veto an MPI based on event record.
173 : // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
174 : // is size of event record before the current MPI.
175 0 : virtual bool doVetoMPIEmission( int, const Event &) { return false; }
176 :
177 : // Possibility to reconnect colours from resonance decay systems.
178 0 : virtual bool canReconnectResonanceSystems() { return false; }
179 :
180 : // Do reconnect colours from resonance decay systems.
181 : // Usage: doVetoFSREmission( oldSizeEvt, event)
182 : // where oldSizeEvent is the event size before resonance decays.
183 : // Should normally return true, while false means serious failure.
184 : // Value of PartonLevel:earlyResDec determines where method is called.
185 0 : virtual bool doReconnectResonanceSystems( int, Event &) {return true;}
186 :
187 : protected:
188 :
189 : // Constructor.
190 : UserHooks() : infoPtr(0), settingsPtr(0), particleDataPtr(0), rndmPtr(0),
191 : beamAPtr(0), beamBPtr(0), beamPomAPtr(0), beamPomBPtr(0), coupSMPtr(0),
192 : partonSystemsPtr(0), sigmaTotPtr(0), selBias(1.) {}
193 :
194 : // Pointer to various information on the generation.
195 : Info* infoPtr;
196 :
197 : // Pointer to the settings database.
198 : Settings* settingsPtr;
199 :
200 : // Pointer to the particle data table.
201 : ParticleData* particleDataPtr;
202 :
203 : // Pointer to the random number generator.
204 : Rndm* rndmPtr;
205 :
206 : // Pointers to the two incoming beams and to Pomeron beam-inside-beam.
207 : BeamParticle* beamAPtr;
208 : BeamParticle* beamBPtr;
209 : BeamParticle* beamPomAPtr;
210 : BeamParticle* beamPomBPtr;
211 :
212 : // Pointers to Standard Model couplings.
213 : CoupSM* coupSMPtr;
214 :
215 : // Pointer to information on subcollision parton locations.
216 : PartonSystems* partonSystemsPtr;
217 :
218 : // Pointer to the total/elastic/diffractive cross sections.
219 : SigmaTotal* sigmaTotPtr;
220 :
221 : // omitResonanceDecays omits resonance decay chains from process record.
222 : void omitResonanceDecays(const Event& process, bool finalOnly = false);
223 :
224 : // subEvent extracts currently resolved partons in the hard process.
225 : void subEvent(const Event& event, bool isHardest = true);
226 :
227 : // Have one event object around as work area.
228 : Event workEvent;
229 :
230 : // User-imposed selection bias.
231 : double selBias;
232 :
233 : };
234 :
235 : //==========================================================================
236 :
237 : // SuppressSmallPT is a derived class for user access to program execution.
238 : // It is a simple example, illustrating how to suppress the cross section
239 : // of 2 -> 2 processes by a factor pT^4 / (pT0^2 + pT^2)^2, with pT0 input,
240 : // and also modify alpha_strong scale similarly.
241 :
242 0 : class SuppressSmallPT : public UserHooks {
243 :
244 : public:
245 :
246 : // Constructor.
247 : SuppressSmallPT( double pT0timesMPIIn = 1., int numberAlphaSIn = 0,
248 : bool useSameAlphaSasMPIIn = true) : pT20(0.) {isInit = false;
249 : pT0timesMPI = pT0timesMPIIn; numberAlphaS = numberAlphaSIn;
250 : useSameAlphaSasMPI = useSameAlphaSasMPIIn;}
251 :
252 : // Possibility to modify cross section of process.
253 0 : virtual bool canModifySigma() {return true;}
254 :
255 : // Multiplicative factor modifying the cross section of a hard process.
256 : // Usage: inEvent is true for event generation, false for initialization.
257 : virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
258 : const PhaseSpace* phaseSpacePtr, bool );
259 :
260 : private:
261 :
262 : // Save input properties and the squared pT0 scale.
263 : bool isInit, useSameAlphaSasMPI;
264 : int numberAlphaS;
265 : double pT0timesMPI, pT20;
266 :
267 : // Alpha_strong calculation.
268 : AlphaStrong alphaS;
269 :
270 : };
271 :
272 : //==========================================================================
273 :
274 : } // end namespace Pythia8
275 :
276 : #endif // Pythia8_UserHooks_H
|