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