1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 ////////////////////////////////////////////////////////////////////////////////
18 // AliGlauberNucleus implementation
19 // support class for Glauber MC
21 // origin: PHOBOS experiment
22 // alification: Mikolaj Krzewicki, Nikhef, mikolaj.krzewicki@cern.ch
24 ////////////////////////////////////////////////////////////////////////////////
26 #include <Riostream.h>
30 #include <TObjArray.h>
33 #include "AliGlauberNucleon.h"
34 #include "AliGlauberNucleus.h"
36 ClassImp(AliGlauberNucleus)
38 //______________________________________________________________________________
39 AliGlauberNucleus::AliGlauberNucleus(Option_t* iname, Int_t iN, Double_t iR, Double_t ia, Double_t iw, TF1* ifunc) :
52 cout << "Setting up nucleus " << iname << endl;
57 //______________________________________________________________________________
58 AliGlauberNucleus::~AliGlauberNucleus()
66 //______________________________________________________________________________
67 AliGlauberNucleus::AliGlauberNucleus(const AliGlauberNucleus& in):
73 fMinDist(in.fMinDist),
76 fFunction(in.fFunction),
81 fNucleons=static_cast<TObjArray*>((in.fNucleons)->Clone());
84 //______________________________________________________________________________
85 AliGlauberNucleus& AliGlauberNucleus::operator=(const AliGlauberNucleus& in)
88 if (&in==this) return *this;
96 fFunction=in.fFunction;
98 fNucleons=static_cast<TObjArray*>((in.fNucleons)->Clone());
99 fNucleons->SetOwner();
103 //______________________________________________________________________________
104 void AliGlauberNucleus::Draw(Double_t xs, Int_t col)
106 Double_t r = 0.5*sqrt(xs/TMath::Pi()/10.);
112 for (Int_t i = 0;i<fNucleons->GetEntries();++i) {
113 AliGlauberNucleon* gn = (AliGlauberNucleon*) fNucleons->UncheckedAt(i);
115 if (gn->IsSpectator()) e.SetLineStyle(3);
116 e.DrawEllipse(gn->GetX(),gn->GetY(),r,r,0,360,0,"");
120 //______________________________________________________________________________
121 void AliGlauberNucleus::Lookup(Option_t* name)
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;}
140 cout << "Could not find nucleus " << name << endl;
147 fFunction = new TF1("prot","x*x*exp(-x/[0])",0,10);
148 fFunction->SetParameter(0,fR);
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);
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);
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);
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);
167 cerr << "Could not find function type " << fF << endl;
172 //______________________________________________________________________________
173 void AliGlauberNucleus::SetR(Double_t ir)
179 fFunction->SetParameter(0,fR);
182 fFunction->SetParameter(0,fR);
185 fFunction->SetParameter(0,fR);
190 //______________________________________________________________________________
191 void AliGlauberNucleus::SetA(Double_t ia)
199 fFunction->SetParameter(1,fA);
202 fFunction->SetParameter(1,fA);
207 //______________________________________________________________________________
208 void AliGlauberNucleus::SetW(Double_t iw)
216 fFunction->SetParameter(2,fW);
219 fFunction->SetParameter(2,fW);
224 //______________________________________________________________________________
225 void AliGlauberNucleus::ThrowNucleons(Double_t xshift)
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);
242 Bool_t hulthen = (TString(GetName())=="dh");
243 if (fN==2 && hulthen) { //special treatmeant for Hulten
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);
250 AliGlauberNucleon *nucleon1=(AliGlauberNucleon*)(fNucleons->UncheckedAt(0));
251 AliGlauberNucleon *nucleon2=(AliGlauberNucleon*)(fNucleons->UncheckedAt(1));
253 nucleon1->SetXYZ(r * stheta * cos(phi) + xshift,
254 r * stheta * sin(phi),
257 nucleon2->SetXYZ(-nucleon1->GetX() + 2*xshift,
264 for (Int_t i = 0; i<fN; i++) {
265 AliGlauberNucleon *nucleon=(AliGlauberNucleon*)(fNucleons->UncheckedAt(i));
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;
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)+
293 if (test) break; //found nucleuon outside of mindist
296 sumx += nucleon->GetX();
297 sumy += nucleon->GetY();
298 sumz += nucleon->GetZ();
301 if(1) { // set the centre-of-mass to be at zero (+xshift)
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);