Line data Source code
1 : /**************************************************************************
2 : * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 : * *
4 : * Author: The ALICE Off-line Project. *
5 : * Contributors are mentioned in the code where appropriate. *
6 : * *
7 : * Permission to use, copy, modify and distribute this software and its *
8 : * documentation strictly for non-commercial purposes is hereby granted *
9 : * without fee, provided that the above copyright notice appears in all *
10 : * copies and that both the copyright notice and this permission notice *
11 : * appear in the supporting documentation. The authors make no claims *
12 : * about the suitability of this software for any purpose. It is *
13 : * provided "as is" without express or implied warranty. *
14 : **************************************************************************/
15 :
16 : /* $Id$ */
17 :
18 : //
19 : // Generator using the TPythia interface (via AliPythia)
20 : // to generate pp collisions.
21 : // Using SetNuclei() also nuclear modifications to the structure functions
22 : // can be taken into account. This makes, of course, only sense for the
23 : // generation of the products of hard processes (heavy flavor, jets ...)
24 : //
25 : // andreas.morsch@cern.ch
26 : //
27 :
28 : #include <TMath.h>
29 : #include <TClonesArray.h>
30 : #include <TDatabasePDG.h>
31 : #include <TParticle.h>
32 : #include <TPDGCode.h>
33 : #include <TObjArray.h>
34 : #include <TSystem.h>
35 : #include <TTree.h>
36 : #include "AliConst.h"
37 : #include "AliDecayerPythia.h"
38 : #include "AliGenPythia.h"
39 : #include "AliFastGlauber.h"
40 : #include "AliHeader.h"
41 : #include "AliGenPythiaEventHeader.h"
42 : #include "AliPythia.h"
43 : #include "AliPythiaRndm.h"
44 : #include "AliRun.h"
45 : #include "AliStack.h"
46 : #include "AliRunLoader.h"
47 : #include "AliMC.h"
48 : #include "AliLog.h"
49 : #include "PyquenCommon.h"
50 :
51 2 : ClassImp(AliGenPythia)
52 :
53 :
54 : AliGenPythia::AliGenPythia():
55 0 : AliGenMC(),
56 0 : fProcess(kPyCharm),
57 0 : fItune(-1),
58 0 : fStrucFunc(kCTEQ5L),
59 0 : fKineBias(0.),
60 0 : fTrials(0),
61 0 : fTrialsRun(0),
62 0 : fQ(0.),
63 0 : fX1(0.),
64 0 : fX2(0.),
65 0 : fEventTime(0.),
66 0 : fInteractionRate(0.),
67 0 : fTimeWindow(0.),
68 0 : fCurSubEvent(0),
69 0 : fEventsTime(0),
70 0 : fNev(0),
71 0 : fFlavorSelect(0),
72 0 : fXsection(0.),
73 0 : fPythia(0),
74 0 : fWeightPower(0.),
75 0 : fPtHardMin(0.),
76 0 : fPtHardMax(1.e4),
77 0 : fYHardMin(-1.e10),
78 0 : fYHardMax(1.e10),
79 0 : fGinit(1),
80 0 : fGfinal(1),
81 0 : fCRoff(0),
82 0 : fHadronisation(1),
83 0 : fPatchOmegaDalitz(0),
84 0 : fDecayerExodus(0),
85 0 : fNpartons(0),
86 0 : fReadFromFile(0),
87 0 : fReadLHEF(0),
88 0 : fQuench(0),
89 0 : fQhat(0.),
90 0 : fLength(0.),
91 0 : fpyquenT(1.),
92 0 : fpyquenTau(0.1),
93 0 : fpyquenNf(0),
94 0 : fpyquenEloss(0),
95 0 : fpyquenAngle(0),
96 0 : fImpact(0.),
97 0 : fPtKick(1.),
98 0 : fFullEvent(kTRUE),
99 0 : fDecayer(new AliDecayerPythia()),
100 0 : fDebugEventFirst(-1),
101 0 : fDebugEventLast(-1),
102 0 : fEtMinJet(0.),
103 0 : fEtMaxJet(1.e4),
104 0 : fEtaMinJet(-20.),
105 0 : fEtaMaxJet(20.),
106 0 : fPhiMinJet(0.),
107 0 : fPhiMaxJet(2.* TMath::Pi()),
108 0 : fJetReconstruction(kCell),
109 0 : fEtaMinGamma(-20.),
110 0 : fEtaMaxGamma(20.),
111 0 : fPhiMinGamma(0.),
112 0 : fPhiMaxGamma(2. * TMath::Pi()),
113 0 : fUseYCutHQ(kFALSE),
114 0 : fYMinHQ(-20.),
115 0 : fYMaxHQ(20.),
116 0 : fPycellEtaMax(2.),
117 0 : fPycellNEta(274),
118 0 : fPycellNPhi(432),
119 0 : fPycellThreshold(0.),
120 0 : fPycellEtSeed(4.),
121 0 : fPycellMinEtJet(10.),
122 0 : fPycellMaxRadius(1.),
123 0 : fStackFillOpt(kFlavorSelection),
124 0 : fFeedDownOpt(kTRUE),
125 0 : fFragmentation(kTRUE),
126 0 : fSetNuclei(kFALSE),
127 0 : fUseNuclearPDF(kFALSE),
128 0 : fUseLorentzBoost(kTRUE),
129 0 : fNewMIS(kFALSE),
130 0 : fHFoff(kFALSE),
131 0 : fNucPdf(0),
132 0 : fTriggerParticle(0),
133 0 : fTriggerEta(0.9),
134 0 : fTriggerY(999.),
135 0 : fTriggerEtaMin(0.9),
136 0 : fTriggerMinPt(-1),
137 0 : fTriggerMaxPt(1000),
138 0 : fTriggerMultiplicity(0),
139 0 : fTriggerMultiplicityEta(0),
140 0 : fTriggerMultiplicityEtaMin(0),
141 0 : fTriggerMultiplicityEtaMax(0),
142 0 : fTriggerMultiplicityPtMin(0),
143 0 : fCountMode(kCountAll),
144 0 : fHeader(0),
145 0 : fRL(0),
146 0 : fkFileName(0),
147 0 : fkNameLHEF(0),
148 0 : fFragPhotonInCalo(kFALSE),
149 0 : fHadronInCalo(kFALSE) ,
150 0 : fPi0InCalo(kFALSE) ,
151 0 : fEtaInCalo(kFALSE) ,
152 0 : fPhotonInCalo(kFALSE), // not in use
153 0 : fDecayPhotonInCalo(kFALSE),
154 0 : fForceNeutralMeson2PhotonDecay(kFALSE),
155 0 : fEleInCalo(kFALSE),
156 0 : fEleInEMCAL(kFALSE), // not in use
157 0 : fCheckBarrel(kFALSE),
158 0 : fCheckEMCAL(kFALSE),
159 0 : fCheckPHOS(kFALSE),
160 0 : fCheckPHOSeta(kFALSE),
161 0 : fPHOSRotateCandidate(-1),
162 0 : fTriggerParticleMinPt(0),
163 0 : fPhotonMinPt(0), // not in use
164 0 : fElectronMinPt(0), // not in use
165 0 : fPHOSMinPhi(219.),
166 0 : fPHOSMaxPhi(321.),
167 0 : fPHOSEta(0.13),
168 0 : fEMCALMinPhi(79.),
169 0 : fEMCALMaxPhi(191.),
170 0 : fEMCALEta(0.71),
171 0 : fkTuneForDiff(0),
172 0 : fProcDiff(0)
173 0 : {
174 : // Default Constructor
175 0 : fEnergyCMS = 5500.;
176 0 : if (!AliPythiaRndm::GetPythiaRandom())
177 0 : AliPythiaRndm::SetPythiaRandom(GetRandom());
178 0 : }
179 :
180 : AliGenPythia::AliGenPythia(Int_t npart)
181 0 : :AliGenMC(npart),
182 0 : fProcess(kPyCharm),
183 0 : fItune(-1),
184 0 : fStrucFunc(kCTEQ5L),
185 0 : fKineBias(0.),
186 0 : fTrials(0),
187 0 : fTrialsRun(0),
188 0 : fQ(0.),
189 0 : fX1(0.),
190 0 : fX2(0.),
191 0 : fEventTime(0.),
192 0 : fInteractionRate(0.),
193 0 : fTimeWindow(0.),
194 0 : fCurSubEvent(0),
195 0 : fEventsTime(0),
196 0 : fNev(0),
197 0 : fFlavorSelect(0),
198 0 : fXsection(0.),
199 0 : fPythia(0),
200 0 : fWeightPower(0.),
201 0 : fPtHardMin(0.),
202 0 : fPtHardMax(1.e4),
203 0 : fYHardMin(-1.e10),
204 0 : fYHardMax(1.e10),
205 0 : fGinit(kTRUE),
206 0 : fGfinal(kTRUE),
207 0 : fCRoff(kFALSE),
208 0 : fHadronisation(kTRUE),
209 0 : fPatchOmegaDalitz(0),
210 0 : fDecayerExodus(0),
211 0 : fNpartons(0),
212 0 : fReadFromFile(kFALSE),
213 0 : fReadLHEF(0),
214 0 : fQuench(kFALSE),
215 0 : fQhat(0.),
216 0 : fLength(0.),
217 0 : fpyquenT(1.),
218 0 : fpyquenTau(0.1),
219 0 : fpyquenNf(0),
220 0 : fpyquenEloss(0),
221 0 : fpyquenAngle(0),
222 0 : fImpact(0.),
223 0 : fPtKick(1.),
224 0 : fFullEvent(kTRUE),
225 0 : fDecayer(new AliDecayerPythia()),
226 0 : fDebugEventFirst(-1),
227 0 : fDebugEventLast(-1),
228 0 : fEtMinJet(0.),
229 0 : fEtMaxJet(1.e4),
230 0 : fEtaMinJet(-20.),
231 0 : fEtaMaxJet(20.),
232 0 : fPhiMinJet(0.),
233 0 : fPhiMaxJet(2.* TMath::Pi()),
234 0 : fJetReconstruction(kCell),
235 0 : fEtaMinGamma(-20.),
236 0 : fEtaMaxGamma(20.),
237 0 : fPhiMinGamma(0.),
238 0 : fPhiMaxGamma(2. * TMath::Pi()),
239 0 : fUseYCutHQ(kFALSE),
240 0 : fYMinHQ(-20.),
241 0 : fYMaxHQ(20.),
242 0 : fPycellEtaMax(2.),
243 0 : fPycellNEta(274),
244 0 : fPycellNPhi(432),
245 0 : fPycellThreshold(0.),
246 0 : fPycellEtSeed(4.),
247 0 : fPycellMinEtJet(10.),
248 0 : fPycellMaxRadius(1.),
249 0 : fStackFillOpt(kFlavorSelection),
250 0 : fFeedDownOpt(kTRUE),
251 0 : fFragmentation(kTRUE),
252 0 : fSetNuclei(kFALSE),
253 0 : fUseNuclearPDF(kFALSE),
254 0 : fUseLorentzBoost(kTRUE),
255 0 : fNewMIS(kFALSE),
256 0 : fHFoff(kFALSE),
257 0 : fNucPdf(0),
258 0 : fTriggerParticle(0),
259 0 : fTriggerEta(0.9),
260 0 : fTriggerY(999.),
261 0 : fTriggerEtaMin(0.9),
262 0 : fTriggerMinPt(-1),
263 0 : fTriggerMaxPt(1000),
264 0 : fTriggerMultiplicity(0),
265 0 : fTriggerMultiplicityEta(0),
266 0 : fTriggerMultiplicityEtaMin(0),
267 0 : fTriggerMultiplicityEtaMax(0),
268 0 : fTriggerMultiplicityPtMin(0),
269 0 : fCountMode(kCountAll),
270 0 : fHeader(0),
271 0 : fRL(0),
272 0 : fkFileName(0),
273 0 : fkNameLHEF(0),
274 0 : fFragPhotonInCalo(kFALSE),
275 0 : fHadronInCalo(kFALSE) ,
276 0 : fPi0InCalo(kFALSE) ,
277 0 : fEtaInCalo(kFALSE) ,
278 0 : fPhotonInCalo(kFALSE), // not in use
279 0 : fDecayPhotonInCalo(kFALSE),
280 0 : fForceNeutralMeson2PhotonDecay(kFALSE),
281 0 : fEleInCalo(kFALSE),
282 0 : fEleInEMCAL(kFALSE), // not in use
283 0 : fCheckBarrel(kFALSE),
284 0 : fCheckEMCAL(kFALSE),
285 0 : fCheckPHOS(kFALSE),
286 0 : fCheckPHOSeta(kFALSE),
287 0 : fPHOSRotateCandidate(-1),
288 0 : fTriggerParticleMinPt(0),
289 0 : fPhotonMinPt(0), // not in use
290 0 : fElectronMinPt(0), // not in use
291 0 : fPHOSMinPhi(219.),
292 0 : fPHOSMaxPhi(321.),
293 0 : fPHOSEta(0.13),
294 0 : fEMCALMinPhi(79.),
295 0 : fEMCALMaxPhi(191.),
296 0 : fEMCALEta(0.71),
297 0 : fkTuneForDiff(0),
298 0 : fProcDiff(0)
299 0 : {
300 : // default charm production at 5. 5 TeV
301 : // semimuonic decay
302 : // structure function GRVHO
303 : //
304 0 : fEnergyCMS = 5500.;
305 0 : fName = "Pythia";
306 0 : fTitle= "Particle Generator using PYTHIA";
307 0 : SetForceDecay();
308 : // Set random number generator
309 0 : if (!AliPythiaRndm::GetPythiaRandom())
310 0 : AliPythiaRndm::SetPythiaRandom(GetRandom());
311 0 : }
312 :
313 0 : AliGenPythia::~AliGenPythia()
314 0 : {
315 : // Destructor
316 0 : if(fEventsTime) delete fEventsTime;
317 0 : }
318 :
319 : void AliGenPythia::SetInteractionRate(Float_t rate,Float_t timewindow)
320 : {
321 : // Generate pileup using user specified rate
322 0 : fInteractionRate = rate;
323 0 : fTimeWindow = timewindow;
324 0 : GeneratePileup();
325 0 : }
326 :
327 : void AliGenPythia::GeneratePileup()
328 : {
329 : // Generate sub events time for pileup
330 0 : fEventsTime = 0;
331 0 : if(fInteractionRate == 0.) {
332 0 : Warning("GeneratePileup","Zero interaction specified. Skipping pileup generation.\n");
333 0 : return;
334 : }
335 :
336 0 : Int_t npart = NumberParticles();
337 0 : if(npart < 0) {
338 0 : Warning("GeneratePileup","Negative number of particles. Skipping pileup generation.\n");
339 0 : return;
340 : }
341 :
342 0 : if(fEventsTime) delete fEventsTime;
343 0 : fEventsTime = new TArrayF(npart);
344 : TArrayF &array = *fEventsTime;
345 0 : for(Int_t ipart = 0; ipart < npart; ipart++)
346 0 : array[ipart] = 0.;
347 :
348 : Float_t eventtime = 0.;
349 0 : while(1)
350 : {
351 0 : eventtime += (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate);
352 0 : if(eventtime > fTimeWindow) break;
353 0 : array.Set(array.GetSize()+1);
354 0 : array[array.GetSize()-1] = eventtime;
355 : }
356 :
357 : eventtime = 0.;
358 0 : while(1)
359 : {
360 0 : eventtime -= (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate);
361 0 : if(TMath::Abs(eventtime) > fTimeWindow) break;
362 0 : array.Set(array.GetSize()+1);
363 0 : array[array.GetSize()-1] = eventtime;
364 : }
365 :
366 0 : SetNumberParticles(fEventsTime->GetSize());
367 0 : }
368 :
369 : void AliGenPythia::SetPycellParameters(Float_t etamax, Int_t neta, Int_t nphi,
370 : Float_t thresh, Float_t etseed, Float_t minet, Float_t r)
371 : {
372 : // Set pycell parameters
373 0 : fPycellEtaMax = etamax;
374 0 : fPycellNEta = neta;
375 0 : fPycellNPhi = nphi;
376 0 : fPycellThreshold = thresh;
377 0 : fPycellEtSeed = etseed;
378 0 : fPycellMinEtJet = minet;
379 0 : fPycellMaxRadius = r;
380 0 : }
381 :
382 :
383 :
384 : void AliGenPythia::SetEventListRange(Int_t eventFirst, Int_t eventLast)
385 : {
386 : // Set a range of event numbers, for which a table
387 : // of generated particle will be printed
388 0 : fDebugEventFirst = eventFirst;
389 0 : fDebugEventLast = eventLast;
390 0 : if (fDebugEventLast==-1) fDebugEventLast=fDebugEventFirst;
391 0 : }
392 :
393 : void AliGenPythia::Init()
394 : {
395 : // Initialisation
396 :
397 0 : SetMC(AliPythia::Instance());
398 0 : fPythia=(AliPythia*) fMCEvGen;
399 :
400 : //
401 0 : fParentWeight=1./Float_t(fNpart);
402 : //
403 0 : if (fWeightPower != 0)
404 0 : fPythia->SetWeightPower(fWeightPower);
405 0 : fPythia->SetCKIN(3,fPtHardMin);
406 0 : fPythia->SetCKIN(4,fPtHardMax);
407 0 : fPythia->SetCKIN(7,fYHardMin);
408 0 : fPythia->SetCKIN(8,fYHardMax);
409 :
410 0 : if (fProjectile != "p" || fTarget != "p") fPythia->SetCollisionSystem(fProjectile,fTarget);
411 :
412 0 : if(fUseNuclearPDF)
413 0 : fPythia->SetNuclei(fAProjectile, fATarget, fNucPdf);
414 : // Fragmentation?
415 0 : if (fFragmentation) {
416 0 : fPythia->SetMSTP(111,1);
417 0 : } else {
418 0 : fPythia->SetMSTP(111,0);
419 : }
420 :
421 :
422 : // initial state radiation
423 0 : fPythia->SetMSTP(61,fGinit);
424 : // final state radiation
425 0 : fPythia->SetMSTP(71,fGfinal);
426 : //color reconnection strength
427 0 : if(fCRoff==1)fPythia->SetMSTP(95,0);
428 : // pt - kick
429 0 : if (fPtKick > 0.) {
430 0 : fPythia->SetMSTP(91,1);
431 0 : fPythia->SetPARP(91,fPtKick);
432 0 : fPythia->SetPARP(93, 4. * fPtKick);
433 0 : } else {
434 0 : fPythia->SetMSTP(91,0);
435 : }
436 :
437 0 : if (fReadLHEF) fPythia->OpenFortranFile(97, const_cast<char*>(fkNameLHEF));
438 :
439 0 : if (fReadFromFile) {
440 0 : fRL = AliRunLoader::Open(fkFileName, "Partons");
441 0 : fRL->LoadKinematics();
442 0 : fRL->LoadHeader();
443 0 : } else {
444 0 : fRL = 0x0;
445 : }
446 : //
447 0 : fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc, fItune);
448 : // Forward Paramters to the AliPythia object
449 0 : fDecayer->SetForceDecay(fForceDecay);
450 : // Switch off Heavy Flavors on request
451 0 : if (fHFoff) {
452 : // Maximum number of quark flavours used in pdf
453 0 : fPythia->SetMSTP(58, 3);
454 : // Maximum number of flavors that can be used in showers
455 0 : fPythia->SetMSTJ(45, 3);
456 : // Switch off g->QQbar splitting in decay table
457 0 : ((AliDecayerPythia*) fDecayer)->HeavyFlavourOff();
458 0 : }
459 :
460 0 : fDecayer->Init();
461 :
462 :
463 : // Parent and Children Selection
464 0 : switch (fProcess)
465 : {
466 : case kPyOldUEQ2ordered:
467 : case kPyOldUEQ2ordered2:
468 : case kPyOldPopcorn:
469 : break;
470 : case kPyCharm:
471 : case kPyCharmUnforced:
472 : case kPyCharmPbPbMNR:
473 : case kPyCharmpPbMNR:
474 : case kPyCharmppMNR:
475 : case kPyCharmppMNRwmi:
476 : case kPyCharmPWHG:
477 0 : fParentSelect[0] = 411;
478 0 : fParentSelect[1] = 421;
479 0 : fParentSelect[2] = 431;
480 0 : fParentSelect[3] = 4122;
481 0 : fParentSelect[4] = 4232;
482 0 : fParentSelect[5] = 4132;
483 0 : fParentSelect[6] = 4332;
484 0 : fFlavorSelect = 4;
485 0 : break;
486 : case kPyD0PbPbMNR:
487 : case kPyD0pPbMNR:
488 : case kPyD0ppMNR:
489 0 : fParentSelect[0] = 421;
490 0 : fFlavorSelect = 4;
491 0 : break;
492 : case kPyDPlusPbPbMNR:
493 : case kPyDPluspPbMNR:
494 : case kPyDPlusppMNR:
495 0 : fParentSelect[0] = 411;
496 0 : fFlavorSelect = 4;
497 0 : break;
498 : case kPyDPlusStrangePbPbMNR:
499 : case kPyDPlusStrangepPbMNR:
500 : case kPyDPlusStrangeppMNR:
501 0 : fParentSelect[0] = 431;
502 0 : fFlavorSelect = 4;
503 0 : break;
504 : case kPyLambdacppMNR:
505 0 : fParentSelect[0] = 4122;
506 0 : fFlavorSelect = 4;
507 0 : break;
508 : case kPyBeauty:
509 : case kPyBeautyJets:
510 : case kPyBeautyPbPbMNR:
511 : case kPyBeautypPbMNR:
512 : case kPyBeautyppMNR:
513 : case kPyBeautyppMNRwmi:
514 : case kPyBeautyPWHG:
515 0 : fParentSelect[0]= 511;
516 0 : fParentSelect[1]= 521;
517 0 : fParentSelect[2]= 531;
518 0 : fParentSelect[3]= 5122;
519 0 : fParentSelect[4]= 5132;
520 0 : fParentSelect[5]= 5232;
521 0 : fParentSelect[6]= 5332;
522 0 : fFlavorSelect = 5;
523 0 : break;
524 : case kPyBeautyUnforced:
525 0 : fParentSelect[0] = 511;
526 0 : fParentSelect[1] = 521;
527 0 : fParentSelect[2] = 531;
528 0 : fParentSelect[3] = 5122;
529 0 : fParentSelect[4] = 5132;
530 0 : fParentSelect[5] = 5232;
531 0 : fParentSelect[6] = 5332;
532 0 : fFlavorSelect = 5;
533 0 : break;
534 : case kPyJpsiChi:
535 : case kPyJpsi:
536 0 : fParentSelect[0] = 443;
537 0 : break;
538 : case kPyMbDefault:
539 : case kPyMbAtlasTuneMC09:
540 : case kPyMb:
541 : case kPyMbWithDirectPhoton:
542 : case kPyMbNonDiffr:
543 : case kPyMbMSEL1:
544 : case kPyJets:
545 : case kPyJetsPWHG:
546 : case kPyDirectGamma:
547 : case kPyLhwgMb:
548 : break;
549 : case kPyWPWHG:
550 : case kPyW:
551 : case kPyZ:
552 : case kPyZgamma:
553 : case kPyMBRSingleDiffraction:
554 : case kPyMBRDoubleDiffraction:
555 : case kPyMBRCentralDiffraction:
556 : break;
557 : }
558 : //
559 : //
560 : // JetFinder for Trigger
561 : //
562 : // Configure detector (EMCAL like)
563 : //
564 0 : fPythia->SetPARU(51, fPycellEtaMax);
565 0 : fPythia->SetMSTU(51, fPycellNEta);
566 0 : fPythia->SetMSTU(52, fPycellNPhi);
567 : //
568 : // Configure Jet Finder
569 : //
570 0 : fPythia->SetPARU(58, fPycellThreshold);
571 0 : fPythia->SetPARU(52, fPycellEtSeed);
572 0 : fPythia->SetPARU(53, fPycellMinEtJet);
573 0 : fPythia->SetPARU(54, fPycellMaxRadius);
574 0 : fPythia->SetMSTU(54, 2);
575 : //
576 : // This counts the total number of calls to Pyevnt() per run.
577 0 : fTrialsRun = 0;
578 0 : fQ = 0.;
579 0 : fX1 = 0.;
580 0 : fX2 = 0.;
581 0 : fNev = 0 ;
582 : //
583 : //
584 : //
585 0 : AliGenMC::Init();
586 :
587 : // Reset Lorentz boost if demanded
588 0 : if(!fUseLorentzBoost) {
589 0 : fDyBoost = 0;
590 0 : Warning("Init","Demand to discard Lorentz boost.\n");
591 0 : }
592 : //
593 : //
594 : //
595 0 : if (fSetNuclei) {
596 0 : fDyBoost = 0;
597 0 : Warning("Init","Deprecated function SetNuclei() used (nPDFs + no boost). Use SetProjectile + SetTarget + SetUseNuclearPDF + SetUseLorentzBoost instead.\n");
598 0 : }
599 0 : fPythia->SetPARJ(200, 0.0);
600 0 : fPythia->SetPARJ(199, 0.0);
601 0 : fPythia->SetPARJ(198, 0.0);
602 0 : fPythia->SetPARJ(197, 0.0);
603 :
604 0 : if (fQuench == 1) {
605 0 : fPythia->InitQuenching(0., 0.1, 0.6e6, 0);
606 0 : }
607 :
608 0 : if(fQuench ==2){fPythia->SetPyquenParameters(fpyquenT,fpyquenTau,fpyquenNf,fpyquenEloss,fpyquenAngle);}
609 :
610 0 : if (fQuench == 3) {
611 : // Nestor's change of the splittings
612 0 : fPythia->SetPARJ(200, 0.8);
613 0 : fPythia->SetMSTJ(41, 1); // QCD radiation only
614 0 : fPythia->SetMSTJ(42, 2); // angular ordering
615 0 : fPythia->SetMSTJ(44, 2); // option to run alpha_s
616 0 : fPythia->SetMSTJ(47, 0); // No correction back to hard scattering element
617 0 : fPythia->SetMSTJ(50, 0); // No coherence in first branching
618 0 : fPythia->SetPARJ(82, 1.); // Cut off for parton showers
619 0 : } else if (fQuench == 4) {
620 : // Armesto-Cunqueiro-Salgado change of the splittings.
621 0 : AliFastGlauber* glauber = AliFastGlauber::Instance();
622 0 : glauber->Init(2);
623 : //read and store transverse almonds corresponding to differnt
624 : //impact parameters.
625 0 : glauber->SetCentralityClass(0.,0.1);
626 0 : fPythia->SetPARJ(200, 1.);
627 0 : fPythia->SetPARJ(198, fQhat);
628 0 : fPythia->SetPARJ(199, fLength);
629 0 : fPythia->SetMSTJ(42, 2); // angular ordering
630 0 : fPythia->SetMSTJ(44, 2); // option to run alpha_s
631 0 : fPythia->SetPARJ(82, 1.); // Cut off for parton showers
632 0 : }
633 :
634 0 : if ( AliLog::GetDebugLevel("","AliGenPythia") >= 1 ) {
635 0 : fPythia->Pystat(4);
636 0 : fPythia->Pystat(5);
637 0 : }
638 0 : }
639 :
640 : void AliGenPythia::Generate()
641 : {
642 : // Generate one event
643 0 : if (!fPythia) fPythia=(AliPythia*) fMCEvGen;
644 0 : fDecayer->ForceDecay();
645 :
646 : Double_t polar[3] = {0,0,0};
647 : Double_t origin[3] = {0,0,0};
648 : Double_t p[4];
649 : // converts from mm/c to s
650 0 : const Float_t kconv=0.001/TMath::C();
651 : //
652 0 : Int_t nt=0;
653 : Int_t jev=0;
654 : Int_t j, kf;
655 0 : fTrials=0;
656 0 : fEventTime = 0.;
657 :
658 :
659 :
660 : // Set collision vertex position
661 0 : if (fVertexSmear == kPerEvent) Vertex();
662 :
663 : // event loop
664 : while(1)
665 : {
666 : //
667 : // Produce event
668 : //
669 : //
670 : // Switch hadronisation off
671 : //
672 0 : fPythia->SetMSTJ(1, 0);
673 :
674 0 : if (fQuench ==4){
675 0 : Double_t bimp;
676 : // Quenching comes through medium-modified splitting functions.
677 0 : AliFastGlauber::Instance()->GetRandomBHard(bimp);
678 0 : fPythia->SetPARJ(197, bimp);
679 0 : fImpact = bimp;
680 0 : fPythia->Qpygin0();
681 0 : }
682 : //
683 : // Either produce new event or read partons from file
684 : //
685 0 : if (!fReadFromFile) {
686 0 : if (!fNewMIS) {
687 0 : fPythia->Pyevnt();
688 0 : } else {
689 0 : fPythia->Pyevnw();
690 : }
691 0 : fNpartons = fPythia->GetN();
692 0 : } else {
693 0 : printf("Loading Event %d\n",AliRunLoader::Instance()->GetEventNumber());
694 0 : fRL->GetEvent(AliRunLoader::Instance()->GetEventNumber());
695 0 : fPythia->SetN(0);
696 0 : LoadEvent(fRL->Stack(), 0 , 1);
697 0 : fPythia->Pyedit(21);
698 : }
699 :
700 : //
701 : // Run quenching routine
702 : //
703 0 : if (fQuench == 1) {
704 0 : fPythia->Quench();
705 0 : } else if (fQuench == 2){
706 0 : fPythia->Pyquen(208., 0, 0.);
707 0 : } else if (fQuench == 3) {
708 : // Quenching is via multiplicative correction of the splittings
709 : }
710 :
711 : //
712 : // Switch hadronisation on
713 : //
714 0 : if (fHadronisation) {
715 0 : fPythia->SetMSTJ(1, 1);
716 : //
717 : // .. and perform hadronisation
718 : // printf("Calling hadronisation %d\n", fPythia->GetN());
719 :
720 0 : if (fPatchOmegaDalitz) {
721 0 : fPythia->SetMDCY(fPythia->Pycomp(111) ,1, 0);
722 0 : fPythia->Pyexec();
723 0 : fPythia->DalitzDecays();
724 0 : fPythia->SetMDCY(fPythia->Pycomp(111) ,1, 1);
725 0 : }
726 :
727 0 : else if (fDecayerExodus) {
728 :
729 0 : fPythia->SetMDCY(fPythia->Pycomp(22) ,1, 0);
730 0 : fPythia->SetMDCY(fPythia->Pycomp(111) ,1, 0);
731 0 : fPythia->SetMDCY(fPythia->Pycomp(221) ,1, 0);
732 0 : fPythia->Pyexec();
733 0 : fPythia->OmegaDalitz();
734 0 : fPythia->SetMDCY(fPythia->Pycomp(111) ,1, 1);
735 0 : fPythia->PizeroDalitz();
736 0 : fPythia->PhiDalitz();
737 0 : fPythia->SetMDCY(fPythia->Pycomp(221) ,1, 1);
738 0 : fPythia->EtaDalitz();
739 0 : fPythia->EtaprimeDalitz();
740 0 : fPythia->SetMDCY(fPythia->Pycomp(22) ,1, 1);
741 0 : fPythia->RhoDirect();
742 0 : fPythia->OmegaDirect();
743 0 : fPythia->PhiDirect();
744 0 : fPythia->JPsiDirect();
745 0 : }
746 :
747 0 : fPythia->Pyexec();
748 0 : }
749 0 : fTrials++;
750 0 : fPythia->ImportParticles(&fParticles,"All");
751 :
752 0 : if (TMath::Abs(fDyBoost) > 1.e-4) Boost();
753 0 : if(TMath::Abs(fXingAngleY) > 1.e-10) BeamCrossAngle();
754 :
755 : //
756 : //
757 : //
758 : Int_t i;
759 :
760 0 : fNprimaries = 0;
761 0 : Int_t np = fParticles.GetEntriesFast();
762 :
763 0 : if (np == 0) continue;
764 : //
765 :
766 : //
767 0 : Int_t* pParent = new Int_t[np];
768 0 : Int_t* pSelected = new Int_t[np];
769 0 : Int_t* trackIt = new Int_t[np];
770 0 : for (i = 0; i < np; i++) {
771 0 : pParent[i] = -1;
772 0 : pSelected[i] = 0;
773 0 : trackIt[i] = 0;
774 : }
775 :
776 : Int_t nc = 0; // Total n. of selected particles
777 : Int_t nParents = 0; // Selected parents
778 : Int_t nTkbles = 0; // Trackable particles
779 0 : if (fProcess != kPyMbDefault &&
780 0 : fProcess != kPyMb &&
781 0 : fProcess != kPyMbAtlasTuneMC09 &&
782 0 : fProcess != kPyMbWithDirectPhoton &&
783 0 : fProcess != kPyJets &&
784 0 : fProcess != kPyDirectGamma &&
785 0 : fProcess != kPyMbNonDiffr &&
786 0 : fProcess != kPyMbMSEL1 &&
787 0 : fProcess != kPyW &&
788 0 : fProcess != kPyZ &&
789 0 : fProcess != kPyJpsi &&
790 0 : fProcess != kPyZgamma &&
791 0 : fProcess != kPyCharmppMNRwmi &&
792 0 : fProcess != kPyBeautyppMNRwmi &&
793 0 : fProcess != kPyBeautyJets &&
794 0 : fProcess != kPyWPWHG &&
795 0 : fProcess != kPyJetsPWHG &&
796 0 : fProcess != kPyCharmPWHG &&
797 0 : fProcess != kPyBeautyPWHG) {
798 :
799 0 : for (i = 0; i < np; i++) {
800 0 : TParticle* iparticle = (TParticle *) fParticles.At(i);
801 0 : Int_t ks = iparticle->GetStatusCode();
802 0 : kf = CheckPDGCode(iparticle->GetPdgCode());
803 : // No initial state partons
804 0 : if (ks==21) continue;
805 : //
806 : // Heavy Flavor Selection
807 : //
808 : // quark ?
809 0 : kf = TMath::Abs(kf);
810 : Int_t kfl = kf;
811 : // Resonance
812 :
813 0 : if (kfl > 100000) kfl %= 100000;
814 0 : if (kfl > 10000) kfl %= 10000;
815 : // meson ?
816 0 : if (kfl > 10) kfl/=100;
817 : // baryon
818 0 : if (kfl > 10) kfl/=10;
819 0 : Int_t ipa = iparticle->GetFirstMother()-1;
820 : Int_t kfMo = 0;
821 : //
822 : // Establish mother daughter relation between heavy quarks and mesons
823 : //
824 0 : if (kf >= fFlavorSelect && kf <= 6) {
825 0 : Int_t idau = iparticle->GetFirstDaughter() - 1;
826 0 : if (idau > -1) {
827 0 : TParticle* daughter = (TParticle *) fParticles.At(idau);
828 0 : Int_t pdgD = daughter->GetPdgCode();
829 0 : if (pdgD == 91 || pdgD == 92) {
830 0 : Int_t jmin = daughter->GetFirstDaughter() - 1;
831 0 : Int_t jmax = daughter->GetLastDaughter() - 1;
832 0 : for (Int_t jp = jmin; jp <= jmax; jp++)
833 0 : ((TParticle *) fParticles.At(jp))->SetFirstMother(i+1);
834 0 : } // is string or cluster
835 0 : } // has daughter
836 0 : } // heavy quark
837 :
838 :
839 0 : if (ipa > -1) {
840 0 : TParticle * mother = (TParticle *) fParticles.At(ipa);
841 0 : kfMo = TMath::Abs(mother->GetPdgCode());
842 0 : }
843 :
844 : // What to keep in Stack?
845 : Bool_t flavorOK = kFALSE;
846 : Bool_t selectOK = kFALSE;
847 0 : if (fFeedDownOpt) {
848 0 : if (kfl >= fFlavorSelect) flavorOK = kTRUE;
849 : } else {
850 0 : if (kfl > fFlavorSelect) {
851 : nc = -1;
852 0 : break;
853 : }
854 0 : if (kfl == fFlavorSelect) flavorOK = kTRUE;
855 : }
856 0 : switch (fStackFillOpt) {
857 : case kHeavyFlavor:
858 : case kFlavorSelection:
859 : selectOK = kTRUE;
860 0 : break;
861 : case kParentSelection:
862 0 : if (ParentSelected(kf) || kf <= 10) selectOK = kTRUE;
863 : break;
864 : }
865 0 : if (flavorOK && selectOK) {
866 : //
867 : // Heavy flavor hadron or quark
868 : //
869 : // Kinematic seletion on final state heavy flavor mesons
870 0 : if (ParentSelected(kf) && !KinematicSelection(iparticle, 0))
871 : {
872 0 : continue;
873 : }
874 0 : pSelected[i] = 1;
875 0 : if (ParentSelected(kf)) ++nParents; // Update parent count
876 : // printf("\n particle (HF) %d %d %d", i, pSelected[i], kf);
877 : } else {
878 : // Kinematic seletion on decay products
879 0 : if (fCutOnChild && ParentSelected(kfMo) && ChildSelected(kf)
880 0 : && !KinematicSelection(iparticle, 1))
881 : {
882 0 : continue;
883 : }
884 : //
885 : // Decay products
886 : // Select if mother was selected and is not tracked
887 :
888 0 : if (pSelected[ipa] &&
889 0 : !trackIt[ipa] && // mother will be tracked ?
890 0 : kfMo != 5 && // mother is b-quark, don't store fragments
891 0 : kfMo != 4 && // mother is c-quark, don't store fragments
892 0 : kf != 92) // don't store string
893 : {
894 : //
895 : // Semi-stable or de-selected: diselect decay products:
896 : //
897 : //
898 0 : if (pSelected[i] == -1 || fDecayer->GetLifetime(kf) > fMaxLifeTime)
899 : {
900 0 : Int_t ipF = iparticle->GetFirstDaughter();
901 0 : Int_t ipL = iparticle->GetLastDaughter();
902 0 : if (ipF > 0) for (j = ipF-1; j < ipL; j++) pSelected[j] = -1;
903 0 : }
904 : // printf("\n particle (decay) %d %d %d", i, pSelected[i], kf);
905 0 : pSelected[i] = (pSelected[i] == -1) ? 0 : 1;
906 0 : }
907 : }
908 0 : if (pSelected[i] == -1) pSelected[i] = 0;
909 0 : if (!pSelected[i]) continue;
910 : // Count quarks only if you did not include fragmentation
911 0 : if (fFragmentation && kf <= 10) continue;
912 :
913 0 : nc++;
914 : // Decision on tracking
915 0 : trackIt[i] = 0;
916 : //
917 : // Track final state particle
918 0 : if (ks == 1) trackIt[i] = 1;
919 : // Track semi-stable particles
920 : // if ((ks == 1) || (fDecayer->GetLifetime(kf) > fMaxLifeTime)) trackIt[i] = 1;
921 : // Track particles selected by process if undecayed.
922 0 : if (fForceDecay == kNoDecay) {
923 0 : if (ParentSelected(kf)) trackIt[i] = 1;
924 : } else {
925 0 : if (ParentSelected(kf)) trackIt[i] = 0;
926 : }
927 0 : if (trackIt[i] == 1) ++nTkbles; // Update trackable counter
928 : //
929 : //
930 :
931 0 : } // particle selection loop
932 0 : if (nc > 0) {
933 0 : for (i = 0; i<np; i++) {
934 0 : if (!pSelected[i]) continue;
935 0 : TParticle * iparticle = (TParticle *) fParticles.At(i);
936 0 : kf = CheckPDGCode(iparticle->GetPdgCode());
937 0 : Int_t ks = iparticle->GetStatusCode();
938 0 : p[0] = iparticle->Px();
939 0 : p[1] = iparticle->Py();
940 0 : p[2] = iparticle->Pz();
941 0 : p[3] = iparticle->Energy();
942 :
943 0 : origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm]
944 0 : origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm]
945 0 : origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm]
946 :
947 0 : Float_t tof = fTime + kconv*iparticle->T();
948 0 : Int_t ipa = iparticle->GetFirstMother()-1;
949 0 : Int_t iparent = (ipa > -1) ? pParent[ipa] : -1;
950 :
951 0 : PushTrack(fTrackIt*trackIt[i], iparent, kf,
952 : p[0], p[1], p[2], p[3],
953 0 : origin[0], origin[1], origin[2], tof,
954 : polar[0], polar[1], polar[2],
955 : kPPrimary, nt, 1., ks);
956 0 : pParent[i] = nt;
957 0 : KeepTrack(nt);
958 0 : fNprimaries++;
959 0 : } // PushTrack loop
960 : }
961 : } else {
962 0 : nc = GenerateMB();
963 : } // mb ?
964 :
965 0 : GetSubEventTime();
966 :
967 0 : delete[] pParent;
968 0 : delete[] pSelected;
969 0 : delete[] trackIt;
970 :
971 0 : if (nc > 0) {
972 0 : switch (fCountMode) {
973 : case kCountAll:
974 : // printf(" Count all \n");
975 0 : jev += nc;
976 0 : break;
977 : case kCountParents:
978 : // printf(" Count parents \n");
979 0 : jev += nParents;
980 0 : break;
981 : case kCountTrackables:
982 : // printf(" Count trackable \n");
983 0 : jev += nTkbles;
984 0 : break;
985 : }
986 0 : if (jev >= fNpart || fNpart == -1) {
987 0 : fKineBias=Float_t(fNpart)/Float_t(fTrials);
988 :
989 0 : fQ += fPythia->GetVINT(51);
990 0 : fX1 += fPythia->GetVINT(41);
991 0 : fX2 += fPythia->GetVINT(42);
992 0 : fTrialsRun += fTrials;
993 0 : fNev++;
994 0 : MakeHeader();
995 0 : break;
996 : }
997 : }
998 0 : } // event loop
999 0 : SetHighWaterMark(nt);
1000 : // adjust weight due to kinematic selection
1001 : // AdjustWeights();
1002 : // get cross-section
1003 0 : fXsection=fPythia->GetPARI(1);
1004 0 : }
1005 :
1006 : Int_t AliGenPythia::GenerateMB()
1007 : {
1008 : //
1009 : // Min Bias selection and other global selections
1010 : //
1011 0 : Int_t i, kf, nt, iparent;
1012 : Int_t nc = 0;
1013 : Double_t p[4];
1014 : Double_t polar[3] = {0,0,0};
1015 : Double_t origin[3] = {0,0,0};
1016 : // converts from mm/c to s
1017 : const Float_t kconv=0.001/2.999792458e8;
1018 :
1019 :
1020 0 : Int_t np = (fHadronisation) ? fParticles.GetEntriesFast() : fNpartons;
1021 :
1022 :
1023 :
1024 0 : Int_t* pParent = new Int_t[np];
1025 0 : for (i=0; i< np; i++) pParent[i] = -1;
1026 0 : if ((fProcess == kPyJets || fProcess == kPyDirectGamma || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi || fProcess == kPyJetsPWHG || fProcess == kPyCharmPWHG || fProcess == kPyBeautyPWHG)
1027 0 : && fEtMinJet > 0.) {
1028 0 : TParticle* jet1 = (TParticle *) fParticles.At(6);
1029 0 : TParticle* jet2 = (TParticle *) fParticles.At(7);
1030 :
1031 0 : if (!jet1 || ! jet2 || !CheckTrigger(jet1, jet2)) {
1032 0 : delete [] pParent;
1033 0 : return 0;
1034 : }
1035 0 : }
1036 :
1037 :
1038 : // Select events with fragmentation photon, decay photon, pi0, eta or other hadrons going to PHOS or EMCAL or central barrel,
1039 : // implemented primaryly for kPyJets, but extended to any kind of process.
1040 0 : if ((fFragPhotonInCalo || fPi0InCalo || fEtaInCalo || fEleInCalo || fHadronInCalo || fDecayPhotonInCalo) &&
1041 0 : (fCheckPHOS || fCheckEMCAL || fCheckBarrel) ) {
1042 0 : Bool_t ok = TriggerOnSelectedParticles(np);
1043 :
1044 0 : if(!ok) {
1045 0 : delete[] pParent;
1046 0 : return 0;
1047 : }
1048 0 : }
1049 :
1050 : // Check for diffraction
1051 0 : if(fkTuneForDiff) {
1052 0 : if(fItune==320 && ( (TMath::Abs(fEnergyCMS - 900) < 1) || (TMath::Abs(fEnergyCMS - 2760) < 1) || (TMath::Abs(fEnergyCMS - 7000) < 1)) ) {
1053 0 : if(!CheckDiffraction()) {
1054 0 : delete [] pParent;
1055 0 : return 0;
1056 : }
1057 : }
1058 : }
1059 :
1060 : // Check for minimum multiplicity
1061 0 : if (fTriggerMultiplicity > 0) {
1062 : Int_t multiplicity = 0;
1063 0 : for (i = 0; i < np; i++) {
1064 0 : TParticle * iparticle = (TParticle *) fParticles.At(i);
1065 :
1066 0 : Int_t statusCode = iparticle->GetStatusCode();
1067 :
1068 : // Initial state particle
1069 0 : if (statusCode != 1)
1070 0 : continue;
1071 : // eta cut
1072 0 : if (fTriggerMultiplicityEta > 0 && TMath::Abs(iparticle->Eta()) > fTriggerMultiplicityEta)
1073 0 : continue;
1074 : //multiplicity check for a given eta range
1075 0 : if ((fTriggerMultiplicityEtaMin != fTriggerMultiplicityEtaMax) &&
1076 0 : (iparticle->Eta() < fTriggerMultiplicityEtaMin || iparticle->Eta() > fTriggerMultiplicityEtaMax))
1077 0 : continue;
1078 : // pt cut
1079 0 : if (iparticle->Pt() < fTriggerMultiplicityPtMin)
1080 0 : continue;
1081 :
1082 0 : TParticlePDG* pdgPart = iparticle->GetPDG();
1083 0 : if (pdgPart && pdgPart->Charge() == 0)
1084 0 : continue;
1085 :
1086 0 : ++multiplicity;
1087 0 : }
1088 :
1089 0 : if (multiplicity < fTriggerMultiplicity) {
1090 0 : delete [] pParent;
1091 0 : return 0;
1092 : }
1093 0 : Printf("Triggered on event with multiplicity of %d >= %d", multiplicity, fTriggerMultiplicity);
1094 0 : }
1095 :
1096 :
1097 0 : if (fTriggerParticle) {
1098 : Bool_t triggered = kFALSE;
1099 0 : for (i = 0; i < np; i++) {
1100 0 : TParticle * iparticle = (TParticle *) fParticles.At(i);
1101 0 : kf = CheckPDGCode(iparticle->GetPdgCode());
1102 0 : if (kf != fTriggerParticle) continue;
1103 0 : if (iparticle->Pt() == 0.) continue;
1104 0 : if (TMath::Abs(iparticle->Y()) > fTriggerY) continue;
1105 0 : if (fTriggerEtaMin == fTriggerEta) {
1106 0 : if (TMath::Abs(iparticle->Eta()) > fTriggerEta) continue;
1107 : } else {
1108 0 : if (iparticle->Eta() > fTriggerEta || iparticle->Eta() < fTriggerEtaMin) continue;
1109 : }
1110 0 : if ( iparticle->Pt() > fTriggerMaxPt || iparticle->Pt() < fTriggerMinPt ) continue;
1111 : triggered = kTRUE;
1112 0 : break;
1113 : }
1114 0 : if (!triggered) {
1115 0 : delete [] pParent;
1116 0 : return 0;
1117 : }
1118 0 : }
1119 :
1120 :
1121 : // Check if there is a ccbar or bbbar pair with at least one of the two
1122 : // in fYMin < y < fYMax
1123 :
1124 0 : if (fProcess == kPyCharmppMNRwmi || fProcess == kPyBeautyppMNRwmi || fProcess == kPyBeautyJets) {
1125 : TParticle *partCheck;
1126 : TParticle *mother;
1127 : Bool_t theQ=kFALSE,theQbar=kFALSE,inYcut=kFALSE;
1128 : Bool_t theChild=kFALSE;
1129 : Bool_t triggered=kFALSE;
1130 : Float_t y;
1131 : Int_t pdg, mpdg, mpdgUpperFamily;
1132 0 : for(i = 0; i < np; i++) {
1133 0 : partCheck = (TParticle*)fParticles.At(i);
1134 0 : pdg = partCheck->GetPdgCode();
1135 0 : if(TMath::Abs(pdg) == fFlavorSelect) { // quark
1136 0 : if(pdg>0) { theQ=kTRUE; } else { theQbar=kTRUE; }
1137 0 : y = 0.5*TMath::Log((partCheck->Energy()+partCheck->Pz()+1.e-13)/
1138 0 : (partCheck->Energy()-partCheck->Pz()+1.e-13));
1139 0 : if(fUseYCutHQ && y>fYMinHQ && y<fYMaxHQ) inYcut=kTRUE;
1140 0 : if(!fUseYCutHQ && y>fYMin && y<fYMax) inYcut=kTRUE;
1141 : }
1142 0 : if(fTriggerParticle) {
1143 0 : if(TMath::Abs(pdg)==fTriggerParticle) triggered=kTRUE;
1144 : }
1145 0 : if(fCutOnChild && TMath::Abs(pdg) == fPdgCodeParticleforAcceptanceCut) {
1146 0 : Int_t mi = partCheck->GetFirstMother() - 1;
1147 0 : if(mi<0) continue;
1148 0 : mother = (TParticle*)fParticles.At(mi);
1149 0 : mpdg=TMath::Abs(mother->GetPdgCode());
1150 0 : mpdgUpperFamily=(mpdg>1000 ? mpdg+1000 : mpdg+100); // keep e from c from b
1151 0 : if ( ParentSelected(mpdg) ||
1152 0 : (fFlavorSelect==5 && ParentSelected(mpdgUpperFamily))) {
1153 0 : if (KinematicSelection(partCheck,1)) {
1154 : theChild=kTRUE;
1155 0 : }
1156 : }
1157 0 : }
1158 : }
1159 0 : if (!theQ || !theQbar || !inYcut) { // one of the c/b conditions not satisfied
1160 0 : delete[] pParent;
1161 0 : return 0;
1162 : }
1163 0 : if (fCutOnChild && !theChild) { // one of the child conditions not satisfied
1164 0 : delete[] pParent;
1165 0 : return 0;
1166 : }
1167 0 : if(fTriggerParticle && !triggered) { // particle requested is not produced
1168 0 : delete[] pParent;
1169 0 : return 0;
1170 : }
1171 :
1172 0 : }
1173 :
1174 : //Introducing child cuts in case kPyW, kPyZ, kPyMb, and kPyMbNonDiff
1175 0 : if ( (
1176 0 : fProcess == kPyWPWHG ||
1177 0 : fProcess == kPyW ||
1178 0 : fProcess == kPyZ ||
1179 0 : fProcess == kPyZgamma ||
1180 0 : fProcess == kPyMbDefault ||
1181 0 : fProcess == kPyMb ||
1182 0 : fProcess == kPyMbAtlasTuneMC09 ||
1183 0 : fProcess == kPyMbWithDirectPhoton ||
1184 0 : fProcess == kPyMbNonDiffr)
1185 0 : && (fCutOnChild == 1) ) {
1186 0 : if ( !CheckKinematicsOnChild() ) {
1187 0 : delete[] pParent;
1188 0 : return 0;
1189 : }
1190 : }
1191 :
1192 :
1193 0 : for (i = 0; i < np; i++) {
1194 : Int_t trackIt = 0;
1195 0 : TParticle * iparticle = (TParticle *) fParticles.At(i);
1196 0 : kf = CheckPDGCode(iparticle->GetPdgCode());
1197 0 : Int_t ks = iparticle->GetStatusCode();
1198 0 : Int_t km = iparticle->GetFirstMother();
1199 : if (
1200 0 : (((ks == 1 && kf!=0 && KinematicSelection(iparticle, 0)) || (ks !=1)) && ((fStackFillOpt != kHeavyFlavor) || IsFromHeavyFlavor(i))) ||
1201 0 : ((fProcess == kPyJets || fProcess == kPyBeautyJets || fProcess == kPyJetsPWHG || fProcess == kPyCharmPWHG || fProcess == kPyBeautyPWHG) && ks == 21 && km == 0 && i>1)
1202 : )
1203 : {
1204 0 : nc++;
1205 0 : if (ks == 1) trackIt = 1;
1206 0 : Int_t ipa = iparticle->GetFirstMother()-1;
1207 :
1208 0 : iparent = (ipa > -1) ? pParent[ipa] : -1;
1209 :
1210 : //
1211 : // store track information
1212 0 : p[0] = iparticle->Px();
1213 0 : p[1] = iparticle->Py();
1214 0 : p[2] = iparticle->Pz();
1215 0 : p[3] = iparticle->Energy();
1216 :
1217 :
1218 0 : origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm]
1219 0 : origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm]
1220 0 : origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm]
1221 :
1222 0 : Float_t tof = fTime + fEventTime + kconv * iparticle->T();
1223 :
1224 0 : PushTrack(fTrackIt*trackIt, iparent, kf,
1225 : p[0], p[1], p[2], p[3],
1226 0 : origin[0], origin[1], origin[2], tof,
1227 : polar[0], polar[1], polar[2],
1228 : kPPrimary, nt, 1., ks);
1229 0 : fNprimaries++;
1230 0 : KeepTrack(nt);
1231 0 : pParent[i] = nt;
1232 0 : SetHighWaterMark(nt);
1233 :
1234 0 : } // select particle
1235 : } // particle loop
1236 :
1237 0 : delete[] pParent;
1238 :
1239 0 : return 1;
1240 0 : }
1241 :
1242 :
1243 : void AliGenPythia::FinishRun()
1244 : {
1245 : // Print x-section summary
1246 0 : fPythia->Pystat(1);
1247 :
1248 0 : if (fNev > 0.) {
1249 0 : fQ /= fNev;
1250 0 : fX1 /= fNev;
1251 0 : fX2 /= fNev;
1252 0 : }
1253 :
1254 0 : printf("\nTotal number of Pyevnt() calls %d\n", fTrialsRun);
1255 0 : printf("\nMean Q, x1, x2: %f %f %f\n", fQ, fX1, fX2);
1256 0 : }
1257 :
1258 : void AliGenPythia::AdjustWeights() const
1259 : {
1260 : // Adjust the weights after generation of all events
1261 : //
1262 0 : if (gAlice) {
1263 : TParticle *part;
1264 0 : Int_t ntrack=gAlice->GetMCApp()->GetNtrack();
1265 0 : for (Int_t i=0; i<ntrack; i++) {
1266 0 : part= gAlice->GetMCApp()->Particle(i);
1267 0 : part->SetWeight(part->GetWeight()*fKineBias);
1268 : }
1269 0 : }
1270 0 : }
1271 :
1272 : void AliGenPythia::SetNuclei(Int_t a1, Int_t a2, Int_t pdfset)
1273 : {
1274 : // Treat protons as inside nuclei with mass numbers a1 and a2
1275 :
1276 0 : fAProjectile = a1;
1277 0 : fATarget = a2;
1278 0 : fNucPdf = pdfset; // 0 EKS98 9 EPS09LO 19 EPS09NLO
1279 0 : fUseNuclearPDF = kTRUE;
1280 0 : fSetNuclei = kTRUE;
1281 0 : }
1282 :
1283 :
1284 : void AliGenPythia::MakeHeader()
1285 : {
1286 : //
1287 : // Make header for the simulated event
1288 : //
1289 0 : if (gAlice) {
1290 0 : if (gAlice->GetEvNumber()>=fDebugEventFirst &&
1291 0 : gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(2);
1292 : }
1293 :
1294 : // Builds the event header, to be called after each event
1295 0 : if (fHeader) delete fHeader;
1296 0 : fHeader = new AliGenPythiaEventHeader("Pythia");
1297 0 : fHeader->SetTitle(GetTitle());
1298 : //
1299 : // Event type
1300 0 : if(fProcDiff>0){
1301 : // if(fProcDiff == 92 || fProcDiff == 93) printf("\n\n\n\n\n");
1302 : // printf("fPythia->GetMSTI(1) = %d fProcDiff = %d\n",fPythia->GetMSTI(1), fProcDiff);
1303 0 : ((AliGenPythiaEventHeader*) fHeader)->SetProcessType(fProcDiff);
1304 0 : }
1305 : else
1306 0 : ((AliGenPythiaEventHeader*) fHeader)->SetProcessType(fPythia->GetMSTI(1));
1307 : //
1308 : // Number of trials
1309 0 : ((AliGenPythiaEventHeader*) fHeader)->SetTrials(fTrials);
1310 : //
1311 : // Event Vertex
1312 0 : fHeader->SetPrimaryVertex(fVertex);
1313 0 : fHeader->SetInteractionTime(fTime+fEventTime);
1314 : //
1315 : // Number of primaries
1316 0 : fHeader->SetNProduced(fNprimaries);
1317 : //
1318 : // Event weight
1319 0 : fHeader->SetEventWeight(fPythia->GetVINT(97));
1320 : //
1321 : // Number of MPI
1322 0 : fHeader->SetNMPI(fPythia->GetNMPI());
1323 : //
1324 : // Jets that have triggered
1325 :
1326 : //Need to store jets for b-jet studies too!
1327 0 : if (fProcess == kPyJets || fProcess == kPyDirectGamma || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi || fProcess == kPyJetsPWHG || fProcess == kPyCharmPWHG || fProcess == kPyBeautyPWHG)
1328 : {
1329 0 : Int_t ntrig, njet;
1330 0 : Float_t jets[4][10];
1331 0 : GetJets(njet, ntrig, jets);
1332 :
1333 :
1334 0 : for (Int_t i = 0; i < ntrig; i++) {
1335 0 : ((AliGenPythiaEventHeader*) fHeader)->AddJet(jets[0][i], jets[1][i], jets[2][i],
1336 0 : jets[3][i]);
1337 : }
1338 0 : }
1339 : //
1340 : // Copy relevant information from external header, if present.
1341 : //
1342 0 : Float_t uqJet[4];
1343 :
1344 0 : if (fRL) {
1345 0 : AliGenPythiaEventHeader* exHeader = (AliGenPythiaEventHeader*) (fRL->GetHeader()->GenEventHeader());
1346 0 : for (Int_t i = 0; i < exHeader->NTriggerJets(); i++)
1347 : {
1348 0 : printf("Adding Jet %d %d \n", i, exHeader->NTriggerJets());
1349 :
1350 :
1351 0 : exHeader->TriggerJet(i, uqJet);
1352 0 : ((AliGenPythiaEventHeader*) fHeader)->AddUQJet(uqJet[0], uqJet[1], uqJet[2], uqJet[3]);
1353 : }
1354 0 : }
1355 : //
1356 : // Store quenching parameters
1357 : //
1358 0 : if (fQuench){
1359 0 : Double_t z[4] = {0.};
1360 0 : Double_t xp = 0.;
1361 0 : Double_t yp = 0.;
1362 :
1363 0 : if (fQuench == 1) {
1364 : // Pythia::Quench()
1365 0 : fPythia->GetQuenchingParameters(xp, yp, z);
1366 0 : } else if (fQuench == 2){
1367 : // Pyquen
1368 0 : Double_t r1 = PARIMP.rb1;
1369 0 : Double_t r2 = PARIMP.rb2;
1370 0 : Double_t b = PARIMP.b1;
1371 0 : Double_t r = 0.5 * TMath::Sqrt(2. * (r1 * r1 + r2 * r2) - b * b);
1372 0 : Double_t phi = PARIMP.psib1;
1373 0 : xp = r * TMath::Cos(phi);
1374 0 : yp = r * TMath::Sin(phi);
1375 :
1376 0 : } else if (fQuench == 4) {
1377 : // QPythia
1378 0 : Double_t xy[2];
1379 0 : Double_t i0i1[2];
1380 0 : AliFastGlauber::Instance()->GetSavedXY(xy);
1381 0 : AliFastGlauber::Instance()->GetSavedI0I1(i0i1);
1382 0 : xp = xy[0];
1383 0 : yp = xy[1];
1384 0 : ((AliGenPythiaEventHeader*) fHeader)->SetImpactParameter(fImpact);
1385 0 : }
1386 :
1387 0 : ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp);
1388 0 : ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z);
1389 0 : }
1390 : //
1391 : // Store pt^hard and cross-section
1392 0 : ((AliGenPythiaEventHeader*) fHeader)->SetPtHard(fPythia->GetVINT(47));
1393 0 : ((AliGenPythiaEventHeader*) fHeader)->SetXsection(fPythia->GetPARI(1));
1394 :
1395 : //
1396 : // Store Event Weight
1397 0 : ((AliGenPythiaEventHeader*) fHeader)->SetEventWeight(fPythia->GetPARI(7)*fPythia->GetPARI(10));
1398 : // PARI(7) is 1 or -1, for weighted generation with accept/reject, e.g. POWHEG
1399 : // PARI(10) is a weight associated with reweighted generation, using Pyevwt
1400 : //
1401 : // Pass header
1402 : //
1403 0 : AddHeader(fHeader);
1404 0 : fHeader = 0x0;
1405 0 : }
1406 :
1407 : Bool_t AliGenPythia::CheckTrigger(const TParticle* jet1, const TParticle* jet2)
1408 : {
1409 : // Check the kinematic trigger condition
1410 : //
1411 0 : Double_t eta[2];
1412 0 : eta[0] = jet1->Eta();
1413 0 : eta[1] = jet2->Eta();
1414 0 : Double_t phi[2];
1415 0 : phi[0] = jet1->Phi();
1416 0 : phi[1] = jet2->Phi();
1417 : Int_t pdg[2];
1418 0 : pdg[0] = jet1->GetPdgCode();
1419 0 : pdg[1] = jet2->GetPdgCode();
1420 : Bool_t triggered = kFALSE;
1421 :
1422 0 : if (fProcess == kPyJets || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi || fProcess == kPyJetsPWHG || fProcess == kPyCharmPWHG || fProcess == kPyBeautyPWHG) {
1423 0 : Int_t njets = 0;
1424 0 : Int_t ntrig = 0;
1425 0 : Float_t jets[4][10];
1426 : //
1427 : // Use Pythia clustering on parton level to determine jet axis
1428 : //
1429 0 : GetJets(njets, ntrig, jets);
1430 :
1431 0 : if (ntrig || fEtMinJet == 0.) triggered = kTRUE;
1432 : //
1433 0 : } else {
1434 : Int_t ij = 0;
1435 : Int_t ig = 1;
1436 0 : if (pdg[0] == kGamma) {
1437 : ij = 1;
1438 : ig = 0;
1439 0 : }
1440 : //Check eta range first...
1441 0 : if ((eta[ij] < fEtaMaxJet && eta[ij] > fEtaMinJet) &&
1442 0 : (eta[ig] < fEtaMaxGamma && eta[ig] > fEtaMinGamma))
1443 : {
1444 : //Eta is okay, now check phi range
1445 0 : if ((phi[ij] < fPhiMaxJet && phi[ij] > fPhiMinJet) &&
1446 0 : (phi[ig] < fPhiMaxGamma && phi[ig] > fPhiMinGamma))
1447 : {
1448 : triggered = kTRUE;
1449 0 : }
1450 : }
1451 : }
1452 0 : return triggered;
1453 0 : }
1454 :
1455 :
1456 :
1457 : Bool_t AliGenPythia::CheckKinematicsOnChild(){
1458 : //
1459 : //Checking Kinematics on Child (status code 1, particle code ?, kin cuts
1460 : //
1461 : Bool_t checking = kFALSE;
1462 : Int_t j, kcode, ks, km;
1463 : Int_t nPartAcc = 0; //number of particles in the acceptance range
1464 : Int_t numberOfAcceptedParticles = 1;
1465 0 : if (fNumberOfAcceptedParticles != 0) { numberOfAcceptedParticles = fNumberOfAcceptedParticles; }
1466 0 : Int_t npart = fParticles.GetEntriesFast();
1467 :
1468 0 : for (j = 0; j<npart; j++) {
1469 0 : TParticle * jparticle = (TParticle *) fParticles.At(j);
1470 0 : kcode = TMath::Abs( CheckPDGCode(jparticle->GetPdgCode()) );
1471 0 : ks = jparticle->GetStatusCode();
1472 0 : km = jparticle->GetFirstMother();
1473 :
1474 0 : if( (ks == 1) && (kcode == fPdgCodeParticleforAcceptanceCut) && (KinematicSelection(jparticle,1)) ){
1475 0 : nPartAcc++;
1476 0 : }
1477 0 : if( numberOfAcceptedParticles <= nPartAcc){
1478 : checking = kTRUE;
1479 0 : break;
1480 : }
1481 0 : }
1482 :
1483 0 : return checking;
1484 : }
1485 :
1486 : void AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
1487 : {
1488 : //
1489 : // Load event into Pythia Common Block
1490 : //
1491 :
1492 0 : Int_t npart = stack -> GetNprimary();
1493 : Int_t n0 = 0;
1494 :
1495 0 : if (!flag) {
1496 0 : (fPythia->GetPyjets())->N = npart;
1497 0 : } else {
1498 0 : n0 = (fPythia->GetPyjets())->N;
1499 0 : (fPythia->GetPyjets())->N = n0 + npart;
1500 : }
1501 :
1502 :
1503 0 : for (Int_t part = 0; part < npart; part++) {
1504 0 : TParticle *mPart = stack->Particle(part);
1505 :
1506 0 : Int_t kf = mPart->GetPdgCode();
1507 0 : Int_t ks = mPart->GetStatusCode();
1508 0 : Int_t idf = mPart->GetFirstDaughter();
1509 0 : Int_t idl = mPart->GetLastDaughter();
1510 :
1511 0 : if (reHadr) {
1512 0 : if (ks == 11 || ks == 12) {
1513 0 : ks -= 10;
1514 : idf = -1;
1515 : idl = -1;
1516 0 : }
1517 : }
1518 :
1519 0 : Float_t px = mPart->Px();
1520 0 : Float_t py = mPart->Py();
1521 0 : Float_t pz = mPart->Pz();
1522 0 : Float_t e = mPart->Energy();
1523 0 : Float_t m = mPart->GetCalcMass();
1524 :
1525 :
1526 0 : (fPythia->GetPyjets())->P[0][part+n0] = px;
1527 0 : (fPythia->GetPyjets())->P[1][part+n0] = py;
1528 0 : (fPythia->GetPyjets())->P[2][part+n0] = pz;
1529 0 : (fPythia->GetPyjets())->P[3][part+n0] = e;
1530 0 : (fPythia->GetPyjets())->P[4][part+n0] = m;
1531 :
1532 0 : (fPythia->GetPyjets())->K[1][part+n0] = kf;
1533 0 : (fPythia->GetPyjets())->K[0][part+n0] = ks;
1534 0 : (fPythia->GetPyjets())->K[3][part+n0] = idf + 1;
1535 0 : (fPythia->GetPyjets())->K[4][part+n0] = idl + 1;
1536 0 : (fPythia->GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1;
1537 : }
1538 0 : }
1539 :
1540 : void AliGenPythia::LoadEvent(const TObjArray* stack, Int_t flag, Int_t reHadr)
1541 : {
1542 : //
1543 : // Load event into Pythia Common Block
1544 : //
1545 :
1546 0 : Int_t npart = stack -> GetEntries();
1547 : Int_t n0 = 0;
1548 :
1549 0 : if (!flag) {
1550 0 : (fPythia->GetPyjets())->N = npart;
1551 0 : } else {
1552 0 : n0 = (fPythia->GetPyjets())->N;
1553 0 : (fPythia->GetPyjets())->N = n0 + npart;
1554 : }
1555 :
1556 :
1557 0 : for (Int_t part = 0; part < npart; part++) {
1558 0 : TParticle *mPart = dynamic_cast<TParticle *>(stack->At(part));
1559 0 : if (!mPart) continue;
1560 :
1561 0 : Int_t kf = mPart->GetPdgCode();
1562 0 : Int_t ks = mPart->GetStatusCode();
1563 0 : Int_t idf = mPart->GetFirstDaughter();
1564 0 : Int_t idl = mPart->GetLastDaughter();
1565 :
1566 0 : if (reHadr) {
1567 0 : if (ks == 11 || ks == 12) {
1568 0 : ks -= 10;
1569 : idf = -1;
1570 : idl = -1;
1571 0 : }
1572 : }
1573 :
1574 0 : Float_t px = mPart->Px();
1575 0 : Float_t py = mPart->Py();
1576 0 : Float_t pz = mPart->Pz();
1577 0 : Float_t e = mPart->Energy();
1578 0 : Float_t m = mPart->GetCalcMass();
1579 :
1580 :
1581 0 : (fPythia->GetPyjets())->P[0][part+n0] = px;
1582 0 : (fPythia->GetPyjets())->P[1][part+n0] = py;
1583 0 : (fPythia->GetPyjets())->P[2][part+n0] = pz;
1584 0 : (fPythia->GetPyjets())->P[3][part+n0] = e;
1585 0 : (fPythia->GetPyjets())->P[4][part+n0] = m;
1586 :
1587 0 : (fPythia->GetPyjets())->K[1][part+n0] = kf;
1588 0 : (fPythia->GetPyjets())->K[0][part+n0] = ks;
1589 0 : (fPythia->GetPyjets())->K[3][part+n0] = idf + 1;
1590 0 : (fPythia->GetPyjets())->K[4][part+n0] = idl + 1;
1591 0 : (fPythia->GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1;
1592 0 : }
1593 0 : }
1594 :
1595 :
1596 : void AliGenPythia::RecJetsUA1(Int_t& njets, Float_t jets [4][50])
1597 : {
1598 : //
1599 : // Calls the Pythia jet finding algorithm to find jets in the current event
1600 : //
1601 : //
1602 : //
1603 : // Save jets
1604 0 : Int_t n = fPythia->GetN();
1605 :
1606 : //
1607 : // Run Jet Finder
1608 0 : fPythia->Pycell(njets);
1609 : Int_t i;
1610 0 : for (i = 0; i < njets; i++) {
1611 0 : Float_t px = (fPythia->GetPyjets())->P[0][n+i];
1612 0 : Float_t py = (fPythia->GetPyjets())->P[1][n+i];
1613 0 : Float_t pz = (fPythia->GetPyjets())->P[2][n+i];
1614 0 : Float_t e = (fPythia->GetPyjets())->P[3][n+i];
1615 :
1616 0 : jets[0][i] = px;
1617 0 : jets[1][i] = py;
1618 0 : jets[2][i] = pz;
1619 0 : jets[3][i] = e;
1620 : }
1621 0 : }
1622 :
1623 :
1624 :
1625 : void AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10])
1626 : {
1627 : //
1628 : // Calls the Pythia clustering algorithm to find jets in the current event
1629 : //
1630 0 : Int_t n = fPythia->GetN();
1631 0 : nJets = 0;
1632 0 : nJetsTrig = 0;
1633 0 : if (fJetReconstruction == kCluster) {
1634 : //
1635 : // Configure cluster algorithm
1636 : //
1637 0 : fPythia->SetPARU(43, 2.);
1638 0 : fPythia->SetMSTU(41, 1);
1639 : //
1640 : // Call cluster algorithm
1641 : //
1642 0 : fPythia->Pyclus(nJets);
1643 : //
1644 : // Loading jets from common block
1645 : //
1646 0 : } else {
1647 :
1648 : //
1649 : // Run Jet Finder
1650 0 : fPythia->Pycell(nJets);
1651 : }
1652 :
1653 : Int_t i;
1654 0 : for (i = 0; i < nJets; i++) {
1655 0 : Float_t px = (fPythia->GetPyjets())->P[0][n+i];
1656 0 : Float_t py = (fPythia->GetPyjets())->P[1][n+i];
1657 0 : Float_t pz = (fPythia->GetPyjets())->P[2][n+i];
1658 0 : Float_t e = (fPythia->GetPyjets())->P[3][n+i];
1659 0 : Float_t pt = TMath::Sqrt(px * px + py * py);
1660 0 : Float_t phi = TMath::Pi() + TMath::ATan2(-py, -px);
1661 0 : Float_t theta = TMath::ATan2(pt,pz);
1662 0 : Float_t et = e * TMath::Sin(theta);
1663 0 : Float_t eta = -TMath::Log(TMath::Tan(theta / 2.));
1664 : if (
1665 0 : eta > fEtaMinJet && eta < fEtaMaxJet &&
1666 0 : phi > fPhiMinJet && phi < fPhiMaxJet &&
1667 0 : et > fEtMinJet && et < fEtMaxJet
1668 : )
1669 : {
1670 0 : jets[0][nJetsTrig] = px;
1671 0 : jets[1][nJetsTrig] = py;
1672 0 : jets[2][nJetsTrig] = pz;
1673 0 : jets[3][nJetsTrig] = e;
1674 0 : nJetsTrig++;
1675 : // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg);
1676 0 : } else {
1677 : // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg);
1678 : }
1679 : }
1680 0 : }
1681 :
1682 : void AliGenPythia::GetSubEventTime()
1683 : {
1684 : // Calculates time of the next subevent
1685 0 : fEventTime = 0.;
1686 0 : if (fEventsTime) {
1687 : TArrayF &array = *fEventsTime;
1688 0 : fEventTime = array[fCurSubEvent++];
1689 0 : }
1690 : // printf(" Event time: %d %f %p",fCurSubEvent,fEventTime,fEventsTime);
1691 0 : return;
1692 : }
1693 :
1694 : Bool_t AliGenPythia::IsInBarrel(Float_t eta) const
1695 : {
1696 : // Is particle in Central Barrel acceptance?
1697 : // etamin=-etamax
1698 0 : if( eta < fTriggerEta )
1699 0 : return kTRUE;
1700 : else
1701 0 : return kFALSE;
1702 0 : }
1703 :
1704 : Bool_t AliGenPythia::IsInEMCAL(Float_t phi, Float_t eta) const
1705 : {
1706 : // Is particle in EMCAL acceptance?
1707 : // phi in degrees, etamin=-etamax
1708 0 : if(phi > fEMCALMinPhi && phi < fEMCALMaxPhi &&
1709 0 : eta < fEMCALEta )
1710 0 : return kTRUE;
1711 : else
1712 0 : return kFALSE;
1713 0 : }
1714 :
1715 : Bool_t AliGenPythia::IsInPHOS(Float_t phi, Float_t eta, Int_t iparticle)
1716 : {
1717 : // Is particle in PHOS acceptance?
1718 : // Acceptance slightly larger considered.
1719 : // phi in degrees, etamin=-etamax
1720 : // iparticle is the index of the particle to be checked, for PHOS rotation case
1721 :
1722 0 : if(phi > fPHOSMinPhi && phi < fPHOSMaxPhi &&
1723 0 : eta < fPHOSEta )
1724 0 : return kTRUE;
1725 : else
1726 : {
1727 0 : if( fCheckPHOSeta && eta < fPHOSEta) fPHOSRotateCandidate = iparticle;
1728 :
1729 0 : return kFALSE;
1730 : }
1731 0 : }
1732 :
1733 : void AliGenPythia::RotatePhi(Bool_t& okdd)
1734 : {
1735 : //Rotate event in phi to enhance events in PHOS acceptance
1736 :
1737 0 : if(fPHOSRotateCandidate < 0) return ;
1738 :
1739 : //calculate the new position random between fPHOSMinPhi and fPHOSMaxPhi
1740 0 : Double_t phiPHOSmin = TMath::Pi()*fPHOSMinPhi/180;
1741 0 : Double_t phiPHOSmax = TMath::Pi()*fPHOSMaxPhi/180;
1742 0 : Double_t phiPHOS = (AliPythiaRndm::GetPythiaRandom())->Uniform(phiPHOSmin,phiPHOSmax);
1743 :
1744 : //calculate deltaphi
1745 0 : TParticle* ph = (TParticle *) fParticles.At(fPHOSRotateCandidate);
1746 0 : Double_t phphi = ph->Phi();
1747 0 : Double_t deltaphi = phiPHOS - phphi;
1748 :
1749 :
1750 :
1751 : //loop for all particles and produce the phi rotation
1752 0 : Int_t np = (fHadronisation) ? fParticles.GetEntriesFast() : fNpartons;
1753 : Double_t oldphi, newphi;
1754 : Double_t newVx, newVy, r, vZ, time;
1755 : Double_t newPx, newPy, pt, pz, e;
1756 0 : for(Int_t i=0; i< np; i++) {
1757 0 : TParticle* iparticle = (TParticle *) fParticles.At(i);
1758 0 : oldphi = iparticle->Phi();
1759 0 : newphi = oldphi + deltaphi;
1760 0 : if(newphi < 0) newphi = 2*TMath::Pi() + newphi; // correct angle
1761 0 : if(newphi > 2*TMath::Pi()) newphi = newphi - 2*TMath::Pi(); // correct angle
1762 :
1763 0 : r = iparticle->R();
1764 0 : newVx = r * TMath::Cos(newphi);
1765 0 : newVy = r * TMath::Sin(newphi);
1766 0 : vZ = iparticle->Vz(); // don't transform
1767 0 : time = iparticle->T(); // don't transform
1768 :
1769 0 : pt = iparticle->Pt();
1770 0 : newPx = pt * TMath::Cos(newphi);
1771 0 : newPy = pt * TMath::Sin(newphi);
1772 0 : pz = iparticle->Pz(); // don't transform
1773 0 : e = iparticle->Energy(); // don't transform
1774 :
1775 : // apply rotation
1776 0 : iparticle->SetProductionVertex(newVx, newVy, vZ, time);
1777 0 : iparticle->SetMomentum(newPx, newPy, pz, e);
1778 :
1779 : } //end particle loop
1780 :
1781 : // now let's check that we put correctly the candidate photon in PHOS
1782 0 : Float_t phi = ph->Phi()*180./TMath::Pi(); //Convert to degrees
1783 0 : Float_t eta =TMath::Abs(ph->Eta());//in calos etamin=-etamax
1784 0 : if(IsInPHOS(phi,eta,-1))
1785 0 : okdd = kTRUE;
1786 :
1787 : // reset the value for next event
1788 0 : fPHOSRotateCandidate = -1;
1789 :
1790 0 : }
1791 :
1792 :
1793 : Bool_t AliGenPythia::CheckDiffraction()
1794 : {
1795 : // use this method only with Perugia-0 tune!
1796 :
1797 : // printf("AAA\n");
1798 :
1799 0 : Int_t np = (fHadronisation) ? fParticles.GetEntriesFast() : fNpartons;
1800 :
1801 : Int_t iPart1=-1;
1802 : Int_t iPart2=-1;
1803 :
1804 : Double_t y1 = 1e10;
1805 : Double_t y2 = -1e10;
1806 :
1807 : const Int_t kNstable=20;
1808 : const Int_t pdgStable[20] = {
1809 : 22, // Photon
1810 : 11, // Electron
1811 : 12, // Electron Neutrino
1812 : 13, // Muon
1813 : 14, // Muon Neutrino
1814 : 15, // Tau
1815 : 16, // Tau Neutrino
1816 : 211, // Pion
1817 : 321, // Kaon
1818 : 311, // K0
1819 : 130, // K0s
1820 : 310, // K0l
1821 : 2212, // Proton
1822 : 2112, // Neutron
1823 : 3122, // Lambda_0
1824 : 3112, // Sigma Minus
1825 : 3222, // Sigma Plus
1826 : 3312, // Xsi Minus
1827 : 3322, // Xsi0
1828 : 3334 // Omega
1829 : };
1830 :
1831 0 : for (Int_t i = 0; i < np; i++) {
1832 0 : TParticle * part = (TParticle *) fParticles.At(i);
1833 :
1834 0 : Int_t statusCode = part->GetStatusCode();
1835 :
1836 : // Initial state particle
1837 0 : if (statusCode != 1)
1838 0 : continue;
1839 :
1840 0 : Int_t pdg = TMath::Abs(part->GetPdgCode());
1841 : Bool_t isStable = kFALSE;
1842 0 : for (Int_t i1 = 0; i1 < kNstable; i1++) {
1843 0 : if (pdg == pdgStable[i1]) {
1844 : isStable = kTRUE;
1845 0 : break;
1846 : }
1847 : }
1848 0 : if(!isStable)
1849 0 : continue;
1850 :
1851 0 : Double_t y = part->Y();
1852 :
1853 0 : if (y < y1)
1854 : {
1855 : y1 = y;
1856 : iPart1 = i;
1857 0 : }
1858 0 : if (y > y2)
1859 : {
1860 : y2 = y;
1861 : iPart2 = i;
1862 0 : }
1863 0 : }
1864 :
1865 0 : if(iPart1<0 || iPart2<0) return kFALSE;
1866 :
1867 0 : y1=TMath::Abs(y1);
1868 0 : y2=TMath::Abs(y2);
1869 :
1870 0 : TParticle * part1 = (TParticle *) fParticles.At(iPart1);
1871 0 : TParticle * part2 = (TParticle *) fParticles.At(iPart2);
1872 :
1873 0 : Int_t pdg1 = part1->GetPdgCode();
1874 0 : Int_t pdg2 = part2->GetPdgCode();
1875 :
1876 :
1877 : Int_t iPart = -1;
1878 0 : if (pdg1 == 2212 && pdg2 == 2212)
1879 : {
1880 0 : if(y1 > y2)
1881 0 : iPart = iPart1;
1882 0 : else if(y1 < y2)
1883 0 : iPart = iPart2;
1884 : else {
1885 : iPart = iPart1;
1886 0 : if((AliPythiaRndm::GetPythiaRandom())->Uniform(0.,1.)>0.5) iPart = iPart2;
1887 : }
1888 : }
1889 0 : else if (pdg1 == 2212)
1890 0 : iPart = iPart1;
1891 0 : else if (pdg2 == 2212)
1892 0 : iPart = iPart2;
1893 :
1894 :
1895 :
1896 :
1897 :
1898 : Double_t M=-1.;
1899 0 : if(iPart>0) {
1900 0 : TParticle * part = (TParticle *) fParticles.At(iPart);
1901 0 : Double_t E= part->Energy();
1902 0 : Double_t P= part->P();
1903 0 : Double_t M2 = (fEnergyCMS-E-P)*(fEnergyCMS-E+P);
1904 0 : if(M2<0) return kFALSE;
1905 0 : M= TMath::Sqrt(M2);
1906 0 : }
1907 :
1908 0 : Double_t Mmin, Mmax, wSD, wDD, wND;
1909 0 : if(!GetWeightsDiffraction(M, Mmin, Mmax, wSD, wDD, wND)) return kFALSE;
1910 :
1911 0 : if(M>-1 && M<Mmin) return kFALSE;
1912 0 : if(M>Mmax) M=-1;
1913 :
1914 0 : Int_t procType=fPythia->GetMSTI(1);
1915 : Int_t proc0=2;
1916 0 : if(procType== 94) proc0=1;
1917 0 : if(procType== 92 || procType== 93) proc0=0;
1918 :
1919 : Int_t proc=2;
1920 0 : if(M>0) proc=0;
1921 0 : else if(proc0==1) proc=1;
1922 :
1923 0 : if(proc==1 && (AliPythiaRndm::GetPythiaRandom())->Uniform(0.,1.) > wDD) return kFALSE;
1924 0 : if(proc==2 && (AliPythiaRndm::GetPythiaRandom())->Uniform(0.,1.) > wND) return kFALSE;
1925 :
1926 :
1927 : // if(proc==1 || proc==2) return kFALSE;
1928 :
1929 0 : if(proc!=0) {
1930 0 : if(proc0!=0) fProcDiff = procType;
1931 0 : else fProcDiff = 95;
1932 0 : return kTRUE;
1933 : }
1934 :
1935 0 : if(wSD<0) AliError("wSD<0 ! \n");
1936 :
1937 0 : if((AliPythiaRndm::GetPythiaRandom())->Uniform(0.,1.)> wSD) return kFALSE;
1938 :
1939 : // printf("iPart = %d\n", iPart);
1940 :
1941 0 : if(iPart==iPart1) fProcDiff=93;
1942 0 : else if(iPart==iPart2) fProcDiff=92;
1943 : else {
1944 0 : printf("EROOR: iPart!=iPart1 && iPart!=iPart2\n");
1945 :
1946 : }
1947 :
1948 0 : return kTRUE;
1949 0 : }
1950 :
1951 :
1952 :
1953 : Bool_t AliGenPythia::GetWeightsDiffraction(Double_t M, Double_t &Mmin, Double_t &Mmax,
1954 : Double_t &wSD, Double_t &wDD, Double_t &wND)
1955 : {
1956 :
1957 : // 900 GeV
1958 0 : if(TMath::Abs(fEnergyCMS-900)<1 ){
1959 :
1960 : const Int_t nbin=400;
1961 0 : Double_t bin[]={
1962 : 1.080000, 1.577300, 2.074600, 2.571900, 3.069200, 3.566500,
1963 : 4.063800, 4.561100, 5.058400, 5.555700, 6.053000, 6.550300,
1964 : 7.047600, 7.544900, 8.042200, 8.539500, 9.036800, 9.534100,
1965 : 10.031400, 10.528700, 11.026000, 11.523300, 12.020600, 12.517900,
1966 : 13.015200, 13.512500, 14.009800, 14.507100, 15.004400, 15.501700,
1967 : 15.999000, 16.496300, 16.993600, 17.490900, 17.988200, 18.485500,
1968 : 18.982800, 19.480100, 19.977400, 20.474700, 20.972000, 21.469300,
1969 : 21.966600, 22.463900, 22.961200, 23.458500, 23.955800, 24.453100,
1970 : 24.950400, 25.447700, 25.945000, 26.442300, 26.939600, 27.436900,
1971 : 27.934200, 28.431500, 28.928800, 29.426100, 29.923400, 30.420700,
1972 : 30.918000, 31.415300, 31.912600, 32.409900, 32.907200, 33.404500,
1973 : 33.901800, 34.399100, 34.896400, 35.393700, 35.891000, 36.388300,
1974 : 36.885600, 37.382900, 37.880200, 38.377500, 38.874800, 39.372100,
1975 : 39.869400, 40.366700, 40.864000, 41.361300, 41.858600, 42.355900,
1976 : 42.853200, 43.350500, 43.847800, 44.345100, 44.842400, 45.339700,
1977 : 45.837000, 46.334300, 46.831600, 47.328900, 47.826200, 48.323500,
1978 : 48.820800, 49.318100, 49.815400, 50.312700, 50.810000, 51.307300,
1979 : 51.804600, 52.301900, 52.799200, 53.296500, 53.793800, 54.291100,
1980 : 54.788400, 55.285700, 55.783000, 56.280300, 56.777600, 57.274900,
1981 : 57.772200, 58.269500, 58.766800, 59.264100, 59.761400, 60.258700,
1982 : 60.756000, 61.253300, 61.750600, 62.247900, 62.745200, 63.242500,
1983 : 63.739800, 64.237100, 64.734400, 65.231700, 65.729000, 66.226300,
1984 : 66.723600, 67.220900, 67.718200, 68.215500, 68.712800, 69.210100,
1985 : 69.707400, 70.204700, 70.702000, 71.199300, 71.696600, 72.193900,
1986 : 72.691200, 73.188500, 73.685800, 74.183100, 74.680400, 75.177700,
1987 : 75.675000, 76.172300, 76.669600, 77.166900, 77.664200, 78.161500,
1988 : 78.658800, 79.156100, 79.653400, 80.150700, 80.648000, 81.145300,
1989 : 81.642600, 82.139900, 82.637200, 83.134500, 83.631800, 84.129100,
1990 : 84.626400, 85.123700, 85.621000, 86.118300, 86.615600, 87.112900,
1991 : 87.610200, 88.107500, 88.604800, 89.102100, 89.599400, 90.096700,
1992 : 90.594000, 91.091300, 91.588600, 92.085900, 92.583200, 93.080500,
1993 : 93.577800, 94.075100, 94.572400, 95.069700, 95.567000, 96.064300,
1994 : 96.561600, 97.058900, 97.556200, 98.053500, 98.550800, 99.048100,
1995 : 99.545400, 100.042700, 100.540000, 101.037300, 101.534600, 102.031900,
1996 : 102.529200, 103.026500, 103.523800, 104.021100, 104.518400, 105.015700,
1997 : 105.513000, 106.010300, 106.507600, 107.004900, 107.502200, 107.999500,
1998 : 108.496800, 108.994100, 109.491400, 109.988700, 110.486000, 110.983300,
1999 : 111.480600, 111.977900, 112.475200, 112.972500, 113.469800, 113.967100,
2000 : 114.464400, 114.961700, 115.459000, 115.956300, 116.453600, 116.950900,
2001 : 117.448200, 117.945500, 118.442800, 118.940100, 119.437400, 119.934700,
2002 : 120.432000, 120.929300, 121.426600, 121.923900, 122.421200, 122.918500,
2003 : 123.415800, 123.913100, 124.410400, 124.907700, 125.405000, 125.902300,
2004 : 126.399600, 126.896900, 127.394200, 127.891500, 128.388800, 128.886100,
2005 : 129.383400, 129.880700, 130.378000, 130.875300, 131.372600, 131.869900,
2006 : 132.367200, 132.864500, 133.361800, 133.859100, 134.356400, 134.853700,
2007 : 135.351000, 135.848300, 136.345600, 136.842900, 137.340200, 137.837500,
2008 : 138.334800, 138.832100, 139.329400, 139.826700, 140.324000, 140.821300,
2009 : 141.318600, 141.815900, 142.313200, 142.810500, 143.307800, 143.805100,
2010 : 144.302400, 144.799700, 145.297000, 145.794300, 146.291600, 146.788900,
2011 : 147.286200, 147.783500, 148.280800, 148.778100, 149.275400, 149.772700,
2012 : 150.270000, 150.767300, 151.264600, 151.761900, 152.259200, 152.756500,
2013 : 153.253800, 153.751100, 154.248400, 154.745700, 155.243000, 155.740300,
2014 : 156.237600, 156.734900, 157.232200, 157.729500, 158.226800, 158.724100,
2015 : 159.221400, 159.718700, 160.216000, 160.713300, 161.210600, 161.707900,
2016 : 162.205200, 162.702500, 163.199800, 163.697100, 164.194400, 164.691700,
2017 : 165.189000, 165.686300, 166.183600, 166.680900, 167.178200, 167.675500,
2018 : 168.172800, 168.670100, 169.167400, 169.664700, 170.162000, 170.659300,
2019 : 171.156600, 171.653900, 172.151200, 172.648500, 173.145800, 173.643100,
2020 : 174.140400, 174.637700, 175.135000, 175.632300, 176.129600, 176.626900,
2021 : 177.124200, 177.621500, 178.118800, 178.616100, 179.113400, 179.610700,
2022 : 180.108000, 180.605300, 181.102600, 181.599900, 182.097200, 182.594500,
2023 : 183.091800, 183.589100, 184.086400, 184.583700, 185.081000, 185.578300,
2024 : 186.075600, 186.572900, 187.070200, 187.567500, 188.064800, 188.562100,
2025 : 189.059400, 189.556700, 190.054000, 190.551300, 191.048600, 191.545900,
2026 : 192.043200, 192.540500, 193.037800, 193.535100, 194.032400, 194.529700,
2027 : 195.027000, 195.524300, 196.021600, 196.518900, 197.016200, 197.513500,
2028 : 198.010800, 198.508100, 199.005400, 199.502700, 200.000000};
2029 0 : Double_t w[]={
2030 : 1.000000, 0.643457, 0.645609, 0.648347, 0.604563, 0.605002,
2031 : 0.602819, 0.611473, 0.576412, 0.562354, 0.550216, 0.529285,
2032 : 0.534558, 0.534364, 0.530358, 0.518475, 0.489253, 0.469754,
2033 : 0.469825, 0.450513, 0.455849, 0.435312, 0.437210, 0.456686,
2034 : 0.413577, 0.427093, 0.426894, 0.418834, 0.409475, 0.388483,
2035 : 0.412208, 0.388912, 0.389611, 0.382572, 0.389220, 0.370964,
2036 : 0.380463, 0.370873, 0.363701, 0.369363, 0.357361, 0.365759,
2037 : 0.348566, 0.337062, 0.348190, 0.332330, 0.359001, 0.335836,
2038 : 0.339154, 0.335599, 0.336035, 0.335204, 0.353440, 0.337836,
2039 : 0.333874, 0.307120, 0.294963, 0.324978, 0.313359, 0.317688,
2040 : 0.323758, 0.319304, 0.335317, 0.301765, 0.317257, 0.356331,
2041 : 0.323119, 0.297732, 0.303188, 0.326102, 0.316467, 0.294728,
2042 : 0.308135, 0.288857, 0.325692, 0.312493, 0.291100, 0.325921,
2043 : 0.313317, 0.295980, 0.308481, 0.328380, 0.313081, 0.296763,
2044 : 0.295431, 0.317325, 0.320462, 0.286918, 0.316035, 0.335208,
2045 : 0.283193, 0.333945, 0.292534, 0.294164, 0.330680, 0.296992,
2046 : 0.285509, 0.317260, 0.298269, 0.311299, 0.312129, 0.286822,
2047 : 0.287442, 0.319139, 0.283314, 0.318454, 0.297727, 0.301597,
2048 : 0.282483, 0.294792, 0.305569, 0.290957, 0.297817, 0.282908,
2049 : 0.272401, 0.305584, 0.300220, 0.297020, 0.298781, 0.278008,
2050 : 0.277727, 0.323777, 0.287419, 0.342074, 0.287259, 0.303658,
2051 : 0.302668, 0.279622, 0.280586, 0.313630, 0.276068, 0.257051,
2052 : 0.309996, 0.265534, 0.297138, 0.281738, 0.294610, 0.292882,
2053 : 0.286860, 0.312686, 0.293244, 0.293744, 0.271375, 0.278734,
2054 : 0.280308, 0.304739, 0.287907, 0.285261, 0.311180, 0.313476,
2055 : 0.289660, 0.289282, 0.319505, 0.271285, 0.272008, 0.289245,
2056 : 0.281038, 0.285284, 0.295836, 0.281416, 0.283501, 0.295417,
2057 : 0.304372, 0.297764, 0.291378, 0.321530, 0.315604, 0.329507,
2058 : 0.282609, 0.275576, 0.283721, 0.311714, 0.283760, 0.273188,
2059 : 0.312193, 0.264347, 0.281532, 0.301226, 0.281718, 0.336408,
2060 : 0.283157, 0.332010, 0.289974, 0.290256, 0.301569, 0.332228,
2061 : 0.288282, 0.326339, 0.313653, 0.300361, 0.289470, 0.264830,
2062 : 0.298659, 0.272359, 0.278878, 0.306001, 0.328168, 0.294991,
2063 : 0.327737, 0.278056, 0.302435, 0.284183, 0.279270, 0.307279,
2064 : 0.307917, 0.315196, 0.283803, 0.313333, 0.315730, 0.304818,
2065 : 0.307171, 0.295223, 0.333741, 0.346911, 0.310143, 0.336686,
2066 : 0.275459, 0.334781, 0.295405, 0.275816, 0.301678, 0.327242,
2067 : 0.320717, 0.309230, 0.292145, 0.294489, 0.305088, 0.300969,
2068 : 0.277438, 0.326159, 0.297065, 0.301177, 0.303843, 0.275382,
2069 : 0.304019, 0.284166, 0.289610, 0.331611, 0.317131, 0.310880,
2070 : 0.360456, 0.294052, 0.342694, 0.327166, 0.336797, 0.298040,
2071 : 0.295767, 0.260053, 0.325544, 0.335310, 0.320182, 0.301072,
2072 : 0.313117, 0.283407, 0.299206, 0.293525, 0.305067, 0.255978,
2073 : 0.327055, 0.316382, 0.317700, 0.278993, 0.283120, 0.314000,
2074 : 0.274396, 0.291208, 0.348813, 0.319603, 0.313076, 0.289155,
2075 : 0.343988, 0.311426, 0.322896, 0.328726, 0.337062, 0.389093,
2076 : 0.284122, 0.312184, 0.304008, 0.319170, 0.320778, 0.288390,
2077 : 0.337272, 0.356273, 0.343310, 0.312209, 0.330709, 0.297977,
2078 : 0.346146, 0.369162, 0.324385, 0.339831, 0.337037, 0.318739,
2079 : 0.343157, 0.277720, 0.368407, 0.321330, 0.338997, 0.314220,
2080 : 0.328861, 0.321824, 0.328013, 0.356925, 0.359144, 0.296314,
2081 : 0.345415, 0.396711, 0.347032, 0.294928, 0.343799, 0.322331,
2082 : 0.328656, 0.326098, 0.337338, 0.337038, 0.300179, 0.351391,
2083 : 0.324337, 0.330896, 0.302842, 0.310522, 0.337052, 0.359989,
2084 : 0.383250, 0.359355, 0.315382, 0.333113, 0.342598, 0.355348,
2085 : 0.320751, 0.320475, 0.351762, 0.351475, 0.338358, 0.326153,
2086 : 0.302507, 0.340048, 0.318685, 0.381646, 0.339320, 0.299453,
2087 : 0.426599, 0.393515, 0.353929, 0.328435, 0.413976, 0.292558,
2088 : 0.379340, 0.358344, 0.409259, 0.313821, 0.336675, 0.324521,
2089 : 0.408382, 0.346273, 0.312939, 0.362453, 0.343152, 0.330577,
2090 : 0.332831, 0.353299, 0.347745, 0.334818, 0.332234, 0.385585,
2091 : 0.395483, 0.395316, 0.326972, 0.349434, 0.317604, 0.328980,
2092 : 0.375056, 0.317290, 0.357083, 0.346165, 0.310444, 0.356873,
2093 : 0.359523, 0.365308, 0.365122, 0.383685, 0.370975, 0.396928,
2094 : 0.407654, 0.307755, 0.323033, 0.350580, 0.345231, 0.342462,
2095 : 0.400000, 0.318309, 0.403570, 0.322856, 0.383053, 0.422252,
2096 : 0.386112, 0.364314, 0.434375, 0.334629};
2097 0 : wDD = 0.379611;
2098 0 : wND = 0.496961;
2099 0 : wSD = -1;
2100 :
2101 0 : Mmin = bin[0];
2102 0 : Mmax = bin[nbin];
2103 0 : if(M<Mmin || M>Mmax) return kTRUE;
2104 :
2105 : Int_t ibin=nbin-1;
2106 0 : for(Int_t i=1; i<=nbin; i++)
2107 0 : if(M<=bin[i]) {
2108 0 : ibin=i-1;
2109 : // printf("Mi> %f && Mi< %f\n", bin[i-1], bin[i]);
2110 0 : break;
2111 : }
2112 0 : wSD=w[ibin];
2113 : return kTRUE;
2114 0 : }
2115 0 : else if(TMath::Abs(fEnergyCMS-2760)<1 ){
2116 :
2117 : const Int_t nbin=400;
2118 0 : Double_t bin[]={
2119 : 1.080000, 1.577300, 2.074600, 2.571900, 3.069200, 3.566500,
2120 : 4.063800, 4.561100, 5.058400, 5.555700, 6.053000, 6.550300,
2121 : 7.047600, 7.544900, 8.042200, 8.539500, 9.036800, 9.534100,
2122 : 10.031400, 10.528700, 11.026000, 11.523300, 12.020600, 12.517900,
2123 : 13.015200, 13.512500, 14.009800, 14.507100, 15.004400, 15.501700,
2124 : 15.999000, 16.496300, 16.993600, 17.490900, 17.988200, 18.485500,
2125 : 18.982800, 19.480100, 19.977400, 20.474700, 20.972000, 21.469300,
2126 : 21.966600, 22.463900, 22.961200, 23.458500, 23.955800, 24.453100,
2127 : 24.950400, 25.447700, 25.945000, 26.442300, 26.939600, 27.436900,
2128 : 27.934200, 28.431500, 28.928800, 29.426100, 29.923400, 30.420700,
2129 : 30.918000, 31.415300, 31.912600, 32.409900, 32.907200, 33.404500,
2130 : 33.901800, 34.399100, 34.896400, 35.393700, 35.891000, 36.388300,
2131 : 36.885600, 37.382900, 37.880200, 38.377500, 38.874800, 39.372100,
2132 : 39.869400, 40.366700, 40.864000, 41.361300, 41.858600, 42.355900,
2133 : 42.853200, 43.350500, 43.847800, 44.345100, 44.842400, 45.339700,
2134 : 45.837000, 46.334300, 46.831600, 47.328900, 47.826200, 48.323500,
2135 : 48.820800, 49.318100, 49.815400, 50.312700, 50.810000, 51.307300,
2136 : 51.804600, 52.301900, 52.799200, 53.296500, 53.793800, 54.291100,
2137 : 54.788400, 55.285700, 55.783000, 56.280300, 56.777600, 57.274900,
2138 : 57.772200, 58.269500, 58.766800, 59.264100, 59.761400, 60.258700,
2139 : 60.756000, 61.253300, 61.750600, 62.247900, 62.745200, 63.242500,
2140 : 63.739800, 64.237100, 64.734400, 65.231700, 65.729000, 66.226300,
2141 : 66.723600, 67.220900, 67.718200, 68.215500, 68.712800, 69.210100,
2142 : 69.707400, 70.204700, 70.702000, 71.199300, 71.696600, 72.193900,
2143 : 72.691200, 73.188500, 73.685800, 74.183100, 74.680400, 75.177700,
2144 : 75.675000, 76.172300, 76.669600, 77.166900, 77.664200, 78.161500,
2145 : 78.658800, 79.156100, 79.653400, 80.150700, 80.648000, 81.145300,
2146 : 81.642600, 82.139900, 82.637200, 83.134500, 83.631800, 84.129100,
2147 : 84.626400, 85.123700, 85.621000, 86.118300, 86.615600, 87.112900,
2148 : 87.610200, 88.107500, 88.604800, 89.102100, 89.599400, 90.096700,
2149 : 90.594000, 91.091300, 91.588600, 92.085900, 92.583200, 93.080500,
2150 : 93.577800, 94.075100, 94.572400, 95.069700, 95.567000, 96.064300,
2151 : 96.561600, 97.058900, 97.556200, 98.053500, 98.550800, 99.048100,
2152 : 99.545400, 100.042700, 100.540000, 101.037300, 101.534600, 102.031900,
2153 : 102.529200, 103.026500, 103.523800, 104.021100, 104.518400, 105.015700,
2154 : 105.513000, 106.010300, 106.507600, 107.004900, 107.502200, 107.999500,
2155 : 108.496800, 108.994100, 109.491400, 109.988700, 110.486000, 110.983300,
2156 : 111.480600, 111.977900, 112.475200, 112.972500, 113.469800, 113.967100,
2157 : 114.464400, 114.961700, 115.459000, 115.956300, 116.453600, 116.950900,
2158 : 117.448200, 117.945500, 118.442800, 118.940100, 119.437400, 119.934700,
2159 : 120.432000, 120.929300, 121.426600, 121.923900, 122.421200, 122.918500,
2160 : 123.415800, 123.913100, 124.410400, 124.907700, 125.405000, 125.902300,
2161 : 126.399600, 126.896900, 127.394200, 127.891500, 128.388800, 128.886100,
2162 : 129.383400, 129.880700, 130.378000, 130.875300, 131.372600, 131.869900,
2163 : 132.367200, 132.864500, 133.361800, 133.859100, 134.356400, 134.853700,
2164 : 135.351000, 135.848300, 136.345600, 136.842900, 137.340200, 137.837500,
2165 : 138.334800, 138.832100, 139.329400, 139.826700, 140.324000, 140.821300,
2166 : 141.318600, 141.815900, 142.313200, 142.810500, 143.307800, 143.805100,
2167 : 144.302400, 144.799700, 145.297000, 145.794300, 146.291600, 146.788900,
2168 : 147.286200, 147.783500, 148.280800, 148.778100, 149.275400, 149.772700,
2169 : 150.270000, 150.767300, 151.264600, 151.761900, 152.259200, 152.756500,
2170 : 153.253800, 153.751100, 154.248400, 154.745700, 155.243000, 155.740300,
2171 : 156.237600, 156.734900, 157.232200, 157.729500, 158.226800, 158.724100,
2172 : 159.221400, 159.718700, 160.216000, 160.713300, 161.210600, 161.707900,
2173 : 162.205200, 162.702500, 163.199800, 163.697100, 164.194400, 164.691700,
2174 : 165.189000, 165.686300, 166.183600, 166.680900, 167.178200, 167.675500,
2175 : 168.172800, 168.670100, 169.167400, 169.664700, 170.162000, 170.659300,
2176 : 171.156600, 171.653900, 172.151200, 172.648500, 173.145800, 173.643100,
2177 : 174.140400, 174.637700, 175.135000, 175.632300, 176.129600, 176.626900,
2178 : 177.124200, 177.621500, 178.118800, 178.616100, 179.113400, 179.610700,
2179 : 180.108000, 180.605300, 181.102600, 181.599900, 182.097200, 182.594500,
2180 : 183.091800, 183.589100, 184.086400, 184.583700, 185.081000, 185.578300,
2181 : 186.075600, 186.572900, 187.070200, 187.567500, 188.064800, 188.562100,
2182 : 189.059400, 189.556700, 190.054000, 190.551300, 191.048600, 191.545900,
2183 : 192.043200, 192.540500, 193.037800, 193.535100, 194.032400, 194.529700,
2184 : 195.027000, 195.524300, 196.021600, 196.518900, 197.016200, 197.513500,
2185 : 198.010800, 198.508100, 199.005400, 199.502700, 200.000000};
2186 0 : Double_t w[]={
2187 : 1.000000, 0.692593, 0.673384, 0.666273, 0.657285, 0.637723,
2188 : 0.625881, 0.643590, 0.606100, 0.589007, 0.567824, 0.578705,
2189 : 0.538530, 0.517937, 0.528278, 0.515918, 0.539461, 0.466186,
2190 : 0.489869, 0.468402, 0.465017, 0.453336, 0.460769, 0.474638,
2191 : 0.456347, 0.434471, 0.427478, 0.435435, 0.410934, 0.366431,
2192 : 0.382240, 0.379513, 0.394249, 0.386837, 0.353103, 0.382138,
2193 : 0.377497, 0.389479, 0.378736, 0.347933, 0.354605, 0.352077,
2194 : 0.324443, 0.358792, 0.339968, 0.359052, 0.330734, 0.318291,
2195 : 0.333703, 0.358644, 0.335819, 0.332213, 0.309051, 0.309975,
2196 : 0.331626, 0.304407, 0.309819, 0.312097, 0.312462, 0.320411,
2197 : 0.280401, 0.302311, 0.315863, 0.281479, 0.310003, 0.296911,
2198 : 0.313676, 0.281071, 0.294163, 0.306500, 0.283462, 0.274867,
2199 : 0.307149, 0.270555, 0.282264, 0.287373, 0.307849, 0.278675,
2200 : 0.286990, 0.278269, 0.300105, 0.286799, 0.265674, 0.275140,
2201 : 0.285702, 0.257352, 0.267714, 0.248204, 0.252220, 0.255678,
2202 : 0.282946, 0.268257, 0.282375, 0.262675, 0.275564, 0.248345,
2203 : 0.236259, 0.291914, 0.259936, 0.241338, 0.267389, 0.285044,
2204 : 0.289419, 0.253594, 0.284568, 0.231840, 0.260008, 0.268527,
2205 : 0.275363, 0.224115, 0.281260, 0.257913, 0.295152, 0.264399,
2206 : 0.232287, 0.282533, 0.223431, 0.255756, 0.244471, 0.221695,
2207 : 0.272450, 0.284244, 0.253682, 0.270717, 0.275403, 0.240323,
2208 : 0.245081, 0.241859, 0.216340, 0.244789, 0.220291, 0.238478,
2209 : 0.224691, 0.244058, 0.266117, 0.271478, 0.242012, 0.267321,
2210 : 0.248494, 0.253343, 0.255606, 0.235458, 0.241079, 0.233223,
2211 : 0.226813, 0.259224, 0.234239, 0.258606, 0.210892, 0.238186,
2212 : 0.243271, 0.222678, 0.213437, 0.273939, 0.247966, 0.232548,
2213 : 0.263438, 0.222089, 0.272111, 0.248818, 0.244295, 0.238368,
2214 : 0.236908, 0.248776, 0.232604, 0.231194, 0.227117, 0.231152,
2215 : 0.282140, 0.229778, 0.232631, 0.261794, 0.216633, 0.253471,
2216 : 0.242157, 0.227406, 0.269335, 0.230547, 0.210618, 0.285872,
2217 : 0.248776, 0.229875, 0.242728, 0.227388, 0.220567, 0.222062,
2218 : 0.235950, 0.224087, 0.228895, 0.208287, 0.235999, 0.208696,
2219 : 0.230367, 0.267667, 0.220484, 0.233402, 0.233815, 0.250455,
2220 : 0.253120, 0.219556, 0.230980, 0.236661, 0.222395, 0.226111,
2221 : 0.198315, 0.210555, 0.202813, 0.208594, 0.235976, 0.221490,
2222 : 0.243059, 0.204901, 0.216987, 0.229039, 0.231466, 0.221975,
2223 : 0.231220, 0.253638, 0.250448, 0.260291, 0.328345, 0.205739,
2224 : 0.222014, 0.251513, 0.279427, 0.270506, 0.248409, 0.222472,
2225 : 0.291632, 0.227796, 0.248769, 0.276896, 0.214742, 0.200139,
2226 : 0.230693, 0.226031, 0.268900, 0.185160, 0.245353, 0.205843,
2227 : 0.231155, 0.219122, 0.214811, 0.199763, 0.274179, 0.217598,
2228 : 0.274988, 0.237244, 0.211820, 0.225459, 0.252799, 0.235948,
2229 : 0.224986, 0.245385, 0.237770, 0.213373, 0.229737, 0.215487,
2230 : 0.234453, 0.249684, 0.239435, 0.250422, 0.257194, 0.214762,
2231 : 0.212266, 0.228988, 0.253798, 0.201607, 0.239946, 0.205245,
2232 : 0.231670, 0.212774, 0.206768, 0.231563, 0.189388, 0.227926,
2233 : 0.227203, 0.237754, 0.221628, 0.211138, 0.203322, 0.200985,
2234 : 0.231780, 0.220294, 0.232686, 0.234243, 0.218264, 0.255870,
2235 : 0.223213, 0.238670, 0.213713, 0.213064, 0.246700, 0.233446,
2236 : 0.221503, 0.206767, 0.200722, 0.226179, 0.237425, 0.239229,
2237 : 0.238611, 0.240419, 0.247806, 0.215923, 0.205298, 0.232778,
2238 : 0.272312, 0.226773, 0.258103, 0.223287, 0.269404, 0.203398,
2239 : 0.223782, 0.204213, 0.229664, 0.234040, 0.228419, 0.203936,
2240 : 0.263686, 0.199141, 0.236127, 0.214058, 0.204611, 0.224324,
2241 : 0.292140, 0.190735, 0.235157, 0.213018, 0.257085, 0.190554,
2242 : 0.203197, 0.213044, 0.237023, 0.214243, 0.193562, 0.262403,
2243 : 0.206256, 0.221396, 0.233588, 0.256611, 0.249731, 0.226683,
2244 : 0.199330, 0.251026, 0.222596, 0.201941, 0.186374, 0.221038,
2245 : 0.196555, 0.222560, 0.299419, 0.231979, 0.242924, 0.198310,
2246 : 0.217628, 0.235458, 0.278595, 0.218624, 0.277305, 0.239109,
2247 : 0.205600, 0.253715, 0.221173, 0.218195, 0.277647, 0.241974,
2248 : 0.268748, 0.268128, 0.292258, 0.249420, 0.191034, 0.219506,
2249 : 0.216502, 0.250677, 0.193386, 0.201310, 0.259464, 0.255351,
2250 : 0.269628, 0.221063, 0.294079, 0.196726, 0.233634, 0.221870,
2251 : 0.216236, 0.197259, 0.247433, 0.272765, 0.294079, 0.236336,
2252 : 0.206396, 0.238524, 0.247846, 0.269519, 0.237141, 0.230611,
2253 : 0.201712, 0.242225, 0.255565, 0.258738};
2254 0 : wDD = 0.512813;
2255 0 : wND = 0.518820;
2256 0 : wSD = -1;
2257 :
2258 0 : Mmin = bin[0];
2259 0 : Mmax = bin[nbin];
2260 0 : if(M<Mmin || M>Mmax) return kTRUE;
2261 :
2262 : Int_t ibin=nbin-1;
2263 0 : for(Int_t i=1; i<=nbin; i++)
2264 0 : if(M<=bin[i]) {
2265 0 : ibin=i-1;
2266 : // printf("Mi> %f && Mi< %f\n", bin[i-1], bin[i]);
2267 0 : break;
2268 : }
2269 0 : wSD=w[ibin];
2270 : return kTRUE;
2271 0 : }
2272 :
2273 :
2274 0 : else if(TMath::Abs(fEnergyCMS-7000)<1 ){
2275 : const Int_t nbin=400;
2276 0 : Double_t bin[]={
2277 : 1.080000, 1.577300, 2.074600, 2.571900, 3.069200, 3.566500,
2278 : 4.063800, 4.561100, 5.058400, 5.555700, 6.053000, 6.550300,
2279 : 7.047600, 7.544900, 8.042200, 8.539500, 9.036800, 9.534100,
2280 : 10.031400, 10.528700, 11.026000, 11.523300, 12.020600, 12.517900,
2281 : 13.015200, 13.512500, 14.009800, 14.507100, 15.004400, 15.501700,
2282 : 15.999000, 16.496300, 16.993600, 17.490900, 17.988200, 18.485500,
2283 : 18.982800, 19.480100, 19.977400, 20.474700, 20.972000, 21.469300,
2284 : 21.966600, 22.463900, 22.961200, 23.458500, 23.955800, 24.453100,
2285 : 24.950400, 25.447700, 25.945000, 26.442300, 26.939600, 27.436900,
2286 : 27.934200, 28.431500, 28.928800, 29.426100, 29.923400, 30.420700,
2287 : 30.918000, 31.415300, 31.912600, 32.409900, 32.907200, 33.404500,
2288 : 33.901800, 34.399100, 34.896400, 35.393700, 35.891000, 36.388300,
2289 : 36.885600, 37.382900, 37.880200, 38.377500, 38.874800, 39.372100,
2290 : 39.869400, 40.366700, 40.864000, 41.361300, 41.858600, 42.355900,
2291 : 42.853200, 43.350500, 43.847800, 44.345100, 44.842400, 45.339700,
2292 : 45.837000, 46.334300, 46.831600, 47.328900, 47.826200, 48.323500,
2293 : 48.820800, 49.318100, 49.815400, 50.312700, 50.810000, 51.307300,
2294 : 51.804600, 52.301900, 52.799200, 53.296500, 53.793800, 54.291100,
2295 : 54.788400, 55.285700, 55.783000, 56.280300, 56.777600, 57.274900,
2296 : 57.772200, 58.269500, 58.766800, 59.264100, 59.761400, 60.258700,
2297 : 60.756000, 61.253300, 61.750600, 62.247900, 62.745200, 63.242500,
2298 : 63.739800, 64.237100, 64.734400, 65.231700, 65.729000, 66.226300,
2299 : 66.723600, 67.220900, 67.718200, 68.215500, 68.712800, 69.210100,
2300 : 69.707400, 70.204700, 70.702000, 71.199300, 71.696600, 72.193900,
2301 : 72.691200, 73.188500, 73.685800, 74.183100, 74.680400, 75.177700,
2302 : 75.675000, 76.172300, 76.669600, 77.166900, 77.664200, 78.161500,
2303 : 78.658800, 79.156100, 79.653400, 80.150700, 80.648000, 81.145300,
2304 : 81.642600, 82.139900, 82.637200, 83.134500, 83.631800, 84.129100,
2305 : 84.626400, 85.123700, 85.621000, 86.118300, 86.615600, 87.112900,
2306 : 87.610200, 88.107500, 88.604800, 89.102100, 89.599400, 90.096700,
2307 : 90.594000, 91.091300, 91.588600, 92.085900, 92.583200, 93.080500,
2308 : 93.577800, 94.075100, 94.572400, 95.069700, 95.567000, 96.064300,
2309 : 96.561600, 97.058900, 97.556200, 98.053500, 98.550800, 99.048100,
2310 : 99.545400, 100.042700, 100.540000, 101.037300, 101.534600, 102.031900,
2311 : 102.529200, 103.026500, 103.523800, 104.021100, 104.518400, 105.015700,
2312 : 105.513000, 106.010300, 106.507600, 107.004900, 107.502200, 107.999500,
2313 : 108.496800, 108.994100, 109.491400, 109.988700, 110.486000, 110.983300,
2314 : 111.480600, 111.977900, 112.475200, 112.972500, 113.469800, 113.967100,
2315 : 114.464400, 114.961700, 115.459000, 115.956300, 116.453600, 116.950900,
2316 : 117.448200, 117.945500, 118.442800, 118.940100, 119.437400, 119.934700,
2317 : 120.432000, 120.929300, 121.426600, 121.923900, 122.421200, 122.918500,
2318 : 123.415800, 123.913100, 124.410400, 124.907700, 125.405000, 125.902300,
2319 : 126.399600, 126.896900, 127.394200, 127.891500, 128.388800, 128.886100,
2320 : 129.383400, 129.880700, 130.378000, 130.875300, 131.372600, 131.869900,
2321 : 132.367200, 132.864500, 133.361800, 133.859100, 134.356400, 134.853700,
2322 : 135.351000, 135.848300, 136.345600, 136.842900, 137.340200, 137.837500,
2323 : 138.334800, 138.832100, 139.329400, 139.826700, 140.324000, 140.821300,
2324 : 141.318600, 141.815900, 142.313200, 142.810500, 143.307800, 143.805100,
2325 : 144.302400, 144.799700, 145.297000, 145.794300, 146.291600, 146.788900,
2326 : 147.286200, 147.783500, 148.280800, 148.778100, 149.275400, 149.772700,
2327 : 150.270000, 150.767300, 151.264600, 151.761900, 152.259200, 152.756500,
2328 : 153.253800, 153.751100, 154.248400, 154.745700, 155.243000, 155.740300,
2329 : 156.237600, 156.734900, 157.232200, 157.729500, 158.226800, 158.724100,
2330 : 159.221400, 159.718700, 160.216000, 160.713300, 161.210600, 161.707900,
2331 : 162.205200, 162.702500, 163.199800, 163.697100, 164.194400, 164.691700,
2332 : 165.189000, 165.686300, 166.183600, 166.680900, 167.178200, 167.675500,
2333 : 168.172800, 168.670100, 169.167400, 169.664700, 170.162000, 170.659300,
2334 : 171.156600, 171.653900, 172.151200, 172.648500, 173.145800, 173.643100,
2335 : 174.140400, 174.637700, 175.135000, 175.632300, 176.129600, 176.626900,
2336 : 177.124200, 177.621500, 178.118800, 178.616100, 179.113400, 179.610700,
2337 : 180.108000, 180.605300, 181.102600, 181.599900, 182.097200, 182.594500,
2338 : 183.091800, 183.589100, 184.086400, 184.583700, 185.081000, 185.578300,
2339 : 186.075600, 186.572900, 187.070200, 187.567500, 188.064800, 188.562100,
2340 : 189.059400, 189.556700, 190.054000, 190.551300, 191.048600, 191.545900,
2341 : 192.043200, 192.540500, 193.037800, 193.535100, 194.032400, 194.529700,
2342 : 195.027000, 195.524300, 196.021600, 196.518900, 197.016200, 197.513500,
2343 : 198.010800, 198.508100, 199.005400, 199.502700, 200.000000};
2344 0 : Double_t w[]={
2345 : 1.000000, 0.798640, 0.770197, 0.717990, 0.699434, 0.692093,
2346 : 0.620940, 0.655294, 0.636496, 0.633483, 0.600518, 0.594548,
2347 : 0.550759, 0.550413, 0.528926, 0.525877, 0.506701, 0.504117,
2348 : 0.486973, 0.492930, 0.461833, 0.488695, 0.438874, 0.460274,
2349 : 0.451551, 0.454424, 0.449270, 0.387571, 0.427556, 0.428740,
2350 : 0.401682, 0.402260, 0.408961, 0.395417, 0.387426, 0.378602,
2351 : 0.357826, 0.359125, 0.348494, 0.363710, 0.356117, 0.340363,
2352 : 0.337637, 0.396084, 0.341777, 0.340551, 0.348838, 0.344014,
2353 : 0.340468, 0.317654, 0.355584, 0.326023, 0.373416, 0.312298,
2354 : 0.326522, 0.290540, 0.335557, 0.318689, 0.327544, 0.319501,
2355 : 0.331754, 0.312728, 0.282263, 0.274937, 0.303867, 0.307820,
2356 : 0.289344, 0.268934, 0.288908, 0.290018, 0.291369, 0.295242,
2357 : 0.289067, 0.277685, 0.267957, 0.267559, 0.320229, 0.265060,
2358 : 0.305931, 0.305352, 0.262064, 0.281879, 0.287780, 0.270033,
2359 : 0.270814, 0.276667, 0.271531, 0.275283, 0.258189, 0.287969,
2360 : 0.251247, 0.301527, 0.267230, 0.245860, 0.293125, 0.253421,
2361 : 0.272396, 0.243637, 0.236206, 0.278452, 0.246544, 0.263165,
2362 : 0.230484, 0.231102, 0.258281, 0.244707, 0.270111, 0.248295,
2363 : 0.246942, 0.245592, 0.272766, 0.254243, 0.199647, 0.262590,
2364 : 0.226710, 0.243836, 0.214153, 0.233535, 0.235207, 0.234245,
2365 : 0.247698, 0.248379, 0.241463, 0.265230, 0.223242, 0.236191,
2366 : 0.252700, 0.231865, 0.228352, 0.229101, 0.237385, 0.246485,
2367 : 0.254706, 0.245565, 0.224841, 0.257301, 0.240968, 0.202379,
2368 : 0.236943, 0.241683, 0.220809, 0.219014, 0.213015, 0.244470,
2369 : 0.221042, 0.198996, 0.225295, 0.264466, 0.200600, 0.228143,
2370 : 0.250138, 0.216784, 0.268002, 0.275734, 0.218144, 0.229866,
2371 : 0.235443, 0.208909, 0.215067, 0.189743, 0.216741, 0.242686,
2372 : 0.200234, 0.218882, 0.245991, 0.222815, 0.219576, 0.209773,
2373 : 0.205247, 0.203855, 0.202534, 0.192536, 0.223447, 0.225810,
2374 : 0.220583, 0.248421, 0.223424, 0.206033, 0.203130, 0.221225,
2375 : 0.223763, 0.239216, 0.252311, 0.206893, 0.228461, 0.233591,
2376 : 0.201512, 0.242382, 0.215147, 0.232578, 0.207116, 0.239213,
2377 : 0.196215, 0.184367, 0.235135, 0.189768, 0.274084, 0.206267,
2378 : 0.203283, 0.248828, 0.250214, 0.217465, 0.232080, 0.215150,
2379 : 0.207936, 0.176789, 0.191338, 0.188655, 0.196181, 0.223691,
2380 : 0.254257, 0.216874, 0.213536, 0.221399, 0.192024, 0.212534,
2381 : 0.218169, 0.226635, 0.201191, 0.212700, 0.211634, 0.232353,
2382 : 0.223636, 0.208605, 0.249132, 0.183681, 0.221842, 0.187136,
2383 : 0.203772, 0.249575, 0.217713, 0.205193, 0.193941, 0.219707,
2384 : 0.264226, 0.182105, 0.207183, 0.220845, 0.236571, 0.182390,
2385 : 0.212914, 0.186266, 0.195361, 0.217665, 0.204527, 0.188804,
2386 : 0.222832, 0.191193, 0.201073, 0.185616, 0.195011, 0.200828,
2387 : 0.181434, 0.233707, 0.202925, 0.211992, 0.173100, 0.205258,
2388 : 0.182695, 0.216520, 0.218422, 0.209490, 0.211257, 0.215801,
2389 : 0.220652, 0.211409, 0.195731, 0.194957, 0.198888, 0.234926,
2390 : 0.221377, 0.229822, 0.176700, 0.172322, 0.212265, 0.206133,
2391 : 0.170355, 0.253305, 0.198688, 0.240043, 0.225384, 0.174729,
2392 : 0.195820, 0.200093, 0.196912, 0.212308, 0.200490, 0.240415,
2393 : 0.159744, 0.173686, 0.223853, 0.213604, 0.193779, 0.179328,
2394 : 0.180873, 0.237481, 0.179640, 0.235857, 0.202847, 0.167690,
2395 : 0.245093, 0.215504, 0.205848, 0.184408, 0.201626, 0.209651,
2396 : 0.236306, 0.217803, 0.188534, 0.187861, 0.161663, 0.221718,
2397 : 0.152044, 0.190412, 0.173505, 0.208995, 0.174575, 0.180734,
2398 : 0.247704, 0.203388, 0.194184, 0.169679, 0.182703, 0.213402,
2399 : 0.191808, 0.178604, 0.178116, 0.198452, 0.174687, 0.169809,
2400 : 0.222851, 0.156811, 0.229170, 0.210359, 0.178557, 0.248570,
2401 : 0.208536, 0.192571, 0.178912, 0.224505, 0.170822, 0.205492,
2402 : 0.330973, 0.160924, 0.203724, 0.213255, 0.205827, 0.187162,
2403 : 0.181252, 0.191723, 0.238106, 0.182398, 0.202358, 0.212066,
2404 : 0.201255, 0.168159, 0.185219, 0.229176, 0.158222, 0.164092,
2405 : 0.215405, 0.200724, 0.234811, 0.184488, 0.213112, 0.198577,
2406 : 0.219622, 0.160671, 0.179349, 0.206681, 0.206091, 0.165618,
2407 : 0.165203, 0.174442, 0.179287, 0.187318, 0.163481, 0.217440,
2408 : 0.160381, 0.177025, 0.204385, 0.163676, 0.210733, 0.186519,
2409 : 0.230701, 0.191764, 0.185119, 0.190770, 0.242987, 0.186812,
2410 : 0.202906, 0.161935, 0.182426, 0.197922, 0.181475, 0.155903,
2411 : 0.175006, 0.223482, 0.202706, 0.218108};
2412 0 : wDD = 0.207705;
2413 0 : wND = 0.289628;
2414 0 : wSD = -1;
2415 :
2416 0 : Mmin = bin[0];
2417 0 : Mmax = bin[nbin];
2418 :
2419 0 : if(M<Mmin || M>Mmax) return kTRUE;
2420 :
2421 : Int_t ibin=nbin-1;
2422 0 : for(Int_t i=1; i<=nbin; i++)
2423 0 : if(M<=bin[i]) {
2424 0 : ibin=i-1;
2425 : // printf("Mi> %f && Mi< %f\n", bin[i-1], bin[i]);
2426 0 : break;
2427 : }
2428 0 : wSD=w[ibin];
2429 : return kTRUE;
2430 0 : }
2431 :
2432 0 : return kFALSE;
2433 0 : }
2434 :
2435 : Bool_t AliGenPythia::IsFromHeavyFlavor(Int_t ipart)
2436 : {
2437 : // Check if this is a heavy flavor decay product
2438 0 : TParticle * part = (TParticle *) fParticles.At(ipart);
2439 0 : Int_t mpdg = TMath::Abs(part->GetPdgCode());
2440 0 : Int_t mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
2441 : //
2442 : // Light hadron
2443 0 : if (mfl >= 4 && mfl < 6) return kTRUE;
2444 :
2445 0 : Int_t imo = part->GetFirstMother()-1;
2446 : TParticle* pm = part;
2447 : //
2448 : // Heavy flavor hadron produced by generator
2449 0 : while (imo > -1) {
2450 0 : pm = (TParticle*)fParticles.At(imo);
2451 0 : mpdg = TMath::Abs(pm->GetPdgCode());
2452 0 : mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
2453 0 : if ((mfl > 3) && (mfl <6) && mpdg > 400) return kTRUE;
2454 0 : imo = pm->GetFirstMother()-1;
2455 : }
2456 0 : return kFALSE;
2457 0 : }
2458 :
2459 : Bool_t AliGenPythia::CheckDetectorAcceptance(Float_t phi, Float_t eta, Int_t iparticle)
2460 : {
2461 : // check the eta/phi correspond to the detectors acceptance
2462 : // iparticle is the index of the particle to be checked, for PHOS rotation case
2463 0 : if (fCheckPHOS && IsInPHOS (phi,eta,iparticle)) return kTRUE;
2464 0 : else if(fCheckEMCAL && IsInEMCAL (phi,eta)) return kTRUE;
2465 0 : else if(fCheckBarrel && IsInBarrel( eta)) return kTRUE;
2466 0 : else return kFALSE;
2467 0 : }
2468 :
2469 : Bool_t AliGenPythia::TriggerOnSelectedParticles(Int_t np)
2470 : {
2471 : // Select events with fragmentation photon, decay photon, pi0, eta or other hadrons going to PHOS or EMCAL or central barrel,
2472 : // implemented primaryly for kPyJets, but extended to any kind of process.
2473 :
2474 : //printf("Check: frag photon %d, pi0 %d, eta %d, electron %d, hadron %d, decay %d, PHOS %d, EMCAL %d, Barrel %d \n",
2475 : // fFragPhotonInCalo,fPi0InCalo, fEtaInCalo,fEleInCalo,fHadronInCalo,fDecayPhotonInCalo,fCheckPHOS,fCheckEMCAL, fCheckBarrel);
2476 :
2477 0 : Bool_t ok = kFALSE;
2478 0 : for (Int_t i=0; i< np; i++) {
2479 :
2480 0 : TParticle* iparticle = (TParticle *) fParticles.At(i);
2481 :
2482 0 : Int_t pdg = iparticle->GetPdgCode();
2483 0 : Int_t status = iparticle->GetStatusCode();
2484 0 : Int_t imother = iparticle->GetFirstMother() - 1;
2485 :
2486 : TParticle* pmother = 0x0;
2487 : Int_t momStatus = -1;
2488 : Int_t momPdg = -1;
2489 0 : if(imother > 0 ){
2490 0 : pmother = (TParticle *) fParticles.At(imother);
2491 0 : momStatus = pmother->GetStatusCode();
2492 0 : momPdg = pmother->GetPdgCode();
2493 0 : }
2494 :
2495 0 : ok = kFALSE;
2496 :
2497 : //
2498 : // Check the particle type: hadron (not pi0 or eta), electron, decay photon (from pi0 or eta or any), pi0 or eta
2499 : //
2500 : // Hadron
2501 0 : if (fHadronInCalo && status == 1)
2502 : {
2503 0 : if(TMath::Abs(pdg) > 23 && pdg !=221 && pdg != 111) // avoid photons, electrons, muons, neutrinos and eta or pi0
2504 : // (in case neutral mesons were declared stable)
2505 0 : ok = kTRUE;
2506 : }
2507 : //Electron
2508 0 : else if (fEleInCalo && status == 1 && TMath::Abs(pdg) == 11)
2509 : {
2510 0 : ok = kTRUE;
2511 0 : }
2512 : //Fragmentation photon
2513 0 : else if (fFragPhotonInCalo && pdg == 22 && status == 1)
2514 : {
2515 0 : if(momStatus != 11) ok = kTRUE ; // No photon from hadron decay
2516 : }
2517 : // Decay photon
2518 0 : else if (fDecayPhotonInCalo && !fForceNeutralMeson2PhotonDecay && pdg == 22) // pi0 status can be 1 or 11 depending on decay settings, work only for 11
2519 : {
2520 0 : if( momStatus == 11)
2521 : {
2522 : //if(iparticle->Pt() > fTriggerParticleMinPt) printf("Decay photon! pdg %d, status %d, pt %2.2f, mom: pdg %d, pt %2.2f\n",
2523 : // pdg,status,iparticle->Pt(),momPdg,pmother->Pt());
2524 0 : ok = kTRUE ; // photon from hadron decay
2525 :
2526 : //In case only decays from pi0 or eta requested
2527 0 : if(fPi0InCalo && momPdg!=111) ok = kFALSE;
2528 0 : if(fEtaInCalo && momPdg!=221) ok = kFALSE;
2529 : }
2530 :
2531 : }
2532 : // Pi0 or Eta particle
2533 0 : else if ((fPi0InCalo || fEtaInCalo))
2534 : {
2535 0 : if(fDecayPhotonInCalo && !fForceNeutralMeson2PhotonDecay ) continue ;
2536 :
2537 0 : if (fPi0InCalo && pdg == 111) // pi0 status can be 1 or 11 depending on decay settings
2538 : {
2539 : //if(iparticle->Pt() > fTriggerParticleMinPt) printf("Pi0! pdg %d, status %d, pt %2.2f\n",pdg,status,iparticle->Pt());
2540 0 : ok = kTRUE;
2541 0 : }
2542 0 : else if (fEtaInCalo && pdg == 221)
2543 : {
2544 : //if(iparticle->Pt() > fTriggerParticleMinPt) printf("Eta! pdg %d, status %d, pt %2.2f\n",pdg,status,iparticle->Pt());
2545 0 : ok = kTRUE;
2546 0 : }
2547 :
2548 : }// pi0 or eta
2549 :
2550 : //
2551 : // Check that the selected particle is in the calorimeter acceptance
2552 : //
2553 0 : if(ok && iparticle->Pt() > fTriggerParticleMinPt)
2554 : {
2555 : //Just check if the selected particle falls in the acceptance
2556 0 : if(!fForceNeutralMeson2PhotonDecay )
2557 : {
2558 : //printf("\t Check acceptance! \n");
2559 0 : Float_t phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees
2560 0 : Float_t eta =TMath::Abs(iparticle->Eta()); //in calos etamin=-etamax
2561 :
2562 0 : if(CheckDetectorAcceptance(phi,eta,i))
2563 : {
2564 0 : ok =kTRUE;
2565 0 : AliDebug(1,Form("Selected trigger pdg %d, status %d, pt %2.2f, eta %2.2f, phi %2.2f\n",pdg,status,iparticle->Pt(), eta, phi));
2566 : //printf("\t Accept \n");
2567 0 : break;
2568 : }
2569 0 : else ok = kFALSE;
2570 0 : }
2571 : //Mesons have several decay modes, select only those decaying into 2 photons
2572 0 : else if(fForceNeutralMeson2PhotonDecay && (fPi0InCalo || fEtaInCalo))
2573 : {
2574 : // In case we want the pi0/eta trigger,
2575 : // check the decay mode (2 photons)
2576 :
2577 : //printf("\t Force decay 2 gamma\n");
2578 :
2579 0 : Int_t ndaughters = iparticle->GetNDaughters();
2580 0 : if(ndaughters != 2){
2581 0 : ok=kFALSE;
2582 0 : continue;
2583 : }
2584 :
2585 0 : TParticle*d1 = (TParticle *) fParticles.At(iparticle->GetDaughter(0)-1);
2586 0 : TParticle*d2 = (TParticle *) fParticles.At(iparticle->GetDaughter(1)-1);
2587 0 : if(!d1 || !d2) {
2588 0 : ok=kFALSE;
2589 0 : continue;
2590 : }
2591 :
2592 : //iparticle->Print();
2593 : //d1->Print();
2594 : //d2->Print();
2595 :
2596 0 : Int_t pdgD1 = d1->GetPdgCode();
2597 0 : Int_t pdgD2 = d2->GetPdgCode();
2598 : //printf("\t \t 2 daughters pdg = %d - %d\n",pdgD1,pdgD2);
2599 : //printf("mother %d - %d\n",d1->GetFirstMother(),d2->GetFirstMother());
2600 :
2601 0 : if(pdgD1 != 22 || pdgD2 != 22){
2602 0 : ok = kFALSE;
2603 0 : continue;
2604 : }
2605 :
2606 : //printf("\t accept decay\n");
2607 :
2608 : //Trigger on the meson, not on the daughter
2609 0 : if(!fDecayPhotonInCalo){
2610 :
2611 0 : Float_t phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees
2612 0 : Float_t eta =TMath::Abs(iparticle->Eta()); //in calos etamin=-etamax
2613 :
2614 0 : if(CheckDetectorAcceptance(phi,eta,i))
2615 : {
2616 : //printf("\t Accept meson pdg %d\n",pdg);
2617 0 : ok =kTRUE;
2618 0 : AliDebug(1,Form("Selected trigger pdg %d (decay), status %d, pt %2.2f, eta %2.2f, phi %2.2f\n",pdg,status,iparticle->Pt(), eta, phi));
2619 0 : break;
2620 : } else {
2621 0 : ok=kFALSE;
2622 0 : continue;
2623 : }
2624 : }
2625 :
2626 : //printf("Check daughters acceptance\n");
2627 :
2628 : //Trigger on the meson daughters
2629 : //Photon 1
2630 0 : Float_t phi = d1->Phi()*180./TMath::Pi(); //Convert to degrees
2631 0 : Float_t eta =TMath::Abs(d1->Eta()); //in calos etamin=-etamax
2632 0 : if(d1->Pt() > fTriggerParticleMinPt)
2633 : {
2634 : //printf("\t Check acceptance photon 1! \n");
2635 0 : if(CheckDetectorAcceptance(phi,eta,i))
2636 : {
2637 : //printf("\t Accept Photon 1\n");
2638 0 : ok =kTRUE;
2639 0 : AliDebug(1,Form("Selected trigger pdg %d (decay), status %d, pt %2.2f, eta %2.2f, phi %2.2f\n",pdg,status,iparticle->Pt(), eta, phi));
2640 0 : break;
2641 : }
2642 0 : else ok = kFALSE;
2643 0 : } // pt cut
2644 0 : else ok = kFALSE;
2645 :
2646 : //Photon 2
2647 0 : phi = d2->Phi()*180./TMath::Pi(); //Convert to degrees
2648 0 : eta =TMath::Abs(d2->Eta()); //in calos etamin=-etamax
2649 :
2650 0 : if(d2->Pt() > fTriggerParticleMinPt)
2651 : {
2652 : //printf("\t Check acceptance photon 2! \n");
2653 0 : if(CheckDetectorAcceptance(phi,eta,i))
2654 : {
2655 : //printf("\t Accept Photon 2\n");
2656 0 : ok =kTRUE;
2657 0 : AliDebug(1,Form("Selected trigger pdg %d (decay), status %d, pt %2.2f, eta %2.2f, phi %2.2f\n",pdg,status,iparticle->Pt(), eta, phi));
2658 0 : break;
2659 : }
2660 0 : else ok = kFALSE;
2661 0 : } // pt cut
2662 0 : else ok = kFALSE;
2663 0 : } // force 2 photon daughters in pi0/eta decays
2664 0 : else ok = kFALSE;
2665 0 : } else ok = kFALSE; // check acceptance
2666 0 : } // primary loop
2667 :
2668 : //
2669 : // If requested, rotate the particles event in phi to enhance/speed PHOS selection
2670 : // A particle passing all trigger conditions except phi position in PHOS, is used as reference
2671 : //
2672 0 : if(fCheckPHOSeta)
2673 : {
2674 0 : RotatePhi(ok);
2675 0 : }
2676 :
2677 0 : return ok;
2678 0 : }
2679 :
2680 : void AliGenPythia::SetSeed(UInt_t seed)
2681 : {
2682 0 : GetRandom()->SetSeed(seed);
2683 0 : }
|