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