]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowTools/glauberMC/AliGlauberNucleus.cxx
cleanup of structure
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowTools / glauberMC / AliGlauberNucleus.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 //
18 //  AliGlauberNucleus implementation
19 //  support class for Glauber MC
20 //
21 //  origin: PHOBOS experiment
22 //  alification: Mikolaj Krzewicki, Nikhef, mikolaj.krzewicki@cern.ch
23 //
24 ////////////////////////////////////////////////////////////////////////////////
25
26 #include <TMath.h>
27 #include <TEllipse.h>
28 #include <TNamed.h>
29 #include <TObjArray.h>
30 #include <TF1.h>
31 #include <TRandom.h>
32 #include "AliGlauberNucleon.h"
33 #include "AliGlauberNucleus.h"
34
35 ClassImp(AliGlauberNucleus)
36
37 //______________________________________________________________________________
38 AliGlauberNucleus::AliGlauberNucleus(Option_t* iname, Int_t iN, Double_t iR, Double_t ia, Double_t iw, TF1* ifunc) : 
39    fN(iN),fR(iR),fA(ia),fW(iw),fMinDist(-1),
40    fF(0),fTrials(0),fFunction(ifunc),
41    fNucleons(0)
42 {
43    if (fN==0) {
44       cout << "Setting up nucleus " << iname << endl;
45       Lookup(iname);
46    }
47 }
48
49 //______________________________________________________________________________
50 AliGlauberNucleus::~AliGlauberNucleus()
51 {
52    if (fNucleons) {
53       delete fNucleons;
54    }
55    delete fFunction;
56 }
57
58 //______________________________________________________________________________
59 void AliGlauberNucleus::Draw(Double_t xs, Int_t col)
60 {
61    Double_t r = 0.5*sqrt(xs/TMath::Pi()/10.);
62    TEllipse e;
63    e.SetLineColor(col);
64    e.SetFillColor(0);
65    e.SetLineWidth(1);
66
67    for (Int_t i = 0;i<fNucleons->GetEntries();++i) {
68       AliGlauberNucleon* gn = (AliGlauberNucleon*) fNucleons->At(i);
69       e.SetLineStyle(1);
70       if (gn->IsSpectator()) e.SetLineStyle(3);
71       e.DrawEllipse(gn->GetX(),gn->GetY(),r,r,0,360,0,"");
72    }
73 }
74
75 //______________________________________________________________________________
76 void AliGlauberNucleus::Lookup(Option_t* name)
77 {
78    SetName(name);
79
80    if      (TString(name) == "p")    {fN = 1;   fR = 0.6;   fA = 0;      fW =  0;      fF = 0;}
81    else if (TString(name) == "d")    {fN = 2;   fR = 0.01;  fA = 0.5882; fW =  0;      fF = 1;}
82    else if (TString(name) == "dh")   {fN = 2;   fR = 0.01;  fA = 0.5882; fW =  0;      fF = 3;}
83    else if (TString(name) == "dhh")  {fN = 2;   fR = 0.01;  fA = 0.5882; fW =  0;      fF = 4;}
84    else if (TString(name) == "O")    {fN = 16;  fR = 2.608; fA = 0.513;  fW = -0.051;  fF = 1;}
85    else if (TString(name) == "Si")   {fN = 28;  fR = 3.34;  fA = 0.580;  fW = -0.233;  fF = 1;}
86    else if (TString(name) == "S")    {fN = 32;  fR = 2.54;  fA = 2.191;  fW =  0.16;   fF = 2;}
87    else if (TString(name) == "Ca")   {fN = 40;  fR = 3.766; fA = 0.586;  fW = -0.161;  fF = 1;}
88    else if (TString(name) == "Ni")   {fN = 58;  fR = 4.309; fA = 0.517;  fW = -0.1308; fF = 1;}
89    else if (TString(name) == "Cu")   {fN = 63;  fR = 4.2;   fA = 0.596;  fW =  0;      fF = 1;}
90    else if (TString(name) == "W")    {fN = 186; fR = 6.58;  fA = 0.480;  fW =  0;      fF = 1;}
91    else if (TString(name) == "Au")   {fN = 197; fR = 6.38;  fA = 0.535;  fW =  0;      fF = 1;}
92    else if (TString(name) == "Pb")   {fN = 208; fR = 6.62;  fA = 0.546;  fW =  0;      fF = 1;}
93    else if (TString(name) == "U")    {fN = 238; fR = 6.81;  fA = 0.6;    fW =  0;      fF = 1;}
94    else {
95       cout << "Could not find nucleus " << name << endl;
96       return;
97    }
98
99    switch (fF)
100    {
101       case 0: // Proton
102          fFunction = new TF1("prot","x*x*exp(-x/[0])",0,10);
103          fFunction->SetParameter(0,fR);
104          break;
105       case 1: // 3pF
106          fFunction = new TF1(name,"x*x*(1+[2]*(x/[0])**2)/(1+exp((x-[0])/[1]))",0,15);
107          fFunction->SetParameters(fR,fA,fW);
108          break;
109       case 2: // 3pG
110          fFunction = new TF1("3pg","x*x*(1+[2]*(x/[0])**2)/(1+exp((x**2-[0]**2)/[1]**2))",0,15);
111          fFunction->SetParameters(fR,fA,fW);
112          break;
113       case 3: // Hulthen
114          fFunction = new TF1("f3","x*x*([0]*[1]*([0]+[1]))/(2*pi*(pow([0]-[1],2)))*pow((exp(-[0]*x)-exp(-[1]*x))/x,2)",0,10);
115          fFunction->SetParameters(1/4.38,1/.85);
116          break;
117       case 4: // Hulthen HIJING
118          fFunction = new TF1("f4","x*x*([0]*[1]*([0]+[1]))/(2*pi*(pow([0]-[1],2)))*pow((exp(-[0]*x)-exp(-[1]*x))/x,2)",0,20);
119          fFunction->SetParameters(2/4.38,2/.85);
120          break;
121       default:
122          cerr << "Could not find function type " << fF << endl;
123          return;
124    }
125 }
126
127 //______________________________________________________________________________
128 void AliGlauberNucleus::SetR(Double_t ir)
129 {
130    fR = ir;
131    switch (fF)
132    {
133       case 0: // Proton
134          fFunction->SetParameter(0,fR);
135          break;
136       case 1: // 3pF
137          fFunction->SetParameter(0,fR);
138          break;
139       case 2: // 3pG
140          fFunction->SetParameter(0,fR);
141          break;
142    }
143 }
144
145 //______________________________________________________________________________
146 void AliGlauberNucleus::SetA(Double_t ia)
147 {
148    fA = ia;
149    switch (fF)
150    {
151       case 0: // Proton
152          break;
153       case 1: // 3pF
154          fFunction->SetParameter(1,fA);
155          break;
156       case 2: // 3pG
157          fFunction->SetParameter(1,fA);
158          break;
159    }
160 }
161
162 //______________________________________________________________________________
163 void AliGlauberNucleus::SetW(Double_t iw)
164 {
165    fW = iw;
166    switch (fF)
167    {
168       case 0: // Proton
169          break;
170       case 1: // 3pF
171          fFunction->SetParameter(2,fW);
172          break;
173       case 2: // 3pG
174          fFunction->SetParameter(2,fW);
175          break;
176    }
177 }
178
179 //______________________________________________________________________________
180 void AliGlauberNucleus::ThrowNucleons(Double_t xshift)
181 {
182    if (fNucleons==0) {
183       fNucleons=new TObjArray(fN);
184       fNucleons->SetOwner();
185       for(Int_t i=0;i<fN;i++) {
186          AliGlauberNucleon *nucleon=new AliGlauberNucleon(); 
187          fNucleons->Add(nucleon); 
188       }
189    } 
190    
191    fTrials = 0;
192
193    Double_t sumx=0;       
194    Double_t sumy=0;       
195    Double_t sumz=0;       
196
197    Bool_t hulthen = (TString(GetName())=="dh");
198    if (fN==2 && hulthen) { //special treatmeant for Hulten
199
200       Double_t r = fFunction->GetRandom()/2;
201       Double_t phi = gRandom->Rndm() * 2 * TMath::Pi() ;
202       Double_t ctheta = 2*gRandom->Rndm() - 1 ;
203       Double_t stheta = sqrt(1-ctheta*ctheta);
204      
205       AliGlauberNucleon *nucleon1=(AliGlauberNucleon*)(fNucleons->At(0));
206       AliGlauberNucleon *nucleon2=(AliGlauberNucleon*)(fNucleons->At(1));
207       nucleon1->Reset();
208       nucleon1->SetXYZ(r * stheta * cos(phi) + xshift,
209                        r * stheta * sin(phi),
210                        r * ctheta);
211       nucleon2->Reset();
212       nucleon2->SetXYZ(-nucleon1->GetX() + 2*xshift,
213                        -nucleon1->GetY(),
214                        -nucleon1->GetZ());
215       fTrials = 1;
216       return;
217    }
218
219    for (Int_t i = 0; i<fN; i++) {
220       AliGlauberNucleon *nucleon=(AliGlauberNucleon*)(fNucleons->At(i));
221       nucleon->Reset();
222       while(1) {
223          fTrials++;
224          Double_t r = fFunction->GetRandom();
225          Double_t phi = gRandom->Rndm() * 2 * TMath::Pi() ;
226          Double_t ctheta = 2*gRandom->Rndm() - 1 ;
227          Double_t stheta = TMath::Sqrt(1-ctheta*ctheta);
228          Double_t x = r * stheta * cos(phi) + xshift;
229          Double_t y = r * stheta * sin(phi);      
230          Double_t z = r * ctheta;      
231          nucleon->SetXYZ(x,y,z);
232          if(fMinDist<0) break;
233          Bool_t test=1;
234          for (Int_t j = 0; j<i; j++) {
235             AliGlauberNucleon *other=(AliGlauberNucleon*)fNucleons->At(j);
236             Double_t xo=other->GetX();
237             Double_t yo=other->GetY();
238             Double_t zo=other->GetZ();
239             Double_t dist = TMath::Sqrt((x-xo)*(x-xo)+
240                                        (y-yo)*(y-yo)+
241                                        (z-zo)*(z-zo));
242                
243             if(dist<fMinDist) {
244                test=0;
245                break;
246             }
247          }
248          if (test) break; //found nucleuon outside of mindist
249       }
250            
251       sumx += nucleon->GetX();
252       sumy += nucleon->GetY();
253       sumz += nucleon->GetZ();
254    }
255       
256    if(1) { // set the centre-of-mass to be at zero (+xshift)
257       sumx = sumx/fN;  
258       sumy = sumy/fN;  
259       sumz = sumz/fN;  
260       for (Int_t i = 0; i<fN; i++) {
261          AliGlauberNucleon *nucleon=(AliGlauberNucleon*)(fNucleons->At(i));
262          nucleon->SetXYZ(nucleon->GetX()-sumx-xshift,
263                          nucleon->GetY()-sumy,
264                          nucleon->GetZ()-sumz);
265       }
266    }
267 }
268