]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenHijing.cxx
Added const where compiler told me to
[u/mrichter/AliRoot.git] / EVGEN / AliGenHijing.cxx
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 /*
17 $Log$
18 Revision 1.34  2001/12/11 16:55:42  morsch
19 Correct initialization for jet phi-range.
20
21 Revision 1.33  2001/12/05 10:18:51  morsch
22 Possibility of kinematic biasing of jet phi range. (J. Klay)
23
24 Revision 1.32  2001/11/28 13:51:11  morsch
25 Introduce kinematic biasing (etamin, etamax) of jet trigger. Bookkeeping
26 (number of trials) done in AliGenHijingEventHeader.
27
28 Revision 1.31  2001/11/06 12:30:34  morsch
29 Add Boost() method to boost all particles to LHC lab frame. Needed for asymmetric collision systems.
30
31 Revision 1.30  2001/10/21 18:35:56  hristov
32 Several pointers were set to zero in the default constructors to avoid memory management problems
33
34 Revision 1.29  2001/10/15 08:12:24  morsch
35 - Vertex smearing with truncated gaussian.
36 - Store triggered jet info before and after final state radiation into mc-heade
37
38 Revision 1.28  2001/10/08 11:55:25  morsch
39 Store 4-momenta of trigegred jets in event header.
40 Possibility to switch of initial and final state radiation.
41
42 Revision 1.27  2001/10/08 07:13:14  morsch
43 Add setter for minimum transverse momentum of triggered jet.
44
45 Revision 1.26  2001/10/04 08:12:24  morsch
46 Redefinition of stable condition.
47
48 Revision 1.25  2001/07/27 17:09:36  morsch
49 Use local SetTrack, KeepTrack and SetHighWaterMark methods
50 to delegate either to local stack or to stack owned by AliRun.
51 (Piotr Skowronski, A.M.)
52
53 Revision 1.24  2001/07/20 09:34:56  morsch
54 Count the number of spectator neutrons and protons and add information
55 to the event header. (Chiara Oppedisano)
56
57 Revision 1.23  2001/07/13 17:30:22  morsch
58 Derive from AliGenMC.
59
60 Revision 1.22  2001/06/11 13:09:23  morsch
61 - Store cross-Section and number of binary collisions as a function of impact parameter
62 - Pass AliGenHijingEventHeader to gAlice.
63
64 Revision 1.21  2001/02/28 17:35:24  morsch
65 Consider elastic interactions (ks = 1 and ks = 11) as spectator (Chiara Oppedisano)
66
67 Revision 1.20  2001/02/14 15:50:40  hristov
68 The last particle in event marked using SetHighWaterMark
69
70 Revision 1.19  2000/12/21 16:24:06  morsch
71 Coding convention clean-up
72
73 Revision 1.18  2000/12/06 17:46:30  morsch
74 Avoid random numbers 1 and 0.
75
76 Revision 1.17  2000/12/04 11:22:03  morsch
77 Init of sRandom as in 1.15
78
79 Revision 1.16  2000/12/02 11:41:39  morsch
80 Use SetRandom() to initialize random number generator in constructor.
81
82 Revision 1.15  2000/11/30 20:29:02  morsch
83 Initialise static variable sRandom in constructor: sRandom = fRandom;
84
85 Revision 1.14  2000/11/30 07:12:50  alibrary
86 Introducing new Rndm and QA classes
87
88 Revision 1.13  2000/11/09 17:40:27  morsch
89 Possibility to select/unselect spectator protons and neutrons.
90 Method SetSpectators(Int_t spect) added. (FCA, Ch. Oppedisano)
91
92 Revision 1.12  2000/10/20 13:38:38  morsch
93 Debug printouts commented.
94
95 Revision 1.11  2000/10/20 13:22:26  morsch
96 - skip particle type 92 (string)
97 - Charmed and beauty baryions (5122, 4122) are considered as stable consistent with
98   mesons.
99
100 Revision 1.10  2000/10/17 15:10:20  morsch
101 Write first all the parent particles to the stack and then the final state particles.
102
103 Revision 1.9  2000/10/17 13:38:59  morsch
104 Protection against division by zero in EvaluateCrossSection() and KinematicSelection(..)     (FCA)
105
106 Revision 1.8  2000/10/17 12:46:31  morsch
107 Protect EvaluateCrossSections() against division by zero.
108
109 Revision 1.7  2000/10/02 21:28:06  fca
110 Removal of useless dependecies via forward declarations
111
112 Revision 1.6  2000/09/11 13:23:37  morsch
113 Write last seed to file (fortran lun 50) and reed back from same lun using calls to
114 luget_hijing and luset_hijing.
115
116 Revision 1.5  2000/09/07 16:55:40  morsch
117 fHijing->Initialize(); after change of parameters. (Dmitri Yurevitch Peressounko)
118
119 Revision 1.4  2000/07/11 18:24:56  fca
120 Coding convention corrections + few minor bug fixes
121
122 Revision 1.3  2000/06/30 12:08:36  morsch
123 In member data: char* replaced by TString, Init takes care of resizing the strings to
124 8 characters required by Hijing.
125
126 Revision 1.2  2000/06/15 14:15:05  morsch
127 Add possibility for heavy flavor selection: charm and beauty.
128
129 Revision 1.1  2000/06/09 20:47:27  morsch
130 AliGenerator interface class to HIJING using THijing (test version)
131
132 */
133
134
135
136 // Generator using HIJING as an external generator
137 // The main HIJING options are accessable for the user through this interface.
138 // Uses the THijing implementation of TGenerator.
139 //
140 // andreas.morsch@cern.ch
141
142 #include "AliGenHijing.h"
143 #include "AliGenHijingEventHeader.h"
144 #include "AliRun.h"
145 #include "AliPDG.h"
146
147 #include <TArrayI.h>
148 #include <TParticle.h>
149 #include <THijing.h>
150 #include <TGraph.h>
151 #include <TLorentzVector.h>
152
153
154  ClassImp(AliGenHijing)
155
156 AliGenHijing::AliGenHijing()
157                  :AliGenMC()
158 {
159 // Constructor
160   fHijing = 0;
161   fDsigmaDb = 0;
162   fDnDb = 0;
163 }
164
165 AliGenHijing::AliGenHijing(Int_t npart)
166     :AliGenMC(npart)
167 {
168 // Default PbPb collisions at 5. 5 TeV
169 //
170     SetEnergyCMS();
171     SetImpactParameterRange();
172     SetTarget();
173     SetProjectile();
174     SetBoostLHC();
175     SetJetEtaRange();
176     SetJetPhiRange();
177     
178     fKeep       =  0;
179     fQuench     =  1;
180     fShadowing  =  1;
181     fTrigger    =  0;
182     fDecaysOff  =  1;
183     fEvaluate   =  0;
184     fSelectAll  =  0;
185     fFlavor     =  0;
186     fSpectators =  1;
187     fDsigmaDb   =  0;
188     fDnDb       =  0;
189     fPtMinJet   = -2.5;         
190     fRadiation  =  1;
191     fEventVertex.Set(3);
192 //
193 // Set random number generator   
194     sRandom = fRandom;
195     fHijing = 0;
196
197 }
198
199 AliGenHijing::AliGenHijing(const AliGenHijing & Hijing)
200 {
201 // copy constructor
202 }
203
204
205 AliGenHijing::~AliGenHijing()
206 {
207 // Destructor
208     if ( fDsigmaDb) delete  fDsigmaDb;  
209     if ( fDnDb)     delete  fDnDb;  
210 }
211
212 void AliGenHijing::Init()
213 {
214 // Initialisation
215     fFrame.Resize(8);
216     fTarget.Resize(8);
217     fProjectile.Resize(8);
218     
219     SetMC(new THijing(fEnergyCMS, fFrame, fProjectile, fTarget, 
220                       fAProjectile, fZProjectile, fATarget, fZTarget, 
221                       fMinImpactParam, fMaxImpactParam));
222
223     fHijing=(THijing*) fgMCEvGen;
224     fHijing->SetIHPR2(2,  fRadiation);
225     fHijing->SetIHPR2(3,  fTrigger);
226     fHijing->SetIHPR2(4,  fQuench);
227     fHijing->SetIHPR2(6,  fShadowing);
228     fHijing->SetIHPR2(12, fDecaysOff);    
229     fHijing->SetIHPR2(21, fKeep);
230     fHijing->SetHIPR1(10, fPtMinJet);   
231     fHijing->Initialize();
232
233     
234 //
235     if (fEvaluate) EvaluateCrossSections();
236 //
237 //
238 //  Initialize random generator
239 }
240
241 void AliGenHijing::Generate()
242 {
243 // Generate one event
244
245   Float_t polar[3]    =   {0,0,0};
246   Float_t origin[3]   =   {0,0,0};
247   Float_t origin0[3]  =   {0,0,0};
248   Float_t p[3], random[6];
249   Float_t tof;
250
251   static TClonesArray *particles;
252 //  converts from mm/c to s
253   const Float_t kconv = 0.001/2.999792458e8;
254 //
255   Int_t nt  = 0;
256   Int_t jev = 0;
257   Int_t j, kf, ks, imo;
258   kf = 0;
259     
260   if(!particles) particles = new TClonesArray("TParticle",10000);
261     
262   fTrials = 0;
263   for (j = 0;j < 3; j++) origin0[j] = fOrigin[j];
264   if(fVertexSmear == kPerEvent) {
265       Float_t dv[3];
266       dv[2] = 1.e10;
267       while(TMath::Abs(dv[2]) > fCutVertexZ*fOsigma[2]) {
268           Rndm(random,6);
269           for (j=0; j < 3; j++) {
270               dv[j] = fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
271                   TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
272           }
273       }
274       for (j=0; j < 3; j++) origin0[j] += dv[j];
275   } else if (fVertexSmear == kPerTrack) {
276 //          fHijing->SetMSTP(151,0);
277       for (j = 0; j < 3; j++) {
278 //            fHijing->SetPARP(151+j, fOsigma[j]*10.);
279       }
280   }
281   while(1)
282   {
283 //    Generate one event
284 // --------------------------------------------------------------------------
285       fSpecn    = 0;  
286       fSpecp    = 0;
287 // --------------------------------------------------------------------------
288       fHijing->GenerateEvent();
289       fTrials++;
290       if (fTrigger != kNoTrigger) {
291           if (!CheckTrigger()) continue;
292       }
293
294       fHijing->ImportParticles(particles,"All");
295       if (fLHC) Boost(particles);
296       
297       
298       Int_t np = particles->GetEntriesFast();
299       printf("\n **************************************************%d\n",np);
300       Int_t nc = 0;
301       if (np == 0 ) continue;
302       Int_t i;
303       Int_t * newPos = new Int_t[np];
304
305       for (i = 0; i < np; i++) *(newPos+i) = i;
306 //      Get event vertex
307 //
308       TParticle *  iparticle = (TParticle *) particles->At(0);
309       fEventVertex[0] = origin0[0];
310       fEventVertex[1] = origin0[1];     
311       fEventVertex[2] = origin0[2];
312       
313 //
314 //      First write parent particles
315 //
316
317       for (i = 0; i < np; i++) {
318           iparticle       = (TParticle *) particles->At(i);
319 // Is this a parent particle ?
320           if (Stable(iparticle)) continue;
321 //
322         Bool_t  hasMother            =  (iparticle->GetFirstMother()   >=0);
323         Bool_t  selected             =  kTRUE;
324         Bool_t  hasSelectedDaughters =  kFALSE;
325             
326             
327         kf        = iparticle->GetPdgCode();
328         ks        = iparticle->GetStatusCode();
329         if (kf == 92) continue;
330             
331         if (!fSelectAll) selected = KinematicSelection(iparticle, 0)&&SelectFlavor(kf);
332         hasSelectedDaughters = DaughtersSelection(iparticle, particles);
333 //
334 // Put particle on the stack if it is either selected or it is the mother of at least one seleted particle
335 //
336         if (selected || hasSelectedDaughters) {
337             nc++;
338             p[0] = iparticle->Px();
339             p[1] = iparticle->Py();
340             p[2] = iparticle->Pz();
341             origin[0] = origin0[0]+iparticle->Vx()/10;
342             origin[1] = origin0[1]+iparticle->Vy()/10;
343             origin[2] = origin0[2]+iparticle->Vz()/10;
344             tof = kconv*iparticle->T();
345             imo = -1;
346             TParticle* mother = 0;
347             if (hasMother) {
348                 imo = iparticle->GetFirstMother();
349                 mother = (TParticle *) particles->At(imo);
350                 imo = (mother->GetPdgCode() != 92) ? imo =* (newPos+imo) : -1;
351                 
352             }
353 // Put particle on the stack ... 
354 //              printf("\n set track mother: %d %d %d %d %d %d ",i,imo, kf, nt+1, selected, hasSelectedDaughters);
355             
356             SetTrack(0,imo,kf,p,origin,polar, tof,kPPrimary,nt);
357 // ... and keep it there
358             KeepTrack(nt);
359 //
360             *(newPos+i)=nt;
361         } // selected
362       } // particle loop parents
363 //
364 // Now write the final state particles
365 //
366
367       for (i = 0; i<np; i++) {
368         TParticle *  iparticle       = (TParticle *) particles->At(i);
369 // Is this a final state particle ?
370         if (!Stable(iparticle)) continue;
371
372         Bool_t  hasMother            =  (iparticle->GetFirstMother()   >=0);
373         Bool_t  selected             =  kTRUE;
374         kf        = iparticle->GetPdgCode();
375         ks        = iparticle->GetStatusCode();
376 // --------------------------------------------------------------------------
377 // Count spectator neutrons and protons
378         if(ks == 0 || ks == 1 || ks == 10 || ks == 11){
379           if(kf == kNeutron) fSpecn += 1;
380           if(kf == kProton)  fSpecp += 1;
381         }
382 // --------------------------------------------------------------------------
383 //          
384         if (!fSelectAll) {
385             selected = KinematicSelection(iparticle,0)&&SelectFlavor(kf);
386             if (!fSpectators && selected) selected = (ks != 0 && ks != 1 && ks != 10
387                                                       && ks != 11);
388         }
389 //
390 // Put particle on the stack if selected
391 //
392         if (selected) {
393             nc++;
394             p[0] = iparticle->Px();
395             p[1] = iparticle->Py();
396             p[2] = iparticle->Pz();
397             origin[0] = origin0[0]+iparticle->Vx()/10;
398             origin[1] = origin0[1]+iparticle->Vy()/10;
399             origin[2] = origin0[2]+iparticle->Vz()/10;
400             tof = kconv*iparticle->T();
401             imo = -1;
402             TParticle* mother = 0;
403             if (hasMother) {
404                 imo = iparticle->GetFirstMother();
405                 mother = (TParticle *) particles->At(imo);
406                 imo = (mother->GetPdgCode() != 92) ? imo=*(newPos+imo) : -1;
407             }   
408 // Put particle on the stack
409             SetTrack(fTrackIt,imo,kf,p,origin,polar,
410                                                      tof,kPNoProcess,nt);
411             KeepTrack(nt);
412             *(newPos+i)=nt;
413         } // selected
414       } // particle loop final state
415  
416       delete[] newPos;
417
418       printf("\n I've put %i particles on the stack \n",nc);
419       if (nc > 0) {
420           jev += nc;
421           if (jev >= fNpart || fNpart == -1) {
422               fKineBias = Float_t(fNpart)/Float_t(fTrials);
423               printf("\n Trials: %i %i %i\n",fTrials, fNpart, jev);
424               break;
425           }
426       }
427   } // event loop
428   MakeHeader();
429   
430   SetHighWaterMark(nt);
431 }
432
433 void AliGenHijing::KeepFullEvent()
434 {
435     fKeep=1;
436 }
437
438 void AliGenHijing::EvaluateCrossSections()
439 {
440 //     Glauber Calculation of geometrical x-section
441 //
442     Float_t xTot       = 0.;          // barn
443     Float_t xTotHard   = 0.;          // barn 
444     Float_t xPart      = 0.;          // barn
445     Float_t xPartHard  = 0.;          // barn 
446     Float_t sigmaHard  = 0.1;         // mbarn
447     Float_t bMin       = 0.;
448     Float_t bMax       = fHijing->GetHIPR1(34)+fHijing->GetHIPR1(35);
449     const Float_t kdib = 0.2;
450     Int_t   kMax       = Int_t((bMax-bMin)/kdib)+1;
451
452
453     printf("\n Projectile Radius (fm): %f \n",fHijing->GetHIPR1(34));
454     printf("\n Target     Radius (fm): %f \n",fHijing->GetHIPR1(35));    
455     Int_t i;
456     Float_t oldvalue= 0.;
457
458     Float_t* b   = new Float_t[kMax];
459     Float_t* si1 = new Float_t[kMax];    
460     Float_t* si2 = new Float_t[kMax];    
461     
462     for (i = 0; i < kMax; i++)
463     {
464         Float_t xb  = bMin+i*kdib;
465         Float_t ov;
466         ov=fHijing->Profile(xb);
467         Float_t gb  =  2.*0.01*fHijing->GetHIPR1(40)*kdib*xb*(1.-TMath::Exp(-fHijing->GetHINT1(12)*ov));
468         Float_t gbh =  2.*0.01*fHijing->GetHIPR1(40)*kdib*xb*sigmaHard*ov;
469         xTot+=gb;
470         xTotHard += gbh;
471         if (xb > fMinImpactParam && xb < fMaxImpactParam)
472         {
473             xPart += gb;
474             xPartHard += gbh;
475         }
476         
477         if(oldvalue) if ((xTot-oldvalue)/oldvalue<0.0001) break;
478         oldvalue = xTot;
479         printf("\n Total cross section (barn): %d %f %f \n",i, xb, xTot);
480         printf("\n Hard  cross section (barn): %d %f %f \n\n",i, xb, xTotHard);
481         if (i>0) {
482             si1[i] = gb/kdib;
483             si2[i] = gbh/gb;
484             b[i]  = xb;
485         }
486     }
487
488     printf("\n Total cross section (barn): %f \n",xTot);
489     printf("\n Hard  cross section (barn): %f \n \n",xTotHard);
490     printf("\n Partial       cross section (barn): %f %f \n",xPart, xPart/xTot*100.);
491     printf("\n Partial  hard cross section (barn): %f %f \n",xPartHard, xPartHard/xTotHard*100.);
492
493 //  Store result as a graph
494     b[0] = 0;
495     si1[0] = 0;
496     si2[0]=si2[1];
497     
498     fDsigmaDb  = new TGraph(i, b, si1);
499     fDnDb      = new TGraph(i, b, si2);
500 }
501
502 Bool_t AliGenHijing::DaughtersSelection(TParticle* iparticle, TClonesArray* particles)
503 {
504 //
505 // Looks recursively if one of the daughters has been selected
506 //
507 //    printf("\n Consider daughters %d:",iparticle->GetPdgCode());
508     Int_t imin = -1;
509     Int_t imax = -1;
510     Int_t i;
511     Bool_t hasDaughters = (iparticle->GetFirstDaughter() >=0);
512     Bool_t selected = kFALSE;
513     if (hasDaughters) {
514         imin = iparticle->GetFirstDaughter();
515         imax = iparticle->GetLastDaughter();       
516         for (i = imin; i <= imax; i++){
517             TParticle *  jparticle = (TParticle *) particles->At(i);    
518             Int_t ip = jparticle->GetPdgCode();
519             if (KinematicSelection(jparticle,0)&&SelectFlavor(ip)) {
520                 selected=kTRUE; break;
521             }
522             if (DaughtersSelection(jparticle, particles)) {selected=kTRUE; break; }
523         }
524     } else {
525         return kFALSE;
526     }
527     return selected;
528 }
529
530
531 Bool_t AliGenHijing::SelectFlavor(Int_t pid)
532 {
533 // Select flavor of particle
534 // 0: all
535 // 4: charm and beauty
536 // 5: beauty
537     if (fFlavor == 0) return kTRUE;
538     
539     Int_t ifl = TMath::Abs(pid/100);
540     if (ifl > 10) ifl/=10;
541     return (fFlavor == ifl);
542 }
543
544 Bool_t AliGenHijing::Stable(TParticle*  particle)
545 {
546 // Return true for a stable particle
547 //
548     if (particle->GetFirstDaughter() < 0 )
549     {
550         return kTRUE;
551     } else {
552         return kFALSE;
553     }
554 }
555
556
557 void AliGenHijing::Boost(TClonesArray* particles)
558 {
559 //
560 // Boost cms into LHC lab frame
561 //
562     Double_t dy    = - 0.5 * TMath::Log(Double_t(fZProjectile) * Double_t(fATarget) / 
563                                       (Double_t(fZTarget)    * Double_t(fAProjectile)));
564     Double_t beta  = TMath::TanH(dy);
565     Double_t gamma = 1./TMath::Sqrt(1.-beta*beta);
566     Double_t gb    = gamma * beta;
567
568     printf("\n Boosting particles to lab frame %f %f %f", dy, beta, gamma);
569     
570     Int_t i;
571     Int_t np = particles->GetEntriesFast();
572     for (i = 0; i < np; i++) 
573     {
574         TParticle* iparticle = (TParticle*) particles->At(i);
575
576         Double_t e   = iparticle->Energy();
577         Double_t px  = iparticle->Px();
578         Double_t py  = iparticle->Py();
579         Double_t pz  = iparticle->Pz();
580
581         Double_t eb  = gamma * e -      gb * pz;
582         Double_t pzb =   -gb * e +   gamma * pz;
583
584         iparticle->SetMomentum(px, py, pzb, eb);
585     }
586 }
587
588
589 void AliGenHijing::MakeHeader()
590 {
591 // Builds the event header, to be called after each event
592     AliGenEventHeader* header = new AliGenHijingEventHeader("Hijing");
593     ((AliGenHijingEventHeader*) header)->SetNProduced(fHijing->GetNATT());
594     ((AliGenHijingEventHeader*) header)->SetImpactParameter(fHijing->GetHINT1(19));
595     ((AliGenHijingEventHeader*) header)->SetTotalEnergy(fHijing->GetEATT());
596     ((AliGenHijingEventHeader*) header)->SetHardScatters(fHijing->GetJATT());
597     ((AliGenHijingEventHeader*) header)->SetParticipants(fHijing->GetNP(), fHijing->GetNT());
598     ((AliGenHijingEventHeader*) header)->SetCollisions(fHijing->GetN0(),
599                                                        fHijing->GetN01(),
600                                                        fHijing->GetN10(),
601                                                        fHijing->GetN11());
602     ((AliGenHijingEventHeader*) header)->SetSpectators(fSpecn, fSpecp);
603
604 // 4-momentum vectors of the triggered jets.
605 //
606 // Before final state gluon radiation.
607     TLorentzVector* jet1 = new TLorentzVector(fHijing->GetHINT1(21), 
608                                               fHijing->GetHINT1(22),
609                                               fHijing->GetHINT1(23),
610                                               fHijing->GetHINT1(24));
611
612     TLorentzVector* jet2 = new TLorentzVector(fHijing->GetHINT1(31), 
613                                               fHijing->GetHINT1(32),
614                                               fHijing->GetHINT1(33),
615                                               fHijing->GetHINT1(34));
616 // After final state gluon radiation.
617     TLorentzVector* jet3 = new TLorentzVector(fHijing->GetHINT1(26), 
618                                               fHijing->GetHINT1(27),
619                                               fHijing->GetHINT1(28),
620                                               fHijing->GetHINT1(29));
621
622     TLorentzVector* jet4 = new TLorentzVector(fHijing->GetHINT1(36), 
623                                               fHijing->GetHINT1(37),
624                                               fHijing->GetHINT1(38),
625                                               fHijing->GetHINT1(39));
626     ((AliGenHijingEventHeader*) header)->SetJets(jet1, jet2, jet3, jet4);
627 // Bookkeeping for kinematic bias
628     ((AliGenHijingEventHeader*) header)->SetTrials(fTrials);
629 // Event Vertex
630     header->SetPrimaryVertex(fEventVertex);
631     gAlice->SetGenEventHeader(header);    
632 }
633
634 Bool_t AliGenHijing::CheckTrigger()
635 {
636 // Check the kinematic trigger condition
637 // 
638     TLorentzVector* jet1 = new TLorentzVector(fHijing->GetHINT1(26), 
639                                               fHijing->GetHINT1(27),
640                                               fHijing->GetHINT1(28),
641                                               fHijing->GetHINT1(29));
642
643     TLorentzVector* jet2 = new TLorentzVector(fHijing->GetHINT1(36), 
644                                               fHijing->GetHINT1(37),
645                                               fHijing->GetHINT1(38),
646                                               fHijing->GetHINT1(39));
647     Double_t eta1      = jet1->Eta();
648     Double_t eta2      = jet2->Eta();
649     Double_t phi1      = jet1->Phi();
650     Double_t phi2      = jet2->Phi();
651     Bool_t   triggered = kFALSE;
652     //Check eta range first...    
653     if ((eta1 < fEtaMaxJet && eta1 > fEtaMinJet) ||
654         (eta2 < fEtaMaxJet && eta2 > fEtaMinJet))
655     {
656         //Eta is okay, now check phi range
657         if ((phi1 < fPhiMaxJet && phi1 > fPhiMinJet) ||
658             (phi2 < fPhiMaxJet && phi2 > fPhiMinJet))
659         {
660             triggered = kTRUE;
661         }
662     }
663     return triggered;
664 }
665
666 AliGenHijing& AliGenHijing::operator=(const  AliGenHijing& rhs)
667 {
668 // Assignment operator
669     return *this;
670 }
671
672 #ifndef WIN32
673 # define rluget_hijing rluget_hijing_
674 # define rluset_hijing rluset_hijing_
675 # define rlu_hijing    rlu_hijing_
676 # define type_of_call
677 #else
678 # define rluget_hijing RLUGET_HIJING
679 # define rluset_hijing RLUSET_HIJING
680 # define rlu_hijing    RLU_HIJING
681 # define type_of_call _stdcall
682 #endif
683
684
685 extern "C" {
686   void type_of_call rluget_hijing(Int_t & /*lfn*/, Int_t & /*move*/)
687   {printf("Dummy version of rluget_hijing reached\n");}
688
689   void type_of_call rluset_hijing(Int_t & /*lfn*/, Int_t & /*move*/)
690   {printf("Dummy version of rluset_hijing reached\n");}
691
692   Double_t type_of_call rlu_hijing(Int_t & /*idum*/) 
693   {
694       Float_t r;
695       do r=sRandom->Rndm(); while(0 >= r || r >= 1);
696       return r;
697   }
698 }