]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenHBTprocessor.cxx
Several pointers were set to zero in the default constructors to avoid memory managem...
[u/mrichter/AliRoot.git] / EVGEN / AliGenHBTprocessor.cxx
CommitLineData
2b9786f4 1/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
2 * See cxx source for full Copyright notice */
3
4/* $Id$ */
5
6// Implementation of the interface for THBTprocessor
7// Author: Piotr Krzysztof Skowronski <Piotr.Skowronski@cern.ch>
8
4c90f2a0 9// 09.10.2001 Piotr Skowronski
10//
11// Theta - Eta cohernecy facilities added:
12// AliGenerator::SetThetaRange method overriden
13// Static methods added
14// EtaToTheta
15// ThetaToEta
16// DegreesToRadians
17// RadiansToDegrees
18//
19// Class description comments put on proper place
20
2b9786f4 21// 27.09.2001 Piotr Skowronski
22// removing of redefinition of defaults velues
23// in method's implementation.
24//
25
4c90f2a0 26//////////////////////////////////////////////////////////////////////////////////
2b9786f4 27//Wrapper class for "hbt processor" after burner
28//The origibal code is written in fortran by Lanny Ray
29//and is put in the directory $ALICE_ROOT/HBTprocessor
30//Detailed description is on the top of the file hbt_event_processor.f
31//
32//This class can be used ONLY with AliGenCocktailAfterBurner wrapper generator
33//i.e. (in Config.C)
34// ....
35// AliGenCocktailAfterBurner = gener = new AliGenCocktailAfterBurner();
36// gener->SetPhiRange(0, 360); //Set global parameters
37// gener->SetThetaRange(35., 145.);
38// AliGenHIJINGpara *hijing = new AliGenHIJINGpara(10000); //Add some generator
39// hijing->SetMomentumRange(0, 999);
40// gener->AddGenerator(hijing,"HIJING PARAMETRIZATION",1);
41//
42// AliGenHBTprocessor *hbtp = new AliGenHBTprocessor(); //create object
43// hbtp->SetRefControl(2); //set parameters
44// hbtp->SetSwitch1D(1);
45// hbtp->SetR0(6);//fm - radius
46// hbtp->SetLambda(0.7);//chaocity parameter
47// gener->AddAfterBurner(hbtp,"HBT PROCESSOR",1); //add to the main generator
48//
49// gener->Init();
50//
51//CAUTIONS:
52// A) HBT PROCESSOR NEEDS MORE THAN ONE EVENT TO WORK
53// AS MORE AS IT BETTER WORKS
54// B) IT IS ABLE TO "ADD" CORRELATIONS ONLY UP TO TWO PARTICLE TYPES AT ONES
4c90f2a0 55//////////////////////////////////////////////////////////////////////////////////
56
57#include "AliGenHBTprocessor.h"
58#include "TROOT.h"
59#include <iostream.h>
60
61#include "AliRun.h"
62#include "AliStack.h"
63#include "TParticle.h"
64#include "AliGenCocktailAfterBurner.h"
65
66
67
68const Int_t AliGenHBTprocessor::fgkHBTPMAXPART = 100000;
69
70ClassImp(AliGenHBTprocessor)
71
2b9786f4 72
73
74AliGenCocktailAfterBurner* GetGenerator();
75/*******************************************************************/
76
77AliGenHBTprocessor::AliGenHBTprocessor() : AliGenerator(-1)
78{
79 //
80 // Standard constructor
81 // Sets default veues of all parameters
2685bf00 82 fHbtPStatCodes = 0;
2b9786f4 83 SetName("AliGenHBTprocessor");
84 SetTitle("AliGenHBTprocessor");
85
86 sRandom = fRandom;
87 fHBTprocessor = new THBTprocessor();
88
89 fNPDGCodes = 0;
90 DefineParticles();
91
92 SetTrackRejectionFactor();
93 SetRefControl();
94 SetPIDs();
95 SetNPIDtypes();
96 SetDeltap();
97 SetMaxIterations();
98 SetDelChi();
99 SetIRand();
100 SetLambda();
101 SetR1d() ;
102 SetRSide();
103 SetROut() ;
104 SetRLong() ;
105 SetRPerp();
106 SetRParallel();
107 SetR0();
108 SetQ0();
109 SetSwitch1D();
110 SetSwitch3D();
111 SetSwitchType();
112 SetSwitchCoherence();
113 SetSwitchCoulomb();
114 SetSwitchFermiBose();
115 //SetMomentumRange();
116 SetPtRange();
117 SetPxRange();
118 SetPyRange();
119 SetPzRange();
120 SetPhiRange();
121 SetEtaRange();
122 SetNPtBins();
123 SetNPhiBins();
124 SetNEtaBins();
125 SetNPxBins();
126 SetNPyBins();
127 SetNPzBins();
128 SetNBins1DFineMesh();
129 SetBinSize1DFineMesh();
130 SetNBins1DCoarseMesh();
131 SetBinSize1DCoarseMesh();
132 SetNBins3DFineMesh();
133 SetBinSize3DFineMesh();
134 SetNBins3DCoarseMesh();
135 SetBinSize3DCoarseMesh();
136 SetNBins3DFineProjectMesh();
137}
138
139/*******************************************************************/
140
141
142/*******************************************************************/
143
144AliGenHBTprocessor::~AliGenHBTprocessor()
145{
146//destructor
147 CleanStatusCodes();
148 if (fHBTprocessor) delete fHBTprocessor; //delete generator
149
150}
151
152/*******************************************************************/
153
154void AliGenHBTprocessor::InitStatusCodes()
155{
156 //creates and inits status codes array to zero
157 AliGenCocktailAfterBurner *cab = GetGenerator();
158
159 if(!cab) Fatal("AliGenHBTprocessor::InitStatusCodes()","Can not find AliGenCocktailAfterBurner generator");
160
161 Int_t nev = cab->GetNumberOfEvents();
162
163 fHbtPStatCodes = new Int_t* [nev];
164 for( Int_t i =0; i<nev;i++)
165 {
166 Int_t nprim = cab->GetStack(i)->GetNprimary();
167 fHbtPStatCodes[i] = new Int_t[nprim];
168 for (Int_t k =0 ;k<nprim;k++)
169 fHbtPStatCodes[i][k] =0;
170
171 }
172
173}
174/*******************************************************************/
175
176void AliGenHBTprocessor::CleanStatusCodes()
177{//Cleans up status codes
178 if (fHbtPStatCodes)
179 {
180 for (Int_t i =0; i<GetGenerator()->GetNumberOfEvents(); i++)
181 delete [] fHbtPStatCodes[i];
182 delete fHbtPStatCodes;
183 fHbtPStatCodes = 0;
184 }
185
186}
187/*******************************************************************/
188
189void AliGenHBTprocessor::Init()
190 {
191 //sets up parameters in generator
192 THBTprocessor *thbtp = fHBTprocessor;
193
4c90f2a0 194
2b9786f4 195 thbtp->SetTrackRejectionFactor(fTrackRejectionFactor);
196 thbtp->SetRefControl(fReferenceControl);
197
198 if ((fPid[0] == fPid[1]) || (fPid[0] == 0) || (fPid[1] == 0))
199 {
200 if (fPid[0] == 0)
201 thbtp->SetPIDs(IdFromPDG(fPid[1]) ,0);
202 else
203 thbtp->SetPIDs(IdFromPDG(fPid[0]) ,0);
204 thbtp->SetNPIDtypes(1);
205
206 if (fSwitch_type !=1)
207 Warning("AliGenHBTprocessor::Init","\nThere is only one particle type set,\n\
208 and Switch_Type differnt then 1,\n which does not make sense.\n\
209 Setting it to 1.\n");
210
211 SetSwitchType(1);
212 }
213 else
214 {
215 thbtp->SetPIDs(IdFromPDG(fPid[0]) ,IdFromPDG(fPid[1]));
216 SetNPIDtypes(2);
217 thbtp->SetSwitchType(fSwitch_type);
218 }
219
220 thbtp->SetMaxIterations(fMaxit);
221 thbtp->SetDelChi(fDelchi);
222 thbtp->SetIRand(fIrand);
223 thbtp->SetLambda(fLambda);
224 thbtp->SetR1d(fR_1d);
225 thbtp->SetRSide(fRside);
226 thbtp->SetROut(fRout);
227 thbtp->SetRLong(fRlong);
228 thbtp->SetRPerp(fRperp);
229 thbtp->SetRParallel(fRparallel);
230 thbtp->SetR0(fR0);
231 thbtp->SetQ0(fQ0);
232 thbtp->SetSwitch1D(fSwitch_1d);
233 thbtp->SetSwitch3D(fSwitch_3d);
234 thbtp->SetSwitchType(fSwitch_type);
235 thbtp->SetSwitchCoherence(fSwitch_coherence);
236 thbtp->SetSwitchCoulomb(fSwitch_coulomb);
237 thbtp->SetSwitchFermiBose(fSwitch_fermi_bose);
238 thbtp->SetPtRange(fPtMin,fPtMax);
239 thbtp->SetPxRange(fPx_min,fPx_max);
240 thbtp->SetPyRange(fPy_min,fPy_max);
241 thbtp->SetPzRange(fPz_min,fPz_max);
4c90f2a0 242 thbtp->SetPhiRange(fPhiMin*180./TMath::Pi(),fPhiMax*180./TMath::Pi());
2b9786f4 243 thbtp->SetEtaRange(fEta_min,fEta_max);
244 thbtp->SetNPtBins(fN_pt_bins);
245 thbtp->SetNPhiBins(fN_phi_bins);
246 thbtp->SetNEtaBins(fN_eta_bins);
247 thbtp->SetNPxBins(fN_px_bins);
248 thbtp->SetNPyBins(fN_py_bins);
249 thbtp->SetNPzBins(fN_pz_bins);
250 thbtp->SetNBins1DFineMesh(fN_1d_fine);
251 thbtp->SetBinSize1DFineMesh(fBinsize_1d_fine);
252 thbtp->SetNBins1DCoarseMesh(fN_1d_coarse);
253 thbtp->SetBinSize1DCoarseMesh(fBinsize_1d_coarse);
254 thbtp->SetNBins3DFineMesh(fN_3d_fine);
255 thbtp->SetBinSize3DFineMesh(fBinsize_3d_fine);
256 thbtp->SetNBins3DCoarseMesh(fN_3d_coarse);
257 thbtp->SetBinSize3DCoarseMesh(fBinsize_3d_coarse);
258 thbtp->SetNBins3DFineProjectMesh(fN_3d_fine_project);
259
260 }
261/*******************************************************************/
262
263void AliGenHBTprocessor::Generate()
264 {
265 //starts processig
266
267 if (gAlice->GetEventsPerRun() <2)
268 {
269 Fatal("AliGenHBTprocessor::Generate()",
270 "HBT Processor needs more than 1 event to work properly,\
271 but there is only %d set", gAlice->GetEventsPerRun());
272 }
4c90f2a0 273
274 fHBTprocessor->Initialize(); //reset all fields of common blocks
275 //in case there are many HBT processors
276 //run one after one (i.e. in AliCocktailAfterBurner)
277 //between Init() called and Generate there might
278 Init(); //be different instance running - be sure that we have our settings
2b9786f4 279
280 InitStatusCodes(); //Init status codes
281
282 fHBTprocessor->GenerateEvent(); //Generates event
283
284 CleanStatusCodes(); //Clean Status codes - thet are not needed anymore
285 }
286
287/*******************************************************************/
288
289
290/*******************************************************************/
291void AliGenHBTprocessor::GetParticles(TClonesArray * particles)
292 {//practically dumm
293 if (!particles)
294 {
295 cout<<"Particles has to be initialized"<<endl;
296 return;
297 }
298 fHBTprocessor->ImportParticles(particles);
299 }
300
301/*******************************************************************/
302
303Int_t AliGenHBTprocessor::GetHbtPStatusCode(Int_t part) const
304{
305//returns the status code of the given particle in the active event
306//see SetActiveEvent in the bottom of AliGenHBTprocessor.cxx
307//and in AliCocktailAfterBurner
308 Int_t ActiveEvent = GetGenerator()->GetActiveEventNumber();
309 return fHbtPStatCodes[ActiveEvent][part];
310}
311
312/*******************************************************************/
313void AliGenHBTprocessor::SetHbtPStatusCode(Int_t hbtstatcode, Int_t part)
314{
315 //Sets the given status code to given particle number (part) in the active event
316 Int_t ActiveEvent = GetGenerator()->GetActiveEventNumber();
317 fHbtPStatCodes[ActiveEvent][part] = hbtstatcode;
318}
319
320/*******************************************************************/
321
322void AliGenHBTprocessor::SetTrackRejectionFactor(Float_t trf) //def 1.0
323 {
324 //Sets the Track Rejection Factor
325 //variates in range 0.0 <-> 1.0
326 //Describes the factor of particles rejected from the output.
327 //Used only in case of low muliplicity particles e.g. lambdas.
328 //Processor generates addisional particles and builds the
329 //correletions on such a statistics.
330 //At the end these particels are left in the event according
331 //to this factor: 1==all particles are left
332 // 0==all are removed
333
334 fTrackRejectionFactor=trf;
335 fHBTprocessor->SetTrackRejectionFactor(trf);
336 }
337/*******************************************************************/
338
339void AliGenHBTprocessor::SetRefControl(Int_t rc) //default 2
340 {
341 //Sets the Refernce Control Switch
342 //switch wether read reference histograms from file =1
343 // compute from input events =2 - default
344 fReferenceControl=rc;
345 fHBTprocessor->SetRefControl(rc);
346 }
347/*******************************************************************/
348
349void AliGenHBTprocessor::SetPIDs(Int_t pid1,Int_t pid2)
350 {
351 //default pi+ pi-
352 //Sets PDG Codes of particles to be processed, default \\Pi^{+} and \\Pi^{-}
353 //This method accepts PDG codes which is
354 //in opposite to THBBProcessor which accepts GEANT PID
355 if ( (pid1 == 0) && (pid2 == 0) )
356 {
357 Error("AliGenHBTprocessor::SetPIDs","Sensless Particle Types setting: 0 0, Ignoring\n");
358 }
359
360 fPid[0]=pid1;
361 fPid[1]=pid2;
362
363 if(pid1 == pid2)
364 {
365 fHBTprocessor->SetPIDs(IdFromPDG(pid1) ,0);
4c90f2a0 366 SetNPIDtypes(1);
367 SetSwitchType(1);
2b9786f4 368 }
369 else
370 {
371 fHBTprocessor->SetPIDs(IdFromPDG(pid1) ,IdFromPDG(pid2));
4c90f2a0 372 SetNPIDtypes(2);
2b9786f4 373 }
374 }
375/*******************************************************************/
376
377void AliGenHBTprocessor::SetNPIDtypes(Int_t npidt)
378 {
379 //Number ofparticle types to be processed - default 2
380 //see AliGenHBTprocessor::SetPIDs
381
382 fNPidTypes = npidt;
383 fHBTprocessor->SetNPIDtypes(npidt);
384 }
385/*******************************************************************/
386
387void AliGenHBTprocessor::SetDeltap(Float_t deltp)
388 {
389 //default = 0.1 GeV
390 //maximum range for random momentum shifts in GeV/c;
391 //px,py,pz independent; Default = 0.1 GeV/c.
392 fDeltap=deltp;
393 fHBTprocessor->SetDeltap(deltp);
394 }
395/*******************************************************************/
396
397void AliGenHBTprocessor::SetMaxIterations(Int_t maxiter)
398 {
399 //maximum # allowed iterations thru all the
400 //tracks for each event. Default = 50.
401 //If maxit=0, then calculate the correlations
402 //for the input set of events without doing the
403 //track adjustment procedure.
404
405 fMaxit=maxiter;
406 fHBTprocessor->SetMaxIterations(maxiter);
407 }
408
409/*******************************************************************/
410void AliGenHBTprocessor::SetDelChi(Float_t dc)
411 {
412 //min % change in total chi-square for which
413 //the track adjustment iterations may stop,
414 //Default = 0.1%.
415
416 fDelchi=dc;
417 fHBTprocessor->SetDelChi(dc);
418 }
419/*******************************************************************/
420
421void AliGenHBTprocessor::SetIRand(Int_t irnd)
422 {
423 //if fact dummy - used only for compatibility
424 //we are using AliRoot built in RNG
425 //dummy in fact since we are using aliroot build-in RNG
426 //Sets renaodom number generator
427 fIrand=irnd;
428 fHBTprocessor->SetIRand(irnd);
429 }
430/*******************************************************************/
431
432void AliGenHBTprocessor::SetLambda(Float_t lam)
433 {
434 //default = 0.6
435 // Sets Chaoticity Parameter
436 fLambda = lam;
437 fHBTprocessor->SetLambda(lam);
438 }
439/*******************************************************************/
440void AliGenHBTprocessor::SetR1d(Float_t r)
441 {
442 //default 7.0
443 //Sets Spherical source model radius (fm)
444 fR_1d = r;
445 fHBTprocessor->SetR1d(r);
446 }
447/*******************************************************************/
448void AliGenHBTprocessor::SetRSide(Float_t rs)
449 {
450 //default rs = 6.0
451 //Rside,Rout,Rlong = Non-spherical Bertsch-Pratt source model (fm)
452
453 fRside = rs;
454 fHBTprocessor->SetRSide(rs);
455 }
456/*******************************************************************/
457void AliGenHBTprocessor::SetROut(Float_t ro)
458 {
459 //default ro = 7.0
460 //Rside,Rout,Rlong = Non-spherical Bertsch-Pratt source model (fm)
461 fRout = ro;
462 fHBTprocessor->SetROut(ro);
463 }
464/*******************************************************************/
465void AliGenHBTprocessor::SetRLong(Float_t rl)
466 {
467 //default rl = 4.0
468 //Rside,Rout,Rlong = Non-spherical Bertsch-Pratt source model (fm)
469 fRlong = rl;
470 fHBTprocessor->SetRLong(rl);
471 }
472/*******************************************************************/
473void AliGenHBTprocessor::SetRPerp(Float_t rp)
474 {
475 //default rp = 6.0
476 //Rperp,Rparallel,R0= Non-spherical Yano-Koonin-Podgoretski source model (fm).
477 fRperp = rp;
478 fHBTprocessor->SetRPerp(rp);
479 }
480/*******************************************************************/
481void AliGenHBTprocessor::SetRParallel(Float_t rprl)
482 {
483 //default rprl = 4.0
484 //Rperp,Rparallel,R0= Non-spherical Yano-Koonin-Podgoretski source model (fm).
485 fRparallel = rprl;
486 fHBTprocessor->SetRParallel(rprl);
487 }
488/*******************************************************************/
489void AliGenHBTprocessor::SetR0(Float_t r0)
490 {
491 //default r0 = 4.0
492 //Rperp,Rparallel,R0= Non-spherical Yano-Koonin-Podgoretski source model (fm).
493 fR0 = r0;
494 fHBTprocessor->SetR0(r0);
495 }
496/*******************************************************************/
497void AliGenHBTprocessor::SetQ0(Float_t q0)
498 {
499 //default q0 = 9.0
500 //Sets Q0 = NA35 Coulomb parameter for finite source size in (GeV/c)
501 // if fSwitchCoulomb = 2
502 // = Spherical Coulomb source radius in (fm)
503 // if switch_coulomb = 3, used to interpolate the
504 // input Pratt/Cramer discrete Coulomb source
505 // radii tables.
506 fQ0 = q0;
507 fHBTprocessor->SetQ0(q0);
508 }
509
510/*******************************************************************/
511void AliGenHBTprocessor::SetSwitch1D(Int_t s1d)
512 {
513//default s1d = 3
4c90f2a0 514// Sets fSwitch_1d
515// = 0 to not compute the 1D two-body //orrelations.
516// = 1 to compute this using Q-invariant
517// = 2 to compute this using Q-total
518// = 3 to compute this using Q-3-ve//tor
2b9786f4 519
520 fSwitch_1d = s1d;
521 fHBTprocessor->SetSwitch1D(s1d);
522 }
523/*******************************************************************/
524void AliGenHBTprocessor::SetSwitch3D(Int_t s3d)
525 {
526//default s3d = 0
4c90f2a0 527// Sets fSwitch_3d
528// = 0 to not compute the 3D two-body correlations.
529// = 1 to compute this using the side-out-long form
530// = 2 to compute this using the Yanno-Koonin-Pogoredskij form.
2b9786f4 531
532 fSwitch_3d = s3d;
533 fHBTprocessor->SetSwitch3D(s3d);
534 }
535/*******************************************************************/
536void AliGenHBTprocessor::SetSwitchType(Int_t st)
537 {
538//default st = 3
539// Sets switch_type = 1 to fit only the like pair correlations
540// = 2 to fit only the unlike pair correlations
541// = 3 to fit both the like and unlike pair correl.
542//See SetPIDs and Init
543//If only one particle type is set, unly==1 makes sens
544
545 fSwitch_type = st;
546 fHBTprocessor->SetSwitchType(st);
547 }
548/*******************************************************************/
549void AliGenHBTprocessor::SetSwitchCoherence(Int_t sc)
550 {
551// default sc = 0
552// switch_coherence = 0 for purely incoherent source (but can have
553// lambda < 1.0)
554// = 1 for mixed incoherent and coherent source
555
556 fSwitch_coherence = sc;
557 fHBTprocessor->SetSwitchCoherence(sc);
558 }
559/*******************************************************************/
560void AliGenHBTprocessor::SetSwitchCoulomb(Int_t scol)
561 {
562//default scol = 2
563// switch_coulomb = 0 no Coulomb correction
564// = 1 Point source Gamow correction only
565// = 2 NA35 finite source size correction
566// = 3 Pratt/Cramer finite source size correction;
567// interpolated from input tables.
568 fSwitch_coulomb =scol;
569 fHBTprocessor->SetSwitchCoulomb(scol);
570 }
571/*******************************************************************/
572void AliGenHBTprocessor::SetSwitchFermiBose(Int_t sfb)
573 {
574//default sfb = 1
575// switch_fermi_bose = 1 Boson pairs
576// = -1 Fermion pairs
577
578 fSwitch_fermi_bose = sfb;
579 fHBTprocessor->SetSwitchFermiBose(sfb);
580 }
581/*******************************************************************/
582void AliGenHBTprocessor::SetPtRange(Float_t ptmin, Float_t ptmax)
583 {
584// default ptmin = 0.1, ptmax = 0.98
585//Sets Pt range (GeV)
586 AliGenerator::SetPtRange(ptmin,ptmax);
587 fHBTprocessor->SetPtRange(ptmin,ptmax);
588 }
589
590/*******************************************************************/
591void AliGenHBTprocessor::SetPxRange(Float_t pxmin, Float_t pxmax)
592 {
593//default pxmin = -1.0, pxmax = 1.0
594//Sets Px range
595 fPx_min =pxmin;
596 fPx_max =pxmax;
597 fHBTprocessor->SetPxRange(pxmin,pxmax);
598 }
599/*******************************************************************/
600void AliGenHBTprocessor::SetPyRange(Float_t pymin, Float_t pymax)
601 {
602//default pymin = -1.0, pymax = 1.0
603//Sets Py range
604 fPy_min =pymin;
605 fPy_max =pymax;
606 fHBTprocessor->SetPyRange(pymin,pymax);
607 }
608/*******************************************************************/
609void AliGenHBTprocessor::SetPzRange(Float_t pzmin, Float_t pzmax)
610 {
611//default pzmin = -3.6, pzmax = 3.6
612//Sets Py range
613 fPz_min =pzmin;
614 fPz_max =pzmax;
615 fHBTprocessor->SetPzRange(pzmin,pzmax);
616 }
617void AliGenHBTprocessor::SetMomentumRange(Float_t pmin, Float_t pmax)
618 {
619 //default pmin=0, pmax=0
620 //Do not use this method!
621 MayNotUse("AliGenHBTprocessor::SetMomentumRange Method is Dummy");
622 }
623
624 /*******************************************************************/
625void AliGenHBTprocessor::SetPhiRange(Float_t phimin, Float_t phimax)
626 {
627//
628//Sets \\Phi range
629 AliGenerator::SetPhiRange(phimin,phimax);
630
631 fHBTprocessor->SetPhiRange(phimin,phimax);
632 }
633/*******************************************************************/
634void AliGenHBTprocessor::SetEtaRange(Float_t etamin, Float_t etamax)
635 {
636//default etamin = -1.5, etamax = 1.5
637//Sets \\Eta range
638 fEta_min= etamin;
639 fEta_max =etamax;
640 fHBTprocessor->SetEtaRange(etamin,etamax);
4c90f2a0 641
642 //set the azimothal angle range in the AliGeneraor -
643 //to keep coherency between azimuthal angle and pseudorapidity
644 //DO NOT CALL this->SetThetaRange, because it calls this method (where we are)
645 //which must cause INFINITE LOOP
646 AliGenerator::SetThetaRange(RadiansToDegrees(EtaToTheta(fEta_min)),
647 RadiansToDegrees(EtaToTheta(fEta_max)));
648
2b9786f4 649 }
4c90f2a0 650/*******************************************************************/
651void AliGenHBTprocessor::SetThetaRange(Float_t thetamin, Float_t thetamax)
652{
653 //default thetamin = 0, thetamax = 180
654 //Azimuthal angle, override AliGenerator method which uses widely (i.e. wrapper generators)
655 //core fortran HBTProcessor uses Eta (pseudorapidity)
656 //so these methods has to be synchronized
657
658 AliGenerator::SetThetaRange(thetamin,thetamax);
659
660 SetEtaRange( ThetaToEta(fThetaMin) , ThetaToEta(fThetaMax) );
661
662}
663
2b9786f4 664/*******************************************************************/
665void AliGenHBTprocessor::SetNPtBins(Int_t nptbin)
666 {
667 //default nptbin = 50
668 //set number of Pt bins
669 fN_pt_bins= nptbin;
670 fHBTprocessor->SetNPtBins(nptbin);
671 }
672/*******************************************************************/
673void AliGenHBTprocessor::SetNPhiBins(Int_t nphibin)
674{
675 //default nphibin = 50
676 //set number of Phi bins
677 fN_phi_bins=nphibin;
678 fHBTprocessor->SetNPhiBins(nphibin);
679}
680/*******************************************************************/
681void AliGenHBTprocessor::SetNEtaBins(Int_t netabin)
682{
683 //default netabin = 50
684 //set number of Eta bins
685 fN_eta_bins = netabin;
686 fHBTprocessor->SetNEtaBins(netabin);
687}
688/*******************************************************************/
689void AliGenHBTprocessor::SetNPxBins(Int_t npxbin)
690{
691 //default npxbin = 20
692 //set number of Px bins
693 fN_px_bins = npxbin;
694 fHBTprocessor->SetNPxBins(npxbin);
695}
696/*******************************************************************/
697void AliGenHBTprocessor::SetNPyBins(Int_t npybin)
698{
699 //default npybin = 20
700 //set number of Py bins
701 fN_py_bins = npybin;
702 fHBTprocessor->SetNPyBins(npybin);
703}
704/*******************************************************************/
705void AliGenHBTprocessor::SetNPzBins(Int_t npzbin)
706{
707 //default npzbin = 70
708 //set number of Pz bins
709 fN_pz_bins = npzbin;
710 fHBTprocessor->SetNPzBins(npzbin);
711}
712/*******************************************************************/
713void AliGenHBTprocessor::SetNBins1DFineMesh(Int_t n)
714{
715//default n = 10
716//Sets the number of bins in the 1D mesh
717 fN_1d_fine =n;
718 fHBTprocessor->SetNBins1DFineMesh(n);
719
720}
721/*******************************************************************/
722void AliGenHBTprocessor::SetBinSize1DFineMesh(Float_t x)
723{
724//default x=0.01
725//Sets the bin size in the 1D mesh
726 fBinsize_1d_fine = x;
727 fHBTprocessor->SetBinSize1DFineMesh(x);
728}
729/*******************************************************************/
730
731void AliGenHBTprocessor::SetNBins1DCoarseMesh(Int_t n)
732{
733//default n =2
734//Sets the number of bins in the coarse 1D mesh
735 fN_1d_coarse =n;
736 fHBTprocessor->SetNBins1DCoarseMesh(n);
737}
738/*******************************************************************/
739void AliGenHBTprocessor::SetBinSize1DCoarseMesh(Float_t x)
740{
741//default x=0.05
742//Sets the bin size in the coarse 1D mesh
743 fBinsize_1d_coarse =x;
744 fHBTprocessor->SetBinSize1DCoarseMesh(x);
745}
746/*******************************************************************/
747
748void AliGenHBTprocessor::SetNBins3DFineMesh(Int_t n)
749{
750//default n = 8
751//Sets the number of bins in the 3D mesh
752 fN_3d_fine =n;
753 fHBTprocessor->SetNBins3DFineMesh(n);
754}
755/*******************************************************************/
756void AliGenHBTprocessor::SetBinSize3DFineMesh(Float_t x)
757{
758//default x=0.01
759//Sets the bin size in the 3D mesh
760 fBinsize_3d_fine =x;
761 fHBTprocessor->SetBinSize3DFineMesh(x);
762}
763/*******************************************************************/
764
765void AliGenHBTprocessor::SetNBins3DCoarseMesh(Int_t n)
766{
767//default n = 2
768//Sets the number of bins in the coarse 3D mesh
769
770 fN_3d_coarse = n;
771 fHBTprocessor->SetNBins3DCoarseMesh(n);
772}
773/*******************************************************************/
774void AliGenHBTprocessor::SetBinSize3DCoarseMesh(Float_t x)
775{
776//default x=0.08
777//Sets the bin size in the coarse 3D mesh
778 fBinsize_3d_coarse = x;
779 fHBTprocessor->SetBinSize3DCoarseMesh(x);
780}
781/*******************************************************************/
782
783void AliGenHBTprocessor::SetNBins3DFineProjectMesh(Int_t n )
784{
785//default n =3
786//Sets the number of bins in the fine project mesh
787 fN_3d_fine_project = n;
788 fHBTprocessor->SetNBins3DFineProjectMesh(n);
789}
790/*******************************************************************/
791
792
793/*******************************************************************/
794
795
796
797
798
799
800void AliGenHBTprocessor::DefineParticles()
801{
802 //
803 // Load standard numbers for GEANT particles and PDG conversion
804 fNPDGCodes = 0; //this is done in the constructor - but in any case ...
805
806 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
807 fPDGCode[fNPDGCodes++]=22; // 1 = photon
808 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
809 fPDGCode[fNPDGCodes++]=11; // 3 = electron
810 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
811 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
812 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
813 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
814 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
815 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
816 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
817 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
818 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
819 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
820 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
821 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
822 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
823 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
824 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
825 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
826 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
827 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
828 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
829 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
830 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
831 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
832 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
833 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
834 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
835 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
836 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
837 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
838 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
839}
840
841/*******************************************************************/
842Int_t AliGenHBTprocessor::IdFromPDG(Int_t pdg) const
843{
844 //
845 // Return Geant3 code from PDG and pseudo ENDF code
846 //
847 for(Int_t i=0;i<fNPDGCodes;++i)
848 if(pdg==fPDGCode[i]) return i;
849 return -1;
850}
851Int_t AliGenHBTprocessor::PDGFromId(Int_t id) const
852{
853 //
854 // Return PDG code and pseudo ENDF code from Geant3 code
855 //
856 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
857 else return -1;
4c90f2a0 858}
859Double_t AliGenHBTprocessor::ThetaToEta(Double_t arg)
860 {
861 //converts etha(azimuthal angle) to Eta (pseudorapidity). Argument in radians
862
863 if(arg>= TMath::Pi()) return 709.0;//This number is the biggest wich not crashes exp(x)
864 if(arg<= 0.0) return -709.0;//
865
866 arg -= TMath::Pi()/2.;
867 if (arg > 0.0)
868 {
869 return -TMath::Log( TMath::Tan(arg/2.)) ;
870 }
871 else
872 {
873 return TMath::Log( TMath::Tan(-arg/2.)) ;
874 }
875 }
876
2b9786f4 877/*******************************************************************/
878/****** ROUTINES USED FOR COMMUNUCATION ********/
879/******************** WITH FORTRAN ********************/
880/*******************************************************************/
881
882#ifndef WIN32
883 # define hbtpran hbtpran_
884 # define alihbtp_puttrack alihbtp_puttrack_
885 # define alihbtp_gettrack alihbtp_gettrack_
886 # define alihbtp_getnumberevents alihbtp_getnumberevents_
887 # define alihbtp_getnumbertracks alihbtp_getnumbertracks_
888 # define alihbtp_initialize alihbtp_initialize_
889 # define alihbtp_setactiveeventnumber alihbtp_setactiveeventnumber_
890 # define alihbtp_setparameters alihbtp_setparameters_
891 # define type_of_call
892
893#else
894 # define hbtpran HBTPRAN
895 # define alihbtp_puttrack ALIHBTP_PUTTRACK
896 # define alihbtp_gettrack ALIHBTP_GETTRACK
897 # define alihbtp_getnumberevents ALIHBTP_GETNUMBEREVENTS
898 # define alihbtp_getnumbertracks ALIHBTP_GETNUMBERTRACKS
899 # define alihbtp_initialize ALIHBTP_INITIALIZE
900 # define alihbtp_setactiveeventnumber ALIHBTP_SETACTIVEEVENTNUMBER
901 # define alihbtp_setparameters ALIHBTP_SETPARAMETERS
902 # define type_of_call _stdcall
903#endif
904
905#include "AliGenCocktailAfterBurner.h"
906#include <string.h>
907/*******************************************************************/
908
909AliGenCocktailAfterBurner* GetGenerator()
910 {
911 // This function has two tasks:
912 // Check if environment is OK (exist gAlice and generator)
913 // Returns pointer to genrator
914 //to be changed with TFolders
915
916 if(!gAlice)
917 {
918 cout<<endl<<"ERROR: There is NO gALICE! Check what you are doing!"<<endl;
919 gAlice->Fatal("AliGenHBTprocessor.cxx: alihbtp_getnumofev()",
920 "\nRunning HBT Processor without gAlice... Exiting \n");
921 return 0;
922 }
923 AliGenerator * gen = gAlice->Generator();
924
925 if (!gen)
926 {
927 gAlice->Fatal("AliGenHBTprocessor.cxx: GetGenerator()",
928 "\nThere is no generator in gAlice, exiting\n");
929 return 0;
930 }
931
932 const Char_t *genname = gen->GetName();
933 Char_t name[30];
934 strcpy(name,"AliGenCocktailAfterBurner");
935
936 if (strcmp(name,genname))
937 {
938 gAlice->Fatal("AliGenHBTprocessor.cxx: GetGenerator()",
939 "\nThe main Generator is not a AliGenCocktailAfterBurner, exiting\n");
940 return 0;
941 }
942
943 AliGenCocktailAfterBurner* CAB= (AliGenCocktailAfterBurner*) gen;
944
945 // cout<<endl<<"Got generator"<<endl;
946 return CAB;
947
948 }
949/*******************************************************************/
950
951AliGenHBTprocessor* GetAliGenHBTprocessor()
952{
953//returns pointer to the current instance of AliGenHBTprocessor in
954//AliGenCocktailAfterBurner (can be more than one)
955//
956 AliGenCocktailAfterBurner* gen = GetGenerator();
957 AliGenerator * g = gen->GetCurrentGenerator();
958
959 const Char_t *genname = g->GetName();
960 Char_t name[30];
961 strcpy(name,"AliGenHBTprocessor");
962
963 if (strcmp(name,genname))
964 {
965 gAlice->Fatal("AliGenHBTprocessor.cxx: GetAliGenHBTprocessor()",
966 "\nCurrernt generator in AliGenCocktailAfterBurner is not a AliGenHBTprocessor, exiting\n");
967 return 0;
968 }
969
970 AliGenHBTprocessor* hbtp = (AliGenHBTprocessor*)g;
971 return hbtp;
972}
973
974/*******************************************************************/
975extern "C" void type_of_call alihbtp_setparameters()
976 {
977 //dummy
978 }
979
980extern "C" void type_of_call alihbtp_initialize()
981 {
982 //dummy
983 }
984
985/*******************************************************************/
986
987extern "C" void type_of_call alihbtp_getnumberevents(Int_t &nev)
988 {
989 //passes number of events to the fortran
990 if(gDebug) cout<<"alihbtp_getnumberevents("<<nev<<") ....";
991 AliGenCocktailAfterBurner* gen = GetGenerator();
992 if(!gen)
993 {
994 nev = -1;
995 return;
996 }
997
998 nev = gen->GetNumberOfEvents();
999
1000 if(gDebug>5) cout<<"EXITED N Ev = "<<nev<<endl;
1001
1002 }
1003
1004/*******************************************************************/
1005
1006extern "C" void type_of_call alihbtp_setactiveeventnumber(Int_t & nev)
1007 {
1008//sets active event in generator (AliGenCocktailAfterBurner)
1009
1010 if(gDebug>5) cout<<"alihbtp_setactiveeventnumber("<<nev<<") ....";
1011 if(gDebug>0) cout<<"Asked for event "<<nev-1<<endl;
1012 AliGenCocktailAfterBurner* gen = GetGenerator();
1013 if(!gen) return;
1014 gen->SetActiveEventNumber(nev - 1); //fortran numerates events from 1 to N
1015
1016 if(gDebug>5) cout<<"EXITED returned "<<nev<<endl;
1017
1018 }
1019/*******************************************************************/
1020
1021extern "C" void type_of_call alihbtp_getnumbertracks(Int_t &ntracks)
1022 {
1023//passes number of particles in active event to the fortran
1024 if(gDebug>5) cout<<"alihbtp_getnumbertracks("<<ntracks<<") ....";
1025
1026 AliGenCocktailAfterBurner* gen = GetGenerator();
1027 if (!gen)
1028 {
1029 ntracks = -1;
1030 return;
1031 }
1032
1033 ntracks = gen->GetActiveStack()->GetNprimary();
1034 if(gDebug>5) cout<<"EXITED Ntracks = "<<ntracks<<endl;
1035 }
1036
1037/*******************************************************************/
1038
1039extern "C" void type_of_call
1040 alihbtp_puttrack(Int_t & n,Int_t& flag, Float_t& px,
1041 Float_t& py, Float_t& pz, Int_t& geantpid)
1042 {
1043//sets new parameters (momenta) in track number n
1044// in the active event
1045// n - number of the track in active event
1046// flag - flag of the track
1047// px,py,pz - momenta
1048// geantpid - type of the particle - Geant Particle ID
1049
1050 if(gDebug>5) cout<<"alihbtp_puttrack("<<n<<") ....";
1051
1052 AliGenCocktailAfterBurner* gen = GetGenerator();
1053 if(!gen) return;
1054
1055 TParticle * track = gen->GetActiveStack()->Particle(n-1);
1056
1057 AliGenHBTprocessor* g = GetAliGenHBTprocessor();
1058
1059 //check to be deleted
1060 if (geantpid != (g->IdFromPDG( track->GetPdgCode() )))
1061 {
1062 cout<<endl<<" AliGenHBTprocessor.cxx: alihbtp_puttrack: SOMETHING IS GOING BAD:\n GEANTPIDS ARE NOT THE SAME"<<endl;
1063 }
1064
1065 if(gDebug>0)
1066 if (px != track->Px())
1067 cout<<"Px diff. = "<<px - track->Px()<<endl;
1068
1069 if(gDebug>3) cout<<" track->GetPdgCode() --> "<<track->GetPdgCode()<<endl;
1070
1071
1072
5925692f 1073 Float_t m = track->GetMass();
1074 Float_t E = TMath::Sqrt(m*m+px*px+py*py+pz*pz);
1075 track->SetMomentum(px,py,pz,E);
2b9786f4 1076
1077 g->SetHbtPStatusCode(flag,n-1);
1078
1079 if(gDebug>5) cout<<"EXITED "<<endl;
1080 }
1081
1082/*******************************************************************/
1083
1084extern "C" void type_of_call
1085 alihbtp_gettrack(Int_t & n,Int_t & flag, Float_t & px,
1086 Float_t & py, Float_t & pz, Int_t & geantpid)
1087
1088 {
1089//passes track parameters to the fortran
1090// n - number of the track in active event
1091// flag - flag of the track
1092// px,py,pz - momenta
1093// geantpid - type of the particle - Geant Particle ID
1094
1095 if(gDebug>5) cout<<"alihbtp_gettrack("<<n<<") ....";
1096 AliGenCocktailAfterBurner* gen = GetGenerator();
1097
1098 if (!gen)
1099 {
1100 n = -1;
1101 flag =-1;
1102 px = py = pz = -1;
1103 geantpid = -1;
1104 return;
1105 }
1106
1107 TParticle * track = gen->GetActiveStack()->Particle(n-1);
1108 AliGenHBTprocessor* g = GetAliGenHBTprocessor();
1109
1110 flag = g->GetHbtPStatusCode(n-1);
1111
1112 px = (Float_t)track->Px();
1113 py = (Float_t)track->Py();
1114 pz = (Float_t)track->Pz();
1115
1116 geantpid = g->IdFromPDG( track->GetPdgCode() );
1117
1118 if(gDebug>5) cout<<"EXITED "<<endl;
1119 }
1120
1121/*******************************************************************/
1122extern "C" Float_t type_of_call hbtpran(Int_t &)
1123{
1124//interface to the random number generator
1125 return sRandom->Rndm();
1126}
1127
1128/*******************************************************************/
1129
1130
1131/*******************************************************************/