12 #ifndef Pythia8_HeavyIons_H
13 #define Pythia8_HeavyIons_H
15 #include "Pythia8/HIUserHooks.h"
52 virtual bool init() = 0;
58 virtual bool next() = 0;
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);
202 bool setupFullCollision(
EventInfo & ei,
const SubCollision & coll,
206 static int getBeam(
Event & ev,
int i);
214 bool colConnect =
false);
218 vector<int> findRecoilers(
const Event & e,
bool tside,
int beam,
int end,
219 const Vec4 & pdiff,
const Vec4 & pbeam);
223 static void addJunctions(
Event & evnt,
Event & sub,
int coloff);
228 const vector<Nucleon> & targ);
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));
244 struct ProcessSelectorHook:
public UserHooks {
246 ProcessSelectorHook(): proc(0), b(-1.0) {}
249 virtual bool canVetoProcessLevel() {
254 virtual bool doVetoProcessLevel(
Event&) {
255 return proc > 0 && infoPtr->code() != proc;
259 virtual bool canSetImpactParameter()
const {
264 virtual double doSetImpactParameter() {
280 HoldProcess(ProcessSelectorHook & hook,
int proc,
double b = -1.0)
281 : saveHook(&hook), saveProc(hook.proc), saveB(hook.b) {
289 saveHook->proc = saveProc;
295 ProcessSelectorHook * saveHook;
306 ProcessSelectorHook selectMB;
309 ProcessSelectorHook selectSASD;
313 static const int MAXTRY = 999;
314 static const int MAXEVSAVE = 999;
317 vector<Nucleon> projectile;
318 vector<Nucleon> target;
321 multiset<SubCollision> subColls;
328 ImpactParameterGenerator * bGenPtr;
332 NucleusModel * projPtr;
333 NucleusModel * targPtr;
337 SubCollisionModel * collPtr;
356 osSave->rdbuf(ros.rdbuf());
377 #endif // Pythia8_HeavyIons_H
static bool getTransforms(Vec4 p1, Vec4 p2, const Vec4 &p1p, pair< RotBstMatrix, RotBstMatrix > &R12, int, int)
HeavyIons(Pythia &mainPythiaIn)
ostream * osSave
The redirected ostream.
Optional object for signal processes (nn).
virtual bool next()
Produce a collision involving heavy ions.
virtual ~HeavyIons()
Destructor.
bool addNucleonExcitation(EventInfo &orig, EventInfo &add, bool colConnect=false)
std::streambuf * bufferSave
The original buffer of the redirected ostream.
bool setHIUserHooksPtr(HIUserHooks *userHooksPtrIn)
Possibility to pass in pointer for special heavy ion user hooks.
static void addSpecialSettings(Settings &settings)
EventInfo getSignal(const SubCollision &coll)
Generate events from the internal Pythia oblects;.
bool nextSASD(int proc)
Generate a single diffractive.
Optional object for signal processes (pp).
PythiaObject
Enumerate the different internal Pythia objects.
EventInfo mkEventInfo(Pythia &pyt, const SubCollision *coll=0)
Setup an EventInfo object from a Pythia instance.
vector< Pythia * > pythia
void clearProcessLevel(Pythia &pyt)
Single diffractive as one side of non-diffractive.
void sumUpMessages(Info &in, string tag, const Info &other)
static void setupSpecials(Settings &settings, string match)
Duplicate setting on the form match: to settings on the form HImatch:
Optional object for signal processes (pn).
virtual bool init()
Initialize Angantyr.
Optional object for signal processes (np).
Angantyr(Pythia &mainPythiaIn)
Redirect(ostream &os, ostream &ros)
Redirect os to ros for the lifetime of this object.
bool setUserHooksPtr(PythiaObject sel, UserHooks *userHooksPtrIn)
Set UserHooks for specific (or ALL) internal Pythia objects.
bool addNucleusRemnants(const vector< Nucleon > &proj, const vector< Nucleon > &targ)
Status
Enum for specifying the status of a nucleon.
vector< string > pythiaNames
The names associated with the secondary pythia objects.
void addSubEvent(Event &evnt, Event &sub)
Add a sub-event to the final event record.
internal class to redirect stdout
The default HeavyIon model in Pythia.
static bool isHeavyIon(Settings &settings)