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 **************************************************************************/
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
26 #include "AliLhcProcessIBS.h"
28 #include "AliLhcIRegion.h"
29 #include "AliLhcBeam.h"
35 #include <TMultiGraph.h>
37 ClassImp(AliLhcProcessIBS)
39 Double_t func(Double_t *x, Double_t *par);
41 AliLhcProcessIBS::AliLhcProcessIBS(AliLHC* lhc, const char* name, const char* title)
42 :AliLhcProcess(lhc,name,title)
48 AliLhcProcessIBS::~AliLhcProcessIBS()
54 void AliLhcProcessIBS::Init()
57 const Float_t r0=1.535e-16;
59 printf("\n Initializing Process %s", GetName());
60 printf("\n ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^");
62 fIRegions = fAccelerator->IRegions();
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();
71 void AliLhcProcessIBS::Evolve(Float_t dt)
74 // Evolve by one time step dt
75 printf("\n Here process %s %f:", GetName(), dt);
76 for (Int_t i=0; i<2; i++) {
78 Float_t sige = fBeam[i]->EnergySpread();
79 Float_t avbeta = fAccelerator->AverageBeta();
80 Float_t avd = fAccelerator->AverageDisp();
82 Float_t gamma = fBeam[i]->Gamma();
83 Float_t epsx = fBeam[i]->Emittance()*gamma;
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);
90 Float_t ssigy = ssigx;
92 Float_t asd = fBeam[i]->N()*fR[i]*fR[i]*fE[i]/
93 (16.*TMath::Pi()*gamma*epsx*epsy*epse);
96 Float_t d = sige*avd/sigx;
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]);
106 TF1 *fct = new TF1("func",func, 0., 1., 3);
108 Float_t f = (Float_t) 8.0*TMath::Pi()*
109 fct->Integral(0., 1., par, 1.e-5);
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);
120 void AliLhcProcessIBS::SetMonitor(Int_t n)
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];
130 void AliLhcProcessIBS::Record()
132 // Record monitor quantities
133 fTauxArray[fAccelerator->Nt()] = fTaux/3600.;
134 fTaueArray[fAccelerator->Nt()] = fTaue/3600.;
138 void AliLhcProcessIBS::DrawPlots()
140 // Draw monitor plots
141 Float_t* t = fAccelerator->TimeA();
143 TH1 *t1 = new TH1F("t1","Hor. IBS growth time",fNmax,0,t[fNmax]);
145 t1->SetMaximum(fTauxArray[fNmax]*1.1);
147 t1->GetXaxis()->SetTitle("t (h)");
148 t1->GetYaxis()->SetTitle("tau_x (t)");
150 TH1 *t2 = new TH1F("t2","Long. IBS growth time",fNmax,0,t[fNmax]);
152 t2->SetMaximum(fTaueArray[fNmax]*1.1);
154 t2->GetXaxis()->SetTitle("t (h)");
155 t2->GetYaxis()->SetTitle("tau_l (t)");
157 TGraph* grTaux = new TGraph(fNmax, fAccelerator->TimeA(), fTauxArray);
158 grTaux->SetHistogram(t1);
160 TGraph* grTaue = new TGraph(fNmax, fAccelerator->TimeA(), fTaueArray);
161 grTaue->SetHistogram(t2);
162 grTaue->SetLineStyle(2);
164 TMultiGraph* mg = new TMultiGraph();
168 TCanvas *c3 = new TCanvas("c3","IBS", 200, 10, 700, 500);
171 mg->GetXaxis()->SetTitle("t (h)");
172 mg->GetYaxis()->SetTitle("IBS Growth Time (h)");
180 AliLhcProcessIBS& AliLhcProcessIBS::operator=(const AliLhcProcessIBS & /*rhs*/)
182 // Assignment operator
186 Double_t func(Double_t *x, Double_t *par)
190 Double_t cc = par[2];
192 const Double_t kbc = 0.5772;
196 Double_t p=1.0/TMath::Sqrt(x2+a*a*x1);
197 Double_t q=1.0/TMath::Sqrt(x2+b*b*x1);
198 return (1.0-3.0*x2)*p*q*(2.0*TMath::Log(0.5*cc*(p+q))-kbc);