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