]> git.uio.no Git - u/mrichter/AliRoot.git/blob - LHC/AliLhcIRegion.cxx
Updated list of classes
[u/mrichter/AliRoot.git] / LHC / AliLhcIRegion.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 //
19 // Realisation of AliLhcMonitor simulating an LHC interaction region.
20 // The interaction region is described by the two LHC beams,
21 // by the beta* and the crossing angle. 
22 // As a monitor it records the luminosity, average luminosity and beta*
23 // time evolution.
24 // Author: Andreas Morsch
25 // andreas.morsch@cern.ch
26 //
27 #include "AliLhcIRegion.h"
28 #include "AliLhcBeam.h"
29 #include "AliLHC.h"
30
31 #include <TMath.h>
32 #include <TCanvas.h>
33 #include <TPad.h>
34 #include <TGraph.h>
35 #include <TMultiGraph.h>
36 #include <TH1.h>
37
38 ClassImp(AliLhcIRegion)
39
40 AliLhcIRegion::AliLhcIRegion(AliLHC* lhc, const char* name, const char* title)
41     :TNamed(name,title)
42 {
43 // Constructor
44     fAccelerator=lhc;   
45     fLumiArray        = 0;  
46     fAverageLumiArray = 0;
47     fBetaStarArray = 0;
48 }
49
50 AliLhcIRegion::AliLhcIRegion(const AliLhcIRegion& region)
51     : TNamed(region), AliLhcMonitor(region) 
52 {
53 // copy constructor
54 }
55
56 AliLhcIRegion::~AliLhcIRegion()
57 {
58 // Destructor
59
60   if (fLumiArray)        delete fLumiArray;
61   if (fAverageLumiArray) delete fAverageLumiArray;
62   if (fBetaStarArray)    delete fBetaStarArray;
63 }
64
65
66 AliLhcIRegion& AliLhcIRegion::operator=(const  AliLhcIRegion & /*rhs*/)
67 {
68 // Assignment operator
69     return *this;
70 }
71
72 void AliLhcIRegion::Init()
73 {
74   // Initialization
75   printf("\n Initializing Interaction Region %s", GetName());
76   printf("\n ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^");
77   // Initial Luminosity
78   fBeam1 = fAccelerator->Beam(0);
79   fBeam2 = fAccelerator->Beam(1);
80   fFrequency = 3.e10/(2.*TMath::Pi()*fAccelerator->Radius());
81
82   Luminosity();
83
84   fLuminosity0 = fLuminosity;  
85   fAverageLumi=0.;
86   fBetaStar0   = fBetaStar;
87   printf("\n IR Beta*                 :%10.3e cm        ", fBetaStar);
88   printf("\n IR Initial Luminosity    :%10.3e cm^-2s^-1 ", fLuminosity);
89 }
90
91 Float_t AliLhcIRegion::Luminosity()
92 {
93   Float_t sigma1 = TMath::Sqrt(fBeam1->Emittance()*fBetaStar);
94   Float_t sigma2 = TMath::Sqrt(fBeam2->Emittance()*fBetaStar);
95   fLuminosity = fFrequency * 
96     fBeam1->N()*fBeam2->N()/(2.*TMath::Pi()*(sigma1*sigma1+sigma2*sigma2));
97   return fLuminosity;
98 }
99
100 void AliLhcIRegion::Update()
101 {
102   Luminosity();
103 }
104
105 void AliLhcIRegion::SetMonitor(Int_t n)
106 {
107   // Initialize monitors
108
109   if (fLumiArray)        delete fLumiArray;
110   if (fAverageLumiArray) delete fAverageLumiArray;
111   if (fBetaStarArray)    delete fBetaStarArray;
112
113   fLumiArray        = new Float_t[n];
114   fAverageLumiArray = new Float_t[n];
115   fBetaStarArray    = new Float_t[n];
116
117   fAverageLumiArray[0] = 0.;
118   fNmax = n;
119 }
120
121 void AliLhcIRegion::Record()
122 {
123   // Record some time dependent quantities
124   //
125   Int_t n = fAccelerator->Nt();
126   // Luminosity
127   
128   fLumiArray[n] = fLuminosity;
129   
130   // Average Luminosity respecting set-up and filling time
131   if (fAccelerator->Time() > fAccelerator->SetUpTime())
132     fAverageLumi+=fLuminosity*fAccelerator->TimeStep();
133   
134   fAverageLumiArray[n] = fAverageLumi/
135     (Float_t(n+1)*fAccelerator->TimeStep()+fAccelerator->FillingTime())/
136     fLuminosity0;
137   
138   // Beta* 
139   fBetaStarArray[n] = fBetaStar;
140 }
141
142
143 void AliLhcIRegion::DrawPlots()
144 {
145   //
146   // Draw the monitor plots
147   //
148   Float_t* t = fAccelerator->TimeA();
149   //
150   char name1[20], name2[20], hname[20];
151   sprintf(name1,"c%s",GetName());
152   sprintf(name2,"b%s",GetName());
153   char title[30];
154   sprintf(title,"Luminosity Lifetime for %s",GetName());
155   
156   //
157   sprintf(hname,"%s%d",name1,0);
158   TH1 *g1 = new TH1F(hname,"Luminosity",fNmax,0,t[fNmax]);
159   g1->SetMinimum(0);
160   g1->SetMaximum(fLumiArray[0]*1.1);
161   g1->SetStats(0);
162   g1->GetXaxis()->SetTitle("t (h)");
163   g1->GetYaxis()->SetTitle("L(t) (cm**-2 s**-1)");
164   sprintf(hname,"%s%d",name1,1);
165   TH1 *g2 = new TH1F(hname,"Luminosity",fNmax,0,t[fNmax]);
166   g2->SetMinimum(0);
167   g2->SetMaximum(1.1);
168   g2->SetStats(0);
169   g2->GetXaxis()->SetTitle("t (h)");
170   g2->GetYaxis()->SetTitle("L(t)/L0");
171   sprintf(hname,"%s%d",name1,3);
172
173   TH1 *g3 = new TH1F(hname,"Average Luminosity",fNmax,0,t[fNmax]);
174   g3->SetMinimum(0);
175   g3->SetMaximum(1.1);
176   g3->SetStats(0);
177   
178   g3->GetXaxis()->SetTitle("t (h)");
179   g3->GetYaxis()->SetTitle("L(t)/L0");
180   sprintf(hname,"%s%d",name1,3);
181
182   TH1 *g4 = new TH1F(hname,"Beta*",fNmax,0,t[fNmax]);
183   g4->SetMinimum(0);
184   g4->SetMaximum(fBetaStarArray[0]*1.1);
185   g4->SetStats(0);
186   g4->GetXaxis()->SetTitle("t (h)");
187   g4->GetYaxis()->SetTitle("beta* (cm)");
188
189   TGraph* grLumi = new TGraph(fNmax, t, fLumiArray);
190   grLumi->SetHistogram(g1);
191
192   for (Int_t i=0; i<fNmax; i++) {
193     fLumiArray[i]=fLumiArray[i]/fLuminosity0;
194   }
195   TGraph* grLumiN = new TGraph(fNmax, t, fLumiArray);
196   grLumiN->SetHistogram(g2);
197
198   TGraph* grLumiA = new TGraph(fNmax, t, fAverageLumiArray);
199   grLumiA->SetHistogram(g3);
200   grLumiA->SetLineStyle(2);
201
202   TGraph* grLumiB = new TGraph(fNmax, t, fBetaStarArray);
203   grLumiB->SetHistogram(g4);
204
205   TMultiGraph* mg = new TMultiGraph();
206   mg->Add(grLumiN);
207   mg->Add(grLumiA);
208
209  
210
211   TCanvas *c1 = new TCanvas(name1,title, 200, 10, 700, 500);
212   
213   c1->SetGrid();
214   mg->Draw("AC");
215   mg->GetXaxis()->SetTitle("t (h)");
216   mg->GetYaxis()->SetTitle("L(t)/L0 and <L>/L0");
217   mg->Draw("AC");
218
219   TCanvas *c2 = new TCanvas(name2,title, 200, 10, 700, 500);
220   c2->SetGrid();
221   grLumiB->Draw("AC");
222
223 }
224
225
226
227