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