]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/Fragment.C
Using AliPHOSLoader instead of AliPHOSGetter
[u/mrichter/AliRoot.git] / ZDC / Fragment.C
1 void Fragment(Int_t nev=0, Int_t debug=0)
2 {
3  gROOT->Reset();
4  gStyle->SetPalette(1);
5  gStyle->SetCanvasColor(33);
6  gStyle->SetFrameFillColor(38);
7  
8  Float_t b;
9  Int_t nspectp, nspectn, nspectpfree, nspectnfree; 
10  Int_t zz[100], nn[100], nalpha, ztot, ntot;
11  void spectator(Float_t, Int_t*, Int_t*);
12  
13  TH2F *hsp = new TH2F("hsp","Spectator protons vs b",100,0.,20.,100,0.,150.);
14  hsp -> SetXTitle("b (fm)");
15  hsp -> SetYTitle("Num. of spectator protons");
16  
17  TH2F *hsn = new TH2F("hsn","Spectator neutrons vs b",100,0.,20.,100,0.,150.);
18  hsn -> SetXTitle("b (fm)");
19  hsn -> SetYTitle("Num. of spectator neutrons");
20  
21  TH2F *hFragp = new TH2F("hFragp","Free spectator protons vs b",100,0.,20.,100,0.,150.);
22  hFragp -> SetXTitle("b (fm)");
23  hFragp -> SetYTitle("Num. of free spectator protons");
24  
25  TH2F *hFragn = new TH2F("hFragn","Free spectator neutrons vs b",100,0.,20.,100,0.,150.);
26  hFragn -> SetXTitle("b (fm)");
27  hFragn -> SetYTitle("Num. of free spectator neutrons");
28
29  TF1 *funb = new  TF1("funb","x",0,15);
30  for(Int_t ievent=0; ievent<nev; ievent++){  
31    if(ievent%1000==0){printf("Processing event %d\n",ievent);}
32    b= Float_t(funb->GetRandom()); 
33    spectator(b, &nspectp, &nspectn);
34    hsp -> Fill(b,nspectp);
35    hsn -> Fill(b,nspectn);
36    AliZDCFragment *gallio = new AliZDCFragment(b);
37    for(Int_t j=0; j<=99; j++){
38       zz[j] =0;
39       nn[j] =0;
40    }
41 //
42 // Generation of fragments
43    gallio->GenerateIMF(zz, nalpha);
44    Int_t NFrag = gallio->GetFragmentNum();
45 //
46 // Attach neutrons
47    ztot=0;
48    ntot=0;
49    gallio->AttachNeutrons(zz, nn, ztot, ntot);
50    nspectpfree = nspectp-ztot-2*nalpha;
51    nspectnfree = nspectn-ntot-2*nalpha;
52    hFragp -> Fill(b,nspectpfree);
53    hFragn -> Fill(b,nspectnfree);
54 //
55 // Print
56    if(debug == 1){
57      printf("\n b = %f\n",b);
58      printf("\n #spect p = %d, #spect n = %d\n",nspectp,nspectn); 
59      printf("\n #spect p free = %d, #spect n free = %d\n",nspectpfree,nspectnfree); 
60      printf("\n #fragments = %f\n",NFrag);
61      for(Int_t i=0; i<NFrag; i++){
62         printf("\n      ZZ[%d] = %d, NN[%d] = %d\n",i,zz[i],i,nn[i]);
63      }
64      printf("\n NAlpha = %d, Ztot = %d, Ntot = %d\n\n",nalpha,ztot,ntot);
65    }
66    delete gallio;
67   } //Event loop
68
69    TCanvas *c1 = new TCanvas("c1","Fragmentation",0,10,580,700);
70    c1->cd();
71    c1->SetFillColor(29);
72    c1->Divide(2,2);
73    c1->cd(1);
74    hsp -> Draw("box");   
75    c1->cd(2);
76    hsn -> Draw("box");    
77    c1->cd(3);
78    hFragp -> Draw("box");   
79    c1->cd(4);
80    hFragn -> Draw("box");    
81
82 }
83
84
85 void spectator(Float_t b, Int_t* NSpectp, Int_t* NSpectn)
86 {
87    Float_t SppvsB[6] = {3.633,-1.518,1.360,-.4899e-1,-.2398e-2,.1066e-3};
88    Float_t SpnvsB[6] = {5.639,-1.685,1.803,-.3129e-1,-.6618e-2,.2352e-3};
89    Float_t Sigmap[4] = {.5668,-.2200e-1,.3657e-3,-.2201e-5};
90    Float_t Sigman[4] = {.4185,-.9798e-2,.1052e-3,-.4238e-6};
91    
92    Float_t rnsp = SppvsB[0]+SppvsB[1]*b+SppvsB[2]*(b*b)+SppvsB[3]*(b*b*b)+
93                   SppvsB[4]*(b*b*b*b)+SppvsB[5]*(b*b*b*b*b);
94    Float_t rnsn = SpnvsB[0]+SpnvsB[1]*b+SpnvsB[2]*(b*b)+SpnvsB[3]*(b*b*b)+
95                   SpnvsB[4]*(b*b*b*b)+SpnvsB[5]*(b*b*b*b*b);
96    Float_t snsp = Sigmap[0]+Sigmap[1]*rnsp+Sigmap[2]*(rnsp*rnsp)+Sigmap[3]*
97                   (rnsp*rnsp*rnsp); 
98    Float_t snsn = Sigman[0]+Sigman[1]*rnsn+Sigman[2]*(rnsn*rnsn)+Sigman[3]*
99                   (rnsn*rnsn*rnsn); 
100    
101    snsp = snsp*rnsp;
102    snsn = snsn*rnsn;
103    
104    Float_t xgaup = gRandom->Gaus(0.0,1.0);
105    snsp = snsp*xgaup;
106    Float_t xgaun = gRandom->Gaus(0.0,1.0);
107    snsn = snsn*xgaun;
108    rnsp=rnsp+snsp;
109    rnsn=rnsn+snsn;
110    
111    *NSpectp = Int_t(rnsp);
112    *NSpectn = Int_t(rnsn);
113
114 }