]>
Commit | Line | Data |
---|---|---|
11141716 | 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 | ||
803d1ab0 | 16 | /* $Id$ */ |
11141716 | 17 | |
21aa51f2 | 18 | // |
19 | // Realisation of AliLhcProcess for the fast simulation of the | |
20 | // Intra Beam Scattering process | |
21 | // in transverse and longitudinal direction. | |
22 | // Author: Andreas Morsch | |
23 | // andreas.morsch@cern.ch | |
24 | // | |
25 | ||
11141716 | 26 | #include "AliLhcProcessIBS.h" |
27 | #include "AliLHC.h" | |
28 | #include "AliLhcIRegion.h" | |
29 | #include "AliLhcBeam.h" | |
30 | ||
31 | #include <TF1.h> | |
32 | #include <TMath.h> | |
33 | #include <TCanvas.h> | |
34 | #include <TGraph.h> | |
35 | #include <TMultiGraph.h> | |
36 | ||
37 | ClassImp(AliLhcProcessIBS) | |
38 | ||
39 | Double_t func(Double_t *x, Double_t *par); | |
40 | ||
41 | AliLhcProcessIBS::AliLhcProcessIBS(AliLHC* lhc, const char* name, const char* title) | |
42 | :AliLhcProcess(lhc,name,title) | |
43 | { | |
44 | // Constructor | |
45 | } | |
46 | ||
47 | ||
48 | AliLhcProcessIBS::~AliLhcProcessIBS() | |
49 | { | |
50 | // Destructor | |
51 | ||
52 | } | |
53 | ||
54 | void AliLhcProcessIBS::Init() | |
55 | { | |
56 | // Initialization | |
57 | const Float_t r0=1.535e-16; | |
58 | ||
59 | printf("\n Initializing Process %s", GetName()); | |
60 | printf("\n ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^"); | |
61 | ||
62 | fIRegions = fAccelerator->IRegions(); | |
63 | ||
64 | for (Int_t i = 0; i < 2; i++) { | |
65 | fBeam[i] = fAccelerator->Beam(i); | |
66 | fR[i] = r0*fBeam[i]->Z()*fBeam[i]->Z()/fBeam[i]->A(); | |
67 | fE[i] = 0.938*fBeam[i]->A(); | |
68 | } | |
69 | } | |
70 | ||
71 | void AliLhcProcessIBS::Evolve(Float_t dt) | |
72 | { | |
21aa51f2 | 73 | // |
74 | // Evolve by one time step dt | |
11141716 | 75 | printf("\n Here process %s %f:", GetName(), dt); |
76 | for (Int_t i=0; i<2; i++) { | |
77 | // Density | |
78 | Float_t sige = fBeam[i]->EnergySpread(); | |
79 | Float_t avbeta = fAccelerator->AverageBeta(); | |
80 | Float_t avd = fAccelerator->AverageDisp(); | |
81 | ||
82 | Float_t gamma = fBeam[i]->Gamma(); | |
83 | Float_t epsx = fBeam[i]->Emittance()*gamma; | |
84 | Float_t epsy = epsx; | |
85 | Float_t epse = fBeam[i]->LongEmittance(); | |
86 | Float_t sigxb = TMath::Sqrt(epsx/gamma*avbeta); | |
87 | Float_t sigx = TMath::Sqrt(sigxb*sigxb+(avd*avd*sige*sige)); | |
88 | Float_t ssigx = TMath::Sqrt(epsx/gamma/avbeta); | |
89 | Float_t sigy = sigx; | |
90 | Float_t ssigy = ssigx; | |
91 | ||
92 | Float_t asd = fBeam[i]->N()*fR[i]*fR[i]*fE[i]/ | |
93 | (16.*TMath::Pi()*gamma*epsx*epsy*epse); | |
94 | ||
95 | // impact parameter | |
96 | Float_t d = sige*avd/sigx; | |
97 | // | |
98 | Float_t at = sige*TMath::Sqrt(1.-d*d)/(gamma*ssigx); | |
99 | Float_t bt = sige*TMath::Sqrt(1.-d*d)/(gamma*ssigy); | |
100 | Float_t ct = sige*TMath::Sqrt(1.-d*d)*TMath::Sqrt(4.*sigy/fR[i]); | |
101 | // | |
102 | Double_t par[3]; | |
103 | par[0] = at; | |
104 | par[1] = bt; | |
105 | par[2] = ct; | |
106 | TF1 *fct = new TF1("func",func, 0., 1., 3); | |
107 | ||
108 | Float_t f = (Float_t) 8.0*TMath::Pi()* | |
109 | fct->Integral(0., 1., par, 1.e-5); | |
110 | // | |
111 | fTaux = 1./(asd*f*(d*d-at*at/2.)); | |
112 | fTaue = 1./(asd*f*(1.-d*d)); | |
113 | // printf("\n taux, taue %f %f", taux, taue); | |
114 | // Float_t tauy = -2./at*at/asd/f; | |
115 | fBeam[i]->IncreaseEmittance(dt/fTaux, dt/fTaue); | |
116 | } | |
117 | ||
118 | } | |
119 | ||
120 | void AliLhcProcessIBS::SetMonitor(Int_t n) | |
121 | { | |
122 | // Initialize Monitor | |
123 | if (fTauxArray) delete fTauxArray; | |
124 | if (fTaueArray) delete fTaueArray; | |
125 | fTauxArray = new Float_t[n]; | |
126 | fTaueArray = new Float_t[n]; | |
127 | fNmax = n; | |
128 | } | |
129 | ||
130 | void AliLhcProcessIBS::Record() | |
131 | { | |
132 | // Record monitor quantities | |
133 | fTauxArray[fAccelerator->Nt()] = fTaux/3600.; | |
134 | fTaueArray[fAccelerator->Nt()] = fTaue/3600.; | |
135 | } | |
136 | ||
137 | ||
138 | void AliLhcProcessIBS::DrawPlots() | |
139 | { | |
140 | // Draw monitor plots | |
141 | Float_t* t = fAccelerator->TimeA(); | |
142 | ||
143 | TH1 *t1 = new TH1F("t1","Hor. IBS growth time",fNmax,0,t[fNmax]); | |
144 | t1->SetMinimum(0); | |
145 | t1->SetMaximum(fTauxArray[fNmax]*1.1); | |
146 | t1->SetStats(0); | |
147 | t1->GetXaxis()->SetTitle("t (h)"); | |
148 | t1->GetYaxis()->SetTitle("tau_x (t)"); | |
149 | ||
150 | TH1 *t2 = new TH1F("t2","Long. IBS growth time",fNmax,0,t[fNmax]); | |
151 | t2->SetMinimum(0); | |
152 | t2->SetMaximum(fTaueArray[fNmax]*1.1); | |
153 | t2->SetStats(0); | |
154 | t2->GetXaxis()->SetTitle("t (h)"); | |
155 | t2->GetYaxis()->SetTitle("tau_l (t)"); | |
156 | ||
157 | TGraph* grTaux = new TGraph(fNmax, fAccelerator->TimeA(), fTauxArray); | |
158 | grTaux->SetHistogram(t1); | |
159 | ||
160 | TGraph* grTaue = new TGraph(fNmax, fAccelerator->TimeA(), fTaueArray); | |
161 | grTaue->SetHistogram(t2); | |
162 | grTaue->SetLineStyle(2); | |
163 | ||
164 | TMultiGraph* mg = new TMultiGraph(); | |
165 | mg->Add(grTaux); | |
166 | mg->Add(grTaue); | |
167 | ||
168 | TCanvas *c3 = new TCanvas("c3","IBS", 200, 10, 700, 500); | |
169 | c3->SetGrid(); | |
170 | mg->Draw("AC"); | |
171 | mg->GetXaxis()->SetTitle("t (h)"); | |
172 | mg->GetYaxis()->SetTitle("IBS Growth Time (h)"); | |
173 | mg->Draw("AC"); | |
174 | } | |
175 | ||
176 | ||
177 | ||
178 | ||
179 | ||
e76f229f | 180 | AliLhcProcessIBS& AliLhcProcessIBS::operator=(const AliLhcProcessIBS & /*rhs*/) |
11141716 | 181 | { |
182 | // Assignment operator | |
183 | return *this; | |
184 | } | |
185 | ||
186 | Double_t func(Double_t *x, Double_t *par) | |
187 | { | |
188 | Double_t a = par[0]; | |
189 | Double_t b = par[1]; | |
190 | Double_t cc = par[2]; | |
191 | ||
21aa51f2 | 192 | const Double_t kbc = 0.5772; |
11141716 | 193 | Double_t xx = x[0]; |
194 | Double_t x2=xx*xx; | |
195 | Double_t x1=1.0-x2; | |
196 | Double_t p=1.0/TMath::Sqrt(x2+a*a*x1); | |
197 | Double_t q=1.0/TMath::Sqrt(x2+b*b*x1); | |
21aa51f2 | 198 | return (1.0-3.0*x2)*p*q*(2.0*TMath::Log(0.5*cc*(p+q))-kbc); |
11141716 | 199 | } |
200 |