Line data Source code
1 : // Info.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 : // This file contains a class that keep track of generic event info.
7 : // Info: contains information on the generation process and errors.
8 :
9 : #ifndef Pythia8_Info_H
10 : #define Pythia8_Info_H
11 :
12 : #include "Pythia8/PythiaStdlib.h"
13 : #include "Pythia8/LHEF3.h"
14 :
15 : namespace Pythia8 {
16 :
17 : //==========================================================================
18 :
19 : // The Info class contains a mixed bag of information on the event
20 : // generation activity, especially on the current subprocess properties,
21 : // and on the number of errors encountered. This is used by the
22 : // generation machinery, but can also be read by the user.
23 : // Note: some methods that maybe should not be accessible to the user
24 : // are still public, to work also for user-written FSR/ISR classes.
25 :
26 0 : class Info {
27 :
28 : public:
29 :
30 : // Constructor.
31 0 : Info() : LHEFversionSave(0), initrwgt(NULL), generators(NULL),
32 0 : weightgroups(NULL), init_weights(NULL), eventAttributes(NULL),
33 0 : weights_detailed(NULL), weights_compressed(NULL), scales(NULL),
34 0 : weights(NULL), rwgt(NULL), eCMSave(0.), lowPTmin(false), a0MPISave(0.),
35 0 : abortPartonLevel(false), weightCKKWLSave(1.), weightFIRSTSave(0.) {
36 0 : for (int i = 0; i < 40; ++i) counters[i] = 0;}
37 :
38 : // Listing of most available information on current event.
39 : void list(ostream& os = cout) const;
40 :
41 : // Beam particles (in rest frame). CM energy of event.
42 0 : int idA() const {return idASave;}
43 0 : int idB() const {return idBSave;}
44 0 : double pzA() const {return pzASave;}
45 0 : double pzB() const {return pzBSave;}
46 0 : double eA() const {return eASave;}
47 0 : double eB() const {return eBSave;}
48 0 : double mA() const {return mASave;}
49 0 : double mB() const {return mBSave;}
50 0 : double eCM() const {return eCMSave;}
51 : double s() const {return sSave;}
52 :
53 : // Warnings from initialization.
54 : bool tooLowPTmin() const {return lowPTmin;}
55 :
56 : // Process name and code, and the number of final-state particles.
57 : string name() const {return nameSave;}
58 0 : int code() const {return codeSave;}
59 0 : int nFinal() const {return nFinalSave;}
60 :
61 : // Are beam particles resolved, with pdf's? Are they diffractive?
62 0 : bool isResolved() const {return isRes;}
63 0 : bool isDiffractiveA() const {return isDiffA;}
64 0 : bool isDiffractiveB() const {return isDiffB;}
65 0 : bool isDiffractiveC() const {return isDiffC;}
66 0 : bool isNonDiffractive() const {return isND;}
67 : // Retained for backwards compatibility.
68 : bool isMinBias() const {return isND;}
69 :
70 : // Information for Les Houches Accord and reading files.
71 : bool isLHA() const {return isLH;}
72 0 : bool atEndOfFile() const {return atEOF;}
73 :
74 : // For nondiffractive and Les Houches Accord identify hardest subprocess.
75 : bool hasSub(int i = 0) const {return hasSubSave[i];}
76 : string nameSub(int i = 0) const {return nameSubSave[i];}
77 : int codeSub(int i = 0) const {return codeSubSave[i];}
78 : int nFinalSub(int i = 0) const {return nFinalSubSave[i];}
79 :
80 : // Incoming parton flavours and x values.
81 0 : int id1(int i = 0) const {return id1Save[i];}
82 0 : int id2(int i = 0) const {return id2Save[i];}
83 0 : double x1(int i = 0) const {return x1Save[i];}
84 0 : double x2(int i = 0) const {return x2Save[i];}
85 : double y(int i = 0) const {return 0.5 * log( x1Save[i] / x2Save[i]);}
86 : double tau(int i = 0) const {return x1Save[i] * x2Save[i];}
87 :
88 : // Hard process flavours, x values, parton densities, couplings, Q2 scales.
89 0 : int id1pdf(int i = 0) const {return id1pdfSave[i];}
90 0 : int id2pdf(int i = 0) const {return id2pdfSave[i];}
91 0 : double x1pdf(int i = 0) const {return x1pdfSave[i];}
92 0 : double x2pdf(int i = 0) const {return x2pdfSave[i];}
93 0 : double pdf1(int i = 0) const {return pdf1Save[i];}
94 0 : double pdf2(int i = 0) const {return pdf2Save[i];}
95 0 : double QFac(int i = 0) const {return sqrtpos(Q2FacSave[i]);}
96 0 : double Q2Fac(int i = 0) const {return Q2FacSave[i];}
97 : bool isValence1() const {return isVal1;}
98 : bool isValence2() const {return isVal2;}
99 0 : double alphaS(int i = 0) const {return alphaSSave[i];}
100 0 : double alphaEM(int i = 0) const {return alphaEMSave[i];}
101 0 : double QRen(int i = 0) const {return sqrtpos(Q2RenSave[i]);}
102 0 : double Q2Ren(int i = 0) const {return Q2RenSave[i];}
103 0 : double scalup(int i = 0) const {return scalupSave[i];}
104 :
105 : // Mandelstam variables (notation as if subcollision).
106 : double mHat(int i = 0) const {return sqrt(sH[i]);}
107 : double sHat(int i = 0) const {return sH[i];}
108 : double tHat(int i = 0) const {return tH[i];}
109 : double uHat(int i = 0) const {return uH[i];}
110 : double pTHat(int i = 0) const {return pTH[i];}
111 : double pT2Hat(int i = 0) const {return pTH[i] * pTH[i];}
112 : double m3Hat(int i = 0) const {return m3H[i];}
113 : double m4Hat(int i = 0) const {return m4H[i];}
114 : double thetaHat(int i = 0) const {return thetaH[i];}
115 : double phiHat(int i = 0) const {return phiH[i];}
116 :
117 : // Weight of current event; normally 1, but used for Les Houches events
118 : // or when reweighting phase space selection. Conversion from mb to pb
119 : // for LHA strategy +-4. Also cumulative sum.
120 : double weight() const;
121 : double weightSum() const;
122 : double lhaStrategy() const {return lhaStrategySave;}
123 :
124 : // Number of times other steps have been carried out.
125 0 : int nISR() const {return nISRSave;}
126 0 : int nFSRinProc() const {return nFSRinProcSave;}
127 : int nFSRinRes() const {return nFSRinResSave;}
128 :
129 : // Maximum pT scales for MPI, ISR and FSR (in hard process).
130 : double pTmaxMPI() const {return pTmaxMPISave;}
131 : double pTmaxISR() const {return pTmaxISRSave;}
132 : double pTmaxFSR() const {return pTmaxFSRSave;}
133 :
134 : // Current evolution scale (for UserHooks).
135 0 : double pTnow() const {return pTnowSave;}
136 :
137 : // Impact parameter picture, global information
138 : double a0MPI() const {return a0MPISave;}
139 :
140 : // Impact parameter picture, as set by hardest interaction.
141 : double bMPI() const {return (bIsSet) ? bMPISave : 1.;}
142 0 : double enhanceMPI() const {return (bIsSet) ? enhanceMPISave : 1.;}
143 : double eMPI(int i) const {return (bIsSet) ? eMPISave[i] : 1.;}
144 :
145 : // Number of multiparton interactions, with code and pT for them.
146 0 : int nMPI() const {return nMPISave;}
147 0 : int codeMPI(int i) const {return codeMPISave[i];}
148 : double pTMPI(int i) const {return pTMPISave[i];}
149 : int iAMPI(int i) const {return iAMPISave[i];}
150 : int iBMPI(int i) const {return iBMPISave[i];}
151 :
152 : // Cross section estimate, optionally process by process.
153 : vector<int> codesHard();
154 : string nameProc(int i = 0) {return (i == 0) ? "sum"
155 : : ( (procNameM[i] == "") ? "unknown process" : procNameM[i] );}
156 : long nTried(int i = 0) {return (i == 0) ? nTry : nTryM[i];}
157 : long nSelected(int i = 0) {return (i == 0) ? nSel : nSelM[i];}
158 : long nAccepted(int i = 0) {return (i == 0) ? nAcc : nAccM[i];}
159 0 : double sigmaGen(int i = 0) {return (i == 0) ? sigGen : sigGenM[i];}
160 0 : double sigmaErr(int i = 0) {return (i == 0) ? sigErr : sigErrM[i];}
161 :
162 : // Counters for number of loops in various places.
163 0 : int getCounter( int i) const {return counters[i];}
164 :
165 : // Set or increase the value stored in a counter.
166 0 : void setCounter( int i, int value = 0) {counters[i] = value;}
167 0 : void addCounter( int i, int value = 1) {counters[i] += value;}
168 :
169 : // Reset to empty map of error messages.
170 0 : void errorReset() {messages.clear();}
171 :
172 : // Print a message the first few times. Insert in database.
173 : void errorMsg(string messageIn, string extraIn = " ",
174 : bool showAlways = false, ostream& os = cout);
175 :
176 : // Provide total number of errors/aborts/warnings experienced to date.
177 : int errorTotalNumber();
178 :
179 : // Print statistics on errors/aborts/warnings.
180 : void errorStatistics(ostream& os = cout);
181 :
182 : // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
183 0 : void setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
184 :
185 : // Set info on valence character of hard collision partons.
186 0 : void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
187 0 : isVal2 = isVal2In;}
188 :
189 : // Set and get some MPI/ISR/FSR properties needed for matching,
190 : // i.e. mainly of internal relevance.
191 0 : void hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
192 0 : bool hasHistory() {return hasHistorySave;}
193 0 : void zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
194 0 : double zNowISR() {return zNowISRSave;}
195 0 : void pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
196 0 : double pT2NowISR() {return pT2NowISRSave;}
197 :
198 : // Return CKKW-L weight.
199 0 : double getWeightCKKWL() const { return weightCKKWLSave;}
200 : // Set CKKW-L weight.
201 0 : void setWeightCKKWL( double weightIn) { weightCKKWLSave = weightIn;}
202 : // Return merging weight.
203 : double mergingWeight() const { return weightCKKWLSave;}
204 :
205 : // Return the complete NLO weight.
206 : double mergingWeightNLO() const {
207 : return (weightCKKWLSave - weightFIRSTSave); }
208 : // Return the O(\alpha_s)-term of the CKKW-L weight.
209 : double getWeightFIRST() const { return weightFIRSTSave;}
210 : // Set the O(\alpha_s)-term of the CKKW-L weight.
211 0 : void setWeightFIRST( double weightIn) { weightFIRSTSave = weightIn;}
212 :
213 : // Return an LHEF header
214 : string header(const string &key) {
215 0 : if (headers.find(key) == headers.end()) return "";
216 0 : else return headers[key];
217 0 : }
218 :
219 : // Return a list of all header key names
220 : vector < string > headerKeys();
221 :
222 : // Return the number of processes in the LHEF.
223 : int nProcessesLHEF() { return int(sigmaLHEFSave.size());}
224 : // Return the cross section information read from LHEF.
225 : double sigmaLHEF(int iProcess) { return sigmaLHEFSave[iProcess];}
226 :
227 : // LHEF3 information: Public for easy access.
228 : int LHEFversionSave;
229 :
230 : // Save process information as read from init block of LHEF.
231 : vector<double> sigmaLHEFSave;
232 :
233 : // Contents of the LHAinitrwgt tag
234 : LHAinitrwgt *initrwgt;
235 :
236 : // Contents of the LHAgenerator tags.
237 : vector<LHAgenerator > *generators;
238 :
239 : // A map of the LHAweightgroup tags, indexed by name.
240 : map<string,LHAweightgroup > *weightgroups;
241 :
242 : // A map of the LHAweight tags, indexed by name.
243 : map<string,LHAweight > *init_weights;
244 :
245 : // Store current-event Les Houches event tags.
246 : map<string, string > *eventAttributes;
247 :
248 : // The weights associated with this event, as given by the LHAwgt tags
249 : map<string,double > *weights_detailed;
250 :
251 : // The weights associated with this event, as given by the LHAweights tags
252 : vector<double > *weights_compressed;
253 :
254 : // Contents of the LHAscales tag.
255 : LHAscales *scales;
256 :
257 : // Contents of the LHAweights tag (compressed format)
258 : LHAweights *weights;
259 :
260 : // Contents of the LHArwgt tag (detailed format)
261 : LHArwgt *rwgt;
262 :
263 : // Set the LHEF3 objects read from the init and header blocks.
264 : void setLHEF3InitInfo();
265 : void setLHEF3InitInfo( int LHEFversionIn, LHAinitrwgt *initrwgtIn,
266 : vector<LHAgenerator> *generatorsIn,
267 : map<string,LHAweightgroup> *weightgroupsIn,
268 : map<string,LHAweight> *init_weightsIn );
269 : // Set the LHEF3 objects read from the event block.
270 : void setLHEF3EventInfo();
271 : void setLHEF3EventInfo( map<string, string> *eventAttributesIn,
272 : map<string, double > *weights_detailedIn,
273 : vector<double > *weights_compressedIn,
274 : LHAscales *scalesIn, LHAweights *weightsIn,
275 : LHArwgt *rwgtIn );
276 :
277 : // Retrieve events tag information.
278 : string getEventAttribute(string key, bool doRemoveWhitespace = false);
279 :
280 : // Retrieve LHEF version
281 : int LHEFversion();
282 :
283 : // Retrieve initrwgt tag information.
284 : unsigned int getInitrwgtSize();
285 :
286 : // Retrieve generator tag information.
287 : unsigned int getGeneratorSize();
288 : string getGeneratorValue(unsigned int n = 0);
289 : string getGeneratorAttribute( unsigned int n, string key,
290 : bool doRemoveWhitespace = false);
291 :
292 : // Retrieve rwgt tag information.
293 : unsigned int getWeightsDetailedSize();
294 : double getWeightsDetailedValue(string n);
295 : string getWeightsDetailedAttribute(string n, string key,
296 : bool doRemoveWhitespace = false);
297 :
298 : // Retrieve weights tag information.
299 : unsigned int getWeightsCompressedSize();
300 : double getWeightsCompressedValue(unsigned int n);
301 : string getWeightsCompressedAttribute(string key,
302 : bool doRemoveWhitespace = false);
303 :
304 : // Retrieve scales tag information.
305 : string getScalesValue(bool doRemoveWhitespace = false);
306 : double getScalesAttribute(string key);
307 :
308 : // Set LHEF headers
309 : void setHeader(const string &key, const string &val)
310 0 : { headers[key] = val; }
311 :
312 : // Set abort in parton level.
313 0 : void setAbortPartonLevel(bool abortIn) {abortPartonLevel = abortIn;}
314 0 : bool getAbortPartonLevel() {return abortPartonLevel;}
315 :
316 : // Get information on hard diffractive events.
317 0 : bool hasUnresolvedBeams() const {return hasUnresBeams;}
318 0 : bool isHardDiffractive() const {return isHardDiffA || isHardDiffB;}
319 : bool isHardDiffractiveA() const {return isHardDiffA;}
320 : bool isHardDiffractiveB() const {return isHardDiffB;}
321 0 : double xPomeronA() const {return xPomA;}
322 0 : double xPomeronB() const {return xPomB;}
323 : double tPomeronA() const {return tPomA;}
324 : double tPomeronB() const {return tPomB;}
325 :
326 : private:
327 :
328 : // Number of times the same error message is repeated, unless overridden.
329 : static const int TIMESTOPRINT;
330 :
331 : // Allow conversion from mb to pb.
332 : static const double CONVERTMB2PB;
333 :
334 : // Store common beam quantities.
335 : int idASave, idBSave;
336 : double pzASave, eASave,mASave, pzBSave, eBSave, mBSave, eCMSave, sSave;
337 :
338 : // Store initialization information.
339 : bool lowPTmin;
340 :
341 : // Store common integrated cross section quantities.
342 : long nTry, nSel, nAcc;
343 : double sigGen, sigErr, wtAccSum;
344 : map<int, string> procNameM;
345 : map<int, long> nTryM, nSelM, nAccM;
346 : map<int, double> sigGenM, sigErrM;
347 : int lhaStrategySave;
348 :
349 : // Store common MPI information.
350 : double a0MPISave;
351 :
352 : // Store current-event quantities.
353 : bool isRes, isDiffA, isDiffB, isDiffC, isND, isLH, hasSubSave[4],
354 : bIsSet, evolIsSet, atEOF, isVal1, isVal2, hasHistorySave,
355 : abortPartonLevel, isHardDiffA, isHardDiffB, hasUnresBeams;
356 : int codeSave, codeSubSave[4], nFinalSave, nFinalSubSave[4], nTotal,
357 : id1Save[4], id2Save[4], id1pdfSave[4], id2pdfSave[4], nMPISave,
358 : nISRSave, nFSRinProcSave, nFSRinResSave;
359 : double x1Save[4], x2Save[4], x1pdfSave[4], x2pdfSave[4], pdf1Save[4],
360 : pdf2Save[4], Q2FacSave[4], alphaEMSave[4], alphaSSave[4],
361 : Q2RenSave[4], scalupSave[4], sH[4], tH[4], uH[4], pTH[4], m3H[4],
362 : m4H[4], thetaH[4], phiH[4], weightSave, bMPISave, enhanceMPISave,
363 : pTmaxMPISave, pTmaxISRSave, pTmaxFSRSave, pTnowSave,
364 : zNowISRSave, pT2NowISRSave, xPomA, xPomB, tPomA, tPomB;
365 : string nameSave, nameSubSave[4];
366 : vector<int> codeMPISave, iAMPISave, iBMPISave;
367 : vector<double> pTMPISave, eMPISave;
368 :
369 : // Vector of various loop counters.
370 : int counters[50];
371 :
372 : // Map for all error messages.
373 : map<string, int> messages;
374 :
375 : // Map for LHEF headers
376 : map<string, string> headers;
377 :
378 : // Map for plugin libraries.
379 : map<string, pair<void*, int> > plugins;
380 :
381 : // Friend classes allowed to set info.
382 : friend class Pythia;
383 : friend class ProcessLevel;
384 : friend class ProcessContainer;
385 : friend class PartonLevel;
386 : friend class MultipartonInteractions;
387 : friend class LHAup;
388 : friend class LHAPDF;
389 : friend class Diffraction;
390 :
391 : // Set info on the two incoming beams: only from Pythia class.
392 : void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
393 0 : idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
394 : void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
395 0 : idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
396 0 : void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
397 :
398 : // Reset info for current event: only from Pythia class.
399 : void clear() {
400 0 : isRes = isDiffA = isDiffB = isDiffC = isND = isLH = bIsSet
401 0 : = evolIsSet = atEOF = isVal1 = isVal2 = hasHistorySave
402 0 : = hasUnresBeams = false;
403 0 : codeSave = nFinalSave = nTotal = nMPISave = nISRSave = nFSRinProcSave
404 0 : = nFSRinResSave = 0;
405 0 : weightSave = bMPISave = enhanceMPISave = weightCKKWLSave = 1.;
406 0 : pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave = zNowISRSave
407 0 : = pT2NowISRSave = weightFIRSTSave = 0.;
408 0 : nameSave = " ";
409 0 : for (int i = 0; i < 4; ++i) {
410 0 : hasSubSave[i] = false;
411 0 : codeSubSave[i] = nFinalSubSave[i] = id1pdfSave[i] = id2pdfSave[i]
412 0 : = id1Save[i] = id2Save[i] = 0;
413 0 : x1pdfSave[i] = x2pdfSave[i] = pdf1Save[i] = pdf2Save[i]
414 0 : = Q2FacSave[i] = alphaEMSave[i] = alphaSSave[i] = Q2RenSave[i]
415 0 : = scalupSave[i] = x1Save[i] = x2Save[i] = sH[i] = tH[i] = uH[i]
416 0 : = pTH[i] = m3H[i] = m4H[i] = thetaH[i] = phiH[i] = 0.;
417 0 : nameSubSave[i] = " ";
418 : }
419 0 : codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
420 0 : pTMPISave.resize(0); eMPISave.resize(0); setHardDiff();}
421 :
422 : // Reset info arrays only for parton and hadron level.
423 0 : int sizeMPIarrays() const { return codeMPISave.size(); }
424 0 : void resizeMPIarrays(int newSize) { codeMPISave.resize(newSize);
425 0 : iAMPISave.resize(newSize); iBMPISave.resize(newSize);
426 0 : pTMPISave.resize(newSize); eMPISave.resize(newSize); }
427 :
428 : // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
429 : // MultipartonInteractions classes.
430 : void setType( string nameIn, int codeIn, int nFinalIn,
431 : bool isNonDiffIn = false, bool isResolvedIn = true,
432 : bool isDiffractiveAin = false, bool isDiffractiveBin = false,
433 : bool isDiffractiveCin = false, bool isLHAin = false) {
434 0 : nameSave = nameIn; codeSave = codeIn; nFinalSave = nFinalIn;
435 0 : isND = isNonDiffIn; isRes = isResolvedIn; isDiffA = isDiffractiveAin;
436 0 : isDiffB = isDiffractiveBin; isDiffC = isDiffractiveCin; isLH = isLHAin;
437 0 : nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave[0] = false;
438 0 : nameSubSave[0] = " "; codeSubSave[0] = 0; nFinalSubSave[0] = 0;
439 0 : evolIsSet = false;}
440 : void setSubType( int iDS, string nameSubIn, int codeSubIn,
441 0 : int nFinalSubIn) { hasSubSave[iDS] = true; nameSubSave[iDS] = nameSubIn;
442 0 : codeSubSave[iDS] = codeSubIn; nFinalSubSave[iDS] = nFinalSubIn;}
443 : void setPDFalpha( int iDS, int id1pdfIn, int id2pdfIn, double x1pdfIn,
444 : double x2pdfIn, double pdf1In, double pdf2In, double Q2FacIn,
445 : double alphaEMIn, double alphaSIn, double Q2RenIn, double scalupIn) {
446 0 : id1pdfSave[iDS] = id1pdfIn; id2pdfSave[iDS] = id2pdfIn;
447 0 : x1pdfSave[iDS] = x1pdfIn; x2pdfSave[iDS] = x2pdfIn;
448 0 : pdf1Save[iDS] = pdf1In; pdf2Save[iDS] = pdf2In; Q2FacSave[iDS] = Q2FacIn;
449 0 : alphaEMSave[iDS] = alphaEMIn; alphaSSave[iDS] = alphaSIn;
450 0 : Q2RenSave[iDS] = Q2RenIn; scalupSave[iDS] = scalupIn;}
451 : void setScalup( int iDS, double scalupIn) {scalupSave[iDS] = scalupIn;}
452 : void setKin( int iDS, int id1In, int id2In, double x1In, double x2In,
453 : double sHatIn, double tHatIn, double uHatIn, double pTHatIn,
454 : double m3HatIn, double m4HatIn, double thetaHatIn, double phiHatIn) {
455 0 : id1Save[iDS] = id1In; id2Save[iDS] = id2In; x1Save[iDS] = x1In;
456 0 : x2Save[iDS] = x2In; sH[iDS] = sHatIn; tH[iDS] = tHatIn;
457 0 : uH[iDS] = uHatIn; pTH[iDS] = pTHatIn; m3H[iDS] = m3HatIn;
458 0 : m4H[iDS] = m4HatIn; thetaH[iDS] = thetaHatIn; phiH[iDS] = phiHatIn;}
459 : void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
460 0 : int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
461 0 : pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
462 0 : iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
463 :
464 : // Set info on cross section: from ProcessLevel (reset in Pythia).
465 0 : void sigmaReset() { nTry = nSel = nAcc = 0; sigGen = sigErr = wtAccSum = 0.;
466 0 : procNameM.clear(); nTryM.clear(); nSelM.clear(); nAccM.clear();
467 0 : sigGenM.clear(); sigErrM.clear();}
468 : void setSigma( int i, string procNameIn, long nTryIn, long nSelIn,
469 : long nAccIn, double sigGenIn, double sigErrIn, double wtAccSumIn) {
470 0 : if (i == 0) {nTry = nTryIn; nSel = nSelIn; nAcc = nAccIn;
471 0 : sigGen = sigGenIn; sigErr = sigErrIn; wtAccSum = wtAccSumIn;}
472 0 : else {procNameM[i] = procNameIn; nTryM[i] = nTryIn; nSelM[i] = nSelIn;
473 0 : nAccM[i] = nAccIn; sigGenM[i] = sigGenIn; sigErrM[i] = sigErrIn;} }
474 : void addSigma( int i, long nTryIn, long nSelIn, long nAccIn, double sigGenIn,
475 0 : double sigErrIn) { nTryM[i] += nTryIn; nSelM[i] += nSelIn;
476 0 : nAccM[i] += nAccIn; sigGenM[i] += sigGenIn;
477 0 : sigErrM[i] = sqrtpos(sigErrM[i]*sigErrM[i] + sigErrIn*sigErrIn); }
478 :
479 : // Set info on impact parameter: from PartonLevel.
480 : void setImpact( double bMPIIn, double enhanceMPIIn, bool bIsSetIn = true) {
481 0 : bMPISave = bMPIIn; enhanceMPISave = eMPISave[0] = enhanceMPIIn,
482 0 : bIsSet = bIsSetIn;}
483 :
484 : // Set info on pTmax scales and number of evolution steps: from PartonLevel.
485 : void setPartEvolved( int nMPIIn, int nISRIn) {
486 0 : nMPISave = nMPIIn; nISRSave = nISRIn;}
487 : void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
488 : int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
489 0 : pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
490 0 : pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
491 0 : nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
492 0 : evolIsSet = true;}
493 :
494 : // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
495 0 : void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
496 :
497 : // Set a0 from MultipartonInteractions.
498 0 : void seta0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
499 :
500 : // Set info whether reading of Les Houches Accord file at end.
501 0 : void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
502 :
503 : // Set event weight; currently only for Les Houches description.
504 : void setWeight( double weightIn, int lhaStrategyIn) {
505 0 : weightSave = weightIn; lhaStrategySave = lhaStrategyIn; }
506 :
507 : // Save merging weight (i.e. CKKW-L-type weight, summed O(\alpha_s) weight)
508 : double weightCKKWLSave, weightFIRSTSave;
509 :
510 : // Set info on resolved processes.
511 : void setIsResolved(bool isResIn) {isRes = isResIn;}
512 :
513 : // Set info on hard diffraction.
514 : void setHardDiff( bool hasUnresBeamsIn = false,
515 : bool isHardDiffAIn = false, bool isHardDiffBIn = false,
516 : double xPomAIn = 0., double xPomBIn = 0., double tPomAIn = 0.,
517 0 : double tPomBIn = 0.) { hasUnresBeams = hasUnresBeamsIn;
518 0 : isHardDiffA = isHardDiffAIn; isHardDiffB = isHardDiffBIn;
519 0 : xPomA = xPomAIn; xPomB = xPomBIn;
520 0 : tPomA = tPomAIn; tPomB = tPomBIn;}
521 :
522 : // Set information in hard diffractive events.
523 : void setHasUnresolvedBeams(bool hasUnresBeamsIn)
524 0 : {hasUnresBeams = hasUnresBeamsIn;}
525 :
526 : };
527 :
528 : //==========================================================================
529 :
530 : } // end namespace Pythia8
531 :
532 : #endif // Pythia8_Info_H
|