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