StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HeavyIons.h
1 // HeavyIons.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This file contains the definition of the HeavyIons class which
7 // Provides Pythia with infrastructure to combine several nucleon
8 // collisions into a single heavy ion collision event. This file also
9 // includes the definition of the Angantyr class which implements the
10 // default model for heavy ion collisions in Pythia.
11 
12 #ifndef Pythia8_HeavyIons_H
13 #define Pythia8_HeavyIons_H
14 
15 #include "Pythia8/HIUserHooks.h"
16 
17 namespace Pythia8 {
18 
20 class Pythia;
21 
22 //==========================================================================
23 
29 
30 class HeavyIons {
31 
32 public:
33 
37  HeavyIons(Pythia & mainPythiaIn)
38  : mainPythiaPtr(&mainPythiaIn), HIHooksPtr(0),
39  pythia(vector<Pythia*>(1, &mainPythiaIn)) {
40  mainPythiaPtr->info.hiinfo = &hiinfo;
41  }
42 
44  virtual ~HeavyIons() {}
45 
52  virtual bool init() = 0;
53 
58  virtual bool next() = 0;
59 
62  static void addSpecialSettings(Settings & settings);
63 
66  static bool isHeavyIon(Settings & settings);
67 
69  bool setHIUserHooksPtr(HIUserHooks * userHooksPtrIn) {
70  HIHooksPtr = userHooksPtrIn; return true;
71  }
72 
73 protected:
74 
78  void sumUpMessages(Info & in, string tag, const Info & other);
79 
82  void updateInfo();
83 
89  void clearProcessLevel(Pythia & pyt);
90 
92  static void setupSpecials(Settings & settings, string match);
93 
96  static void setupSpecials(Pythia & p, string match);
97 
101 
105 
109 
112  vector<Pythia *> pythia;
113 
115  vector<string> pythiaNames;
116 
117 public:
118 
123 
124  // Print out statistics.
125  virtual void stat();
126 
127 };
128 
129 //==========================================================================
130 
132 
133 class Angantyr: public HeavyIons {
134 
135 public:
136 
139  HADRON = 0,
140  MBIAS = 1,
141  SASD = 2,
142  SIGPP = 3,
143  SIGPN = 4,
144  SIGNP = 5,
145  SIGNN = 6,
146  ALL = 7
147  };
148 
149 public:
150 
154  Angantyr(Pythia & mainPythiaIn);
155 
156  virtual ~Angantyr();
157 
159  virtual bool init();
160 
162  virtual bool next();
163 
164 
166  bool setUserHooksPtr(PythiaObject sel, UserHooks * userHooksPtrIn);
167 
168 protected:
169 
171  EventInfo mkEventInfo(Pythia & pyt, const SubCollision * coll = 0);
172 
174  EventInfo getSignal(const SubCollision & coll);
175  EventInfo getND() { return getMBIAS(0, 101); }
176  EventInfo getND(const SubCollision & coll) { return getMBIAS(&coll, 101); }
177  EventInfo getEl(const SubCollision & coll) { return getMBIAS(&coll, 102); }
178  EventInfo getSDP(const SubCollision & coll) { return getMBIAS(&coll, 103); }
179  EventInfo getSDT(const SubCollision & coll) { return getMBIAS(&coll, 104); }
180  EventInfo getDD(const SubCollision & coll) { return getMBIAS(&coll, 105); }
181  EventInfo getCD(const SubCollision & coll) { return getMBIAS(&coll, 106); }
182  EventInfo getSDabsP(const SubCollision & coll)
183  { return getSASD(&coll, 103); }
184  EventInfo getSDabsT(const SubCollision & coll)
185  { return getSASD(&coll, 104); }
186  EventInfo getMBIAS(const SubCollision * coll, int procid);
187  EventInfo getSASD(const SubCollision * coll, int procid);
188  bool genAbs(const multiset<SubCollision> & coll,
189  list<EventInfo> & subevents);
190  void addSASD(const multiset<SubCollision> & coll);
191  bool addDD(const multiset<SubCollision> & coll, list<EventInfo> & subevents);
192  bool addSD(const multiset<SubCollision> & coll, list<EventInfo> & subevents);
193  void addSDsecond(const multiset<SubCollision> & coll);
194  bool addCD(const multiset<SubCollision> & coll, list<EventInfo> & subevents);
195  void addCDsecond(const multiset<SubCollision> & coll);
196  bool addEL(const multiset<SubCollision> & coll, list<EventInfo> & subevents);
197  void addELsecond(const multiset<SubCollision> & coll);
198  bool buildEvent(list<EventInfo> & subevents,
199  const vector<Nucleon> & proj,
200  const vector<Nucleon> & targ);
201 
202  bool setupFullCollision(EventInfo & ei, const SubCollision & coll,
203  Nucleon::Status ptype, Nucleon::Status ttype);
204  bool fixIsoSpin(EventInfo & ei);
205  EventInfo & shiftEvent(EventInfo & ei);
206  static int getBeam(Event & ev, int i);
207 
209  bool nextSASD(int proc);
210 
213  bool addNucleonExcitation(EventInfo & orig, EventInfo & add,
214  bool colConnect = false);
215 
216  // Find the recoilers in the current event to conserve energy and
217  // momentum in addNucleonExcitation.
218  vector<int> findRecoilers(const Event & e, bool tside, int beam, int end,
219  const Vec4 & pdiff, const Vec4 & pbeam);
220 
222  void addSubEvent(Event & evnt, Event & sub);
223  static void addJunctions(Event & evnt, Event & sub, int coloff);
224 
227  bool addNucleusRemnants(const vector<Nucleon> & proj,
228  const vector<Nucleon> & targ);
229 
230 public:
231 
234  static bool
235  getTransforms(Vec4 p1, Vec4 p2, const Vec4 & p1p,
236  pair<RotBstMatrix,RotBstMatrix> & R12, int, int);
237  static double mT2(const Vec4 & p) { return p.pPos()*p.pNeg(); }
238  static double mT(const Vec4 & p) { return sqrt(max(mT2(p), 0.0));
239  }
240 
241 private:
242 
243  // Private UserHooks class to select a specific process.
244  struct ProcessSelectorHook: public UserHooks {
245 
246  ProcessSelectorHook(): proc(0), b(-1.0) {}
247 
248  // Yes we can veto event after process-level selection.
249  virtual bool canVetoProcessLevel() {
250  return true;
251  }
252 
253  // Veto any unwanted process.
254  virtual bool doVetoProcessLevel(Event&) {
255  return proc > 0 && infoPtr->code() != proc;
256  }
257 
258  // Can set the overall impact parameter for the MPI treatment.
259  virtual bool canSetImpactParameter() const {
260  return b >= 0.0;
261  }
262 
263  // Set the overall impact parameter for the MPI treatment.
264  virtual double doSetImpactParameter() {
265  return b;
266  }
267 
268  // The wanted process;
269  int proc;
270 
271  // The selected b-value.
272  double b;
273 
274  };
275 
276  // Holder class to temporarily select a specific process
277  struct HoldProcess {
278 
279  // Set the given process for the given hook object.
280  HoldProcess(ProcessSelectorHook & hook, int proc, double b = -1.0)
281  : saveHook(&hook), saveProc(hook.proc), saveB(hook.b) {
282  hook.proc = proc;
283  hook.b = b;
284  }
285 
286  // Reset the process of the hook object given in the constructor.
287  ~HoldProcess() {
288  if ( saveHook ) {
289  saveHook->proc = saveProc;
290  saveHook->b = saveB;
291  }
292  }
293 
294  // The hook object.
295  ProcessSelectorHook * saveHook;
296 
297  // The previous process of the hook object.
298  int saveProc;
299 
300  // The previous b-value of the hook object.
301  double saveB;
302 
303  };
304 
305  // The process selector for standard minimum bias processes.
306  ProcessSelectorHook selectMB;
307 
308  // The process selector for the SASD object.
309  ProcessSelectorHook selectSASD;
310 
311 private:
312 
313  static const int MAXTRY = 999;
314  static const int MAXEVSAVE = 999;
315 
317  vector<Nucleon> projectile;
318  vector<Nucleon> target;
319 
321  multiset<SubCollision> subColls;
322 
325  bool hasSignal;
326 
328  ImpactParameterGenerator * bGenPtr;
329 
332  NucleusModel * projPtr;
333  NucleusModel * targPtr;
334 
337  SubCollisionModel * collPtr;
338 
341  int recoilerMode;
342 
344  int bMode;
345 
346 public:
347 
349  class Redirect {
350 
351  public:
352 
354  Redirect(ostream & os, ostream & ros)
355  : osSave(&os), bufferSave(os.rdbuf()) {
356  osSave->rdbuf(ros.rdbuf());
357  }
358 
359  ~Redirect() {
360  osSave->rdbuf(bufferSave);
361  }
362 
364  ostream * osSave;
365 
367  std::streambuf * bufferSave;
368 
369  };
370 
371 };
372 
373 //==========================================================================
374 
375 } // end namespace Pythia8
376 
377 #endif // Pythia8_HeavyIons_H
static bool getTransforms(Vec4 p1, Vec4 p2, const Vec4 &p1p, pair< RotBstMatrix, RotBstMatrix > &R12, int, int)
Definition: HeavyIons.cc:1293
HeavyIons(Pythia &mainPythiaIn)
Definition: HeavyIons.h:37
HIUserHooks * HIHooksPtr
Definition: HeavyIons.h:108
SigmaTotal sigtot
Definition: HeavyIons.h:104
ostream * osSave
The redirected ostream.
Definition: HeavyIons.h:364
Optional object for signal processes (nn).
Definition: HeavyIons.h:145
virtual bool next()
Produce a collision involving heavy ions.
Definition: HeavyIons.cc:1582
virtual ~HeavyIons()
Destructor.
Definition: HeavyIons.h:44
bool addNucleonExcitation(EventInfo &orig, EventInfo &add, bool colConnect=false)
Definition: HeavyIons.cc:1158
std::streambuf * bufferSave
The original buffer of the redirected ostream.
Definition: HeavyIons.h:367
bool setHIUserHooksPtr(HIUserHooks *userHooksPtrIn)
Possibility to pass in pointer for special heavy ion user hooks.
Definition: HeavyIons.h:69
static void addSpecialSettings(Settings &settings)
Definition: HeavyIons.cc:24
EventInfo getSignal(const SubCollision &coll)
Generate events from the internal Pythia oblects;.
Definition: HeavyIons.cc:594
bool nextSASD(int proc)
Generate a single diffractive.
Definition: HeavyIons.cc:1426
Optional object for signal processes (pp).
Definition: HeavyIons.h:142
PythiaObject
Enumerate the different internal Pythia objects.
Definition: HeavyIons.h:138
Minimum Bias processed.
Definition: HeavyIons.h:140
EventInfo mkEventInfo(Pythia &pyt, const SubCollision *coll=0)
Setup an EventInfo object from a Pythia instance.
Definition: HeavyIons.cc:319
Indicates all objects.
Definition: HeavyIons.h:146
For hadronization only.
Definition: HeavyIons.h:139
vector< Pythia * > pythia
Definition: HeavyIons.h:112
void clearProcessLevel(Pythia &pyt)
Definition: HeavyIons.cc:126
Single diffractive as one side of non-diffractive.
Definition: HeavyIons.h:141
void sumUpMessages(Info &in, string tag, const Info &other)
Definition: HeavyIons.cc:115
static void setupSpecials(Settings &settings, string match)
Duplicate setting on the form match: to settings on the form HImatch:
Definition: HeavyIons.cc:32
Optional object for signal processes (pn).
Definition: HeavyIons.h:143
virtual bool init()
Initialize Angantyr.
Definition: HeavyIons.cc:340
Pythia * mainPythiaPtr
Definition: HeavyIons.h:100
Optional object for signal processes (np).
Definition: HeavyIons.h:144
Angantyr(Pythia &mainPythiaIn)
Definition: HeavyIons.cc:272
Redirect(ostream &os, ostream &ros)
Redirect os to ros for the lifetime of this object.
Definition: HeavyIons.h:354
Definition: beam.h:43
Definition: AgUStep.h:26
bool setUserHooksPtr(PythiaObject sel, UserHooks *userHooksPtrIn)
Set UserHooks for specific (or ALL) internal Pythia objects.
Definition: HeavyIons.cc:307
bool addNucleusRemnants(const vector< Nucleon > &proj, const vector< Nucleon > &targ)
Definition: HeavyIons.cc:1497
virtual bool init()=0
Status
Enum for specifying the status of a nucleon.
Definition: HIUserHooks.h:69
virtual bool next()=0
vector< string > pythiaNames
The names associated with the secondary pythia objects.
Definition: HeavyIons.h:115
void addSubEvent(Event &evnt, Event &sub)
Add a sub-event to the final event record.
Definition: HeavyIons.cc:1372
internal class to redirect stdout
Definition: HeavyIons.h:349
The default HeavyIon model in Pythia.
Definition: HeavyIons.h:133
static bool isHeavyIon(Settings &settings)
Definition: HeavyIons.cc:258