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