]>
Commit | Line | Data |
---|---|---|
ff6aba1e | 1 | // Implementation of the interface for THBTprocessor |
2 | // which is a wrapper itself to Fortran | |
3 | // program "HBT processor" written by Lanny Ray | |
4 | // Author: Piotr Krzysztof Skowronski <Piotr.Skowronski@cern.ch> | |
5 | // | |
2b9786f4 | 6 | /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * |
7 | * See cxx source for full Copyright notice */ | |
8 | ||
9 | /* $Id$ */ | |
10 | ||
ff6aba1e | 11 | #ifndef ALIGENHBTPROCESSOR_H |
12 | #define ALIGENHBTPROCESSOR_H | |
2b9786f4 | 13 | |
14 | #include "AliGenerator.h" | |
2b9786f4 | 15 | #include <AliPDG.h> |
20dddfab | 16 | |
17 | class THBTprocessor; | |
18 | class TClonesArray; | |
2b9786f4 | 19 | |
20 | enum {kHBTPMaxParticleTypes = 50}; | |
21 | ||
20dddfab | 22 | class AliGenHBTprocessor : public AliGenerator |
23 | { | |
24 | //Wrapper class for THBTProcessor | |
ff6aba1e | 25 | //which is a wrapper itself to Fortran |
26 | //program "HBT processor" written by Lanny Ray | |
20dddfab | 27 | // |
28 | //Piotr.Skowronski@cern.ch | |
2b9786f4 | 29 | |
30 | public: | |
31 | AliGenHBTprocessor(); | |
20dddfab | 32 | AliGenHBTprocessor(const AliGenHBTprocessor& in); |
2b9786f4 | 33 | virtual ~AliGenHBTprocessor(); |
34 | ||
35 | virtual void Init(); | |
36 | virtual void Generate(); | |
ff6aba1e | 37 | virtual void GetParticles(TClonesArray * particles) const; |
20dddfab | 38 | Int_t IdFromPDG(Int_t pdg) const; |
39 | Int_t PDGFromId(Int_t id) const; | |
2b9786f4 | 40 | |
20dddfab | 41 | Int_t GetHbtPStatusCode(Int_t part) const; |
42 | void SetHbtPStatusCode(Int_t hbtstatcode, Int_t part); | |
2b9786f4 | 43 | /************* S E T T E R S ******************/ |
44 | ||
45 | virtual void SetTrackRejectionFactor(Float_t trf = 1.0); | |
46 | ||
47 | virtual void SetRefControl(Int_t rc =2); | |
48 | virtual void SetPIDs(Int_t pid1 = kPiPlus,Int_t pid2 = kPiMinus); //PDG Codes of particles to be processed, default \\Pi^{+} and \\Pi^{-} | |
49 | virtual void SetNPIDtypes(Int_t npidt = 2); //Number ofparticle types to be processed | |
50 | virtual void SetDeltap(Float_t deltp = 0.1); //maximum range for random momentum shifts in GeV/c; | |
51 | //px,py,pz independent; Default = 0.1 GeV/c. | |
20dddfab | 52 | virtual void SetMaxIterations(Int_t maxiter = 50);// |
2b9786f4 | 53 | virtual void SetDelChi(Float_t dc = 0.1); |
54 | virtual void SetIRand(Int_t irnd = 76564) ; | |
55 | ||
56 | virtual void SetLambda(Float_t lam = 0.6); | |
57 | virtual void SetR1d(Float_t r = 7.0) ; | |
58 | virtual void SetRSide(Float_t rs = 6.0); | |
59 | virtual void SetROut(Float_t ro = 7.0) ; | |
60 | virtual void SetRLong(Float_t rl = 4.0) ; | |
61 | virtual void SetRPerp(Float_t rp = 6.0); | |
62 | virtual void SetRParallel(Float_t rprl = 4.0); | |
63 | virtual void SetR0(Float_t r0 = 4.0) ; | |
64 | virtual void SetQ0(Float_t q0 = 9.0) ; | |
65 | virtual void SetSwitch1D(Int_t s1d = 3); | |
66 | virtual void SetSwitch3D(Int_t s3d = 0) ; | |
67 | virtual void SetSwitchType(Int_t st = 3); | |
68 | virtual void SetSwitchCoherence(Int_t sc = 0); | |
69 | virtual void SetSwitchCoulomb(Int_t scol = 2); | |
70 | virtual void SetSwitchFermiBose(Int_t sfb = 1); | |
71 | ||
72 | virtual void SetMomentumRange(Float_t pmin=0, Float_t pmax=0); //Dummy method | |
73 | virtual void SetPtRange(Float_t ptmin = 0.1, Float_t ptmax = 0.98); | |
74 | virtual void SetPxRange(Float_t pxmin = -1.0, Float_t pxmax = 1.0); | |
75 | virtual void SetPyRange(Float_t pymin = -1.0, Float_t pymax = 1.0); | |
76 | virtual void SetPzRange(Float_t pzmin = -3.6, Float_t pzmax = 3.6); | |
77 | ||
78 | virtual void SetPhiRange(Float_t phimin = 0.0, Float_t phimax = 360.0);//Angle in degrees | |
79 | //coherent with AliGenCocktail | |
80 | //incohernet with AliGenerator | |
81 | virtual void SetEtaRange(Float_t etamin = -1.5, Float_t etamax = 1.5);//Pseudorapidity | |
4c90f2a0 | 82 | void SetThetaRange(Float_t thetamin = 0, Float_t thetamax = 180); //Azimuthal angle, override AliGenerator method |
83 | //which uses this, core fortran HBTProcessor uses Eta (pseudorapidity) | |
84 | //so these methods has to be synchronized | |
85 | ||
2b9786f4 | 86 | virtual void SetNPtBins(Int_t nptbin = 50); |
87 | virtual void SetNPhiBins(Int_t nphibin = 50); | |
88 | virtual void SetNEtaBins(Int_t netabin = 50); | |
89 | virtual void SetNPxBins(Int_t npxbin = 20); | |
90 | virtual void SetNPyBins(Int_t npybin = 20); | |
91 | virtual void SetNPzBins(Int_t npzbin = 70); | |
92 | ||
93 | ||
94 | virtual void SetNBins1DFineMesh(Int_t n = 10); | |
95 | virtual void SetBinSize1DFineMesh(Float_t x=0.01); | |
96 | ||
97 | virtual void SetNBins1DCoarseMesh(Int_t n =2 ); | |
98 | virtual void SetBinSize1DCoarseMesh(Float_t x=0.05); | |
99 | ||
100 | virtual void SetNBins3DFineMesh(Int_t n = 8); | |
101 | virtual void SetBinSize3DFineMesh(Float_t x=0.01); | |
102 | ||
103 | virtual void SetNBins3DCoarseMesh(Int_t n = 2); | |
104 | virtual void SetBinSize3DCoarseMesh(Float_t x=0.08); | |
105 | ||
106 | virtual void SetNBins3DFineProjectMesh(Int_t n =3 ); | |
107 | /***********************************************************************/ | |
108 | /* * * * * * * P R O T E C T E D A R E A * * * * * * * * * * * */ | |
109 | /***********************************************************************/ | |
110 | protected: | |
111 | ||
112 | THBTprocessor * fHBTprocessor; //pointer to generator (TGenerator) | |
113 | Int_t **fHbtPStatCodes; //! hbtp status codes of particles | |
114 | Int_t fNPDGCodes; //! Number of defined particles | |
115 | Int_t fPDGCode[kHBTPMaxParticleTypes]; //! PDG codes (for conversion PDG<->Geant) | |
116 | void DefineParticles(); //initiates array with PDG codes | |
117 | void InitStatusCodes(); //Initiates status codes (allocates memory and sets everything to zero) | |
118 | void CleanStatusCodes(); //deletes array with status codes | |
119 | /********** P A R A M E T E R S OF THE GENERATOR****************/ | |
120 | ||
121 | Float_t fTrackRejectionFactor; //variates in range 0.0 <-> 1.0 | |
122 | //Describes the factor of particles rejected from the output. | |
123 | //Used only in case of low muliplicity particles e.g. lambdas. | |
124 | //Processor generates addisional particles and builds the | |
125 | //correletions on such a statistics. | |
126 | //At the end these particels are left in the event according | |
127 | //to this factor: 1==all particles are left | |
128 | // 0==all are removed | |
129 | Int_t fReferenceControl; //switch wether read reference histograms from file =1 | |
130 | // compute from input events =2 - default | |
131 | Int_t fPrintFull; // Full print out option - each event | |
132 | Int_t fPrintSectorData; // Print sector overflow diagnostics | |
133 | Int_t fNPidTypes; // # particle ID types to correlate | |
134 | Int_t fPid[2]; // Geant particle ID #s, max of 2 types | |
135 | Int_t fNevents ; // # events in input event text file | |
20dddfab | 136 | Int_t fSwitch1d; // Include 1D correlations |
137 | Int_t fSwitch3d; // Include 3D correlations | |
138 | Int_t fSwitchType ; // For like, unlike or both PID pairs | |
139 | Int_t fSwitchCoherence; // To include incoh/coher mixed source | |
140 | Int_t fSwitchCoulomb; // Coulomb correction selection options | |
141 | Int_t fSwitchFermiBose; // For fermions or bosons | |
2b9786f4 | 142 | |
2b9786f4 | 143 | // Counters: |
144 | ||
20dddfab | 145 | Int_t fEventLineCounter; // Input event text file line counter |
2b9786f4 | 146 | Int_t fMaxit; // Max # iterations in track adjustment |
147 | Int_t fIrand; // Random # starting seed (Def=12345) | |
2b9786f4 | 148 | // // line counter |
149 | ||
150 | // Correlation Model Parameters: | |
151 | ||
152 | Float_t fLambda; // Chaoticity parameter | |
20dddfab | 153 | Float_t fR1d; // Spherical source radius (fm) |
2b9786f4 | 154 | Float_t fRside; // 3D Bertsch-Pratt source 'side' R (fm) |
155 | Float_t fRout; // 3D Bertsch-Pratt source 'out' R (fm) | |
156 | Float_t fRlong; // 3D Bertsch-Pratt source 'long' R (fm) | |
157 | Float_t fRperp; // 3D YKP source transverse radius (fm) | |
158 | Float_t fRparallel; // 3D YKP source longitudinal radius(fm) | |
159 | Float_t fR0; // 3D YKP source emission time durat(fm) | |
160 | Float_t fQ0; // NA35 Coulomb parameter (GeV/c) or | |
161 | // // Coul radius for Pratt finite src (fm) | |
162 | ||
163 | // Search Control Parameters: | |
164 | ||
165 | ||
166 | Float_t fDeltap; // Max limit for x,y,z momt shifts(GeV/c) | |
167 | Float_t fDelchi; // Min% change in Chi-Sq to stop iterat. | |
168 | ||
169 | ||
2b9786f4 | 170 | // Particle Masses: |
171 | ||
2b9786f4 | 172 | |
173 | /********** M E S H ****************/ | |
174 | ||
175 | ||
20dddfab | 176 | Int_t fNPtBins; // # one-body pt bins |
177 | Int_t fNPhiBins; // # one-body phi bins | |
178 | Int_t fNEtaBins; // # one-body eta bins | |
2b9786f4 | 179 | |
20dddfab | 180 | Int_t fN1dFine; // # bins for 1D, Fine Mesh |
181 | Int_t fN1dCoarse; // # bins for 1D, Coarse Mesh | |
182 | Int_t fN1dTotal; // Total # bins for 1D | |
183 | Int_t fN3dFine ; // # bins for 3D, Fine Mesh | |
184 | Int_t fN3dCoarse; // # bins for 3D, Coarse Mesh | |
185 | Int_t fN3dTotal; // Total # bins for 3D | |
186 | Int_t fN3dFineProject; // # 3D fine mesh bins to sum over for | |
2b9786f4 | 187 | |
188 | // Momentum Space Sectors for Track Sorting: | |
189 | ||
20dddfab | 190 | Int_t fNPxBins; // # sector bins in px |
191 | Int_t fNPyBins; // # sector bins in py | |
192 | Int_t fNPzBins; // # sector bins in pz | |
193 | Int_t fNSectors; // Total # sectors in 3D momentum space | |
2b9786f4 | 194 | |
2b9786f4 | 195 | |
20dddfab | 196 | Float_t fPtBinSize ; // One-body pt bin size in (GeV/c) |
2b9786f4 | 197 | |
198 | ||
20dddfab | 199 | Float_t fPhiBinSize; // One-body phi bin size in (degrees) |
2b9786f4 | 200 | |
20dddfab | 201 | Float_t fEtaBinSize ; // One-body eta bin size |
ff6aba1e | 202 | Float_t fEtaMin; // One-body eta min |
203 | Float_t fEtaMax; // One-body eta max | |
2b9786f4 | 204 | // Two-Body Histograms and Correlation Mesh for 1D and 3D distributions: |
205 | // // projections onto single axis. | |
206 | ||
20dddfab | 207 | Float_t fBinsize1dFine; // Bin Size - 1D, Fine Mesh in (GeV/c) |
208 | Float_t fBinsize1dCoarse; // Bin Size - 1D, Coarse Mesh in (GeV/c) | |
209 | Float_t fQmid1d; // q (GeV/c) at fine-coarse mesh boundary | |
210 | Float_t fQmax1d; // Max q (GeV/c) for 1D distributions | |
211 | Float_t fBinsize3dFine; // Bin Size - 3D, Fine Mesh in (GeV/c) | |
212 | Float_t fBinsize3dCoarse; // Bin Size - 3D, Coarse Mesh in (GeV/c) | |
213 | Float_t fQmid3d; // q (GeV/c) at fine-coarse mesh boundary | |
214 | Float_t fQmax3d; // Max q (GeV/c) for 3D distributions | |
215 | ||
216 | Float_t fPxMin; // Sector range in px in GeV/c | |
217 | Float_t fPxMax; //--//-- | |
2b9786f4 | 218 | Float_t fDelpx; // Mom. space sector cell size - px(GeV/c) |
219 | ||
20dddfab | 220 | Float_t fPyMin; // Sector range in py in GeV/c |
221 | Float_t fPyMax; // --//-- | |
2b9786f4 | 222 | Float_t fDelpy; // Mom. space sector cell size - py(GeV/c) |
223 | ||
20dddfab | 224 | Float_t fPzMin; // Sector range in pz in GeV/c min |
225 | Float_t fPzMax; // Sector range in pz in GeV/c max | |
2b9786f4 | 226 | Float_t fDelpz; // Mom. space sector cell size - pz(GeV/c) |
227 | ||
4c90f2a0 | 228 | |
229 | /******* P R O T E C T E D M E T H O D S *****/ | |
ff6aba1e | 230 | private: |
4c90f2a0 | 231 | public: |
232 | //conveerts Eta (pseudorapidity) to etha(azimuthal angle). Returns radians | |
233 | static Double_t EtaToTheta(Double_t arg){return 2.*TMath::ATan(TMath::Exp(-arg));} | |
234 | //converts etha(azimuthal angle) to Eta (pseudorapidity). Argument in radians | |
235 | static Double_t ThetaToEta(Double_t arg); | |
236 | //converts Degrees To Radians | |
237 | static Double_t DegreesToRadians(Double_t arg){return arg*TMath::Pi()/180.;} | |
238 | //converts Radians To Degrees | |
239 | static Double_t RadiansToDegrees(Double_t arg){return arg*180./TMath::Pi();} | |
240 | ||
2b9786f4 | 241 | ClassDef(AliGenHBTprocessor,1) // Interface class for AliMevsim |
242 | ||
243 | }; | |
4c90f2a0 | 244 | #include <iostream.h> |
2b9786f4 | 245 | #endif |