]> git.uio.no Git - u/mrichter/AliRoot.git/blame - LHC/AliLhcBeam.cxx
Bug fix. Removed delete statement
[u/mrichter/AliRoot.git] / LHC / AliLhcBeam.cxx
CommitLineData
11141716 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
803d1ab0 16/* $Id$ */
21aa51f2 17//
18// Class that holds all parameters about an LHC beam.
19// The parameters can change with time.
20// A monitor can be set that stores the time distribution of
21// emittance and number of particles per bunch.
22// Author: Andreas Morsch
23// andreas.morsch@cern.ch
24//
11141716 25
11141716 26#include <TCanvas.h>
27#include <TGraph.h>
0b28fd57 28#include <TH1F.h>
29#include <TMath.h>
11141716 30#include <TMultiGraph.h>
31
0b28fd57 32#include "AliLhcBeam.h"
33#include "AliLHC.h"
34
11141716 35ClassImp(AliLhcBeam)
36
64eb707f 37AliLhcBeam::AliLhcBeam(AliLHC* lhc):
38 fAccelerator(lhc),
39 fN(0),
40 fN0(0),
41 fNEmittance(0.),
42 fEmittance(0.),
43 fEmittance0(0.),
44 fEmittanceL(0.),
45 fEmittanceL0(0.),
46 fEnergySpread(0.),
47 fA(0),
48 fZ(0),
49 fEnergy(0.),
50 fGamma(0.),
51 fTimeArray(0),
52 fEmittanceArray(0),
53 fEmittanceLArray(0)
11141716 54{
55// Constructor
11141716 56}
57
64eb707f 58AliLhcBeam::AliLhcBeam(const AliLhcBeam& beam):
59 TNamed(beam), AliLhcMonitor(beam),
60 fAccelerator(0),
61 fN(0),
62 fN0(0),
63 fNEmittance(0.),
64 fEmittance(0.),
65 fEmittance0(0.),
66 fEmittanceL(0.),
67 fEmittanceL0(0.),
68 fEnergySpread(0.),
69 fA(0),
70 fZ(0),
71 fEnergy(0.),
72 fGamma(0.),
73 fTimeArray(0),
74 fEmittanceArray(0),
75 fEmittanceLArray(0)
11141716 76{
77// copy constructor
78}
79
80AliLhcBeam::~AliLhcBeam()
81{
82// Destructor
83
84}
85
86void AliLhcBeam::Init()
87{
88 // Initialization
89 printf("\n Initializing Beam");
90 printf("\n ^^^^^^^^^^^^^^^^^");
91 // Scale energy with regidity
92 fEnergy *= fZ/fA;
93 fGamma = fEnergy/0.938272;
94 fEmittance = fNEmittance/fGamma;
95 fEmittance0 = fEmittance;
96 fEmittanceL *= fZ;
97 fEmittanceL0 = fEmittanceL;
98 fN0=fN;
99
100 printf("\n Beam Energy :%10.3e GeV", fEnergy);
101 printf("\n Beam Normalized Emittance :%10.3e cm ", fNEmittance);
102 printf("\n Beam Particles per Bunch :%10.3e ", fN);
103}
104
105void AliLhcBeam::RemoveParticles(Float_t loss)
106{
107 fN-=loss;
108}
109
110void AliLhcBeam::IncreaseEmittance(Float_t de, Float_t del)
111{
21aa51f2 112//
113// Increase the emittance
11141716 114 fEmittance *= (1.+de);
115 fEmittanceL *= (1.+del);
116 fEnergySpread *= (1.+del);
117}
118
e76f229f 119AliLhcBeam& AliLhcBeam::operator=(const AliLhcBeam & /*rhs*/)
11141716 120{
121// Assignment operator
122 return *this;
123}
124void AliLhcBeam::SetMonitor(Int_t n)
125{
21aa51f2 126//
127// Initialize a monitor with n time bins
11141716 128 fNmax = n;
129 if (fEmittanceArray) delete fEmittanceArray;
130 if (fEmittanceLArray) delete fEmittanceLArray;
131
132
133 fEmittanceArray = new Float_t[n];
134 fEmittanceLArray = new Float_t[n];
135}
136
137void AliLhcBeam::Record()
138{
139 fEmittanceArray [fAccelerator->Nt()] = fEmittance/fEmittance0;
140 fEmittanceLArray[fAccelerator->Nt()] = fEmittanceL/fEmittanceL0;
141}
142
143
144void AliLhcBeam::DrawPlots()
145{
146 // Draw monitor plots
147 Float_t* t = fAccelerator->TimeA();
148
149 TH1 *e1 = new TH1F("e1","Hor. Emittance",fNmax,0,t[fNmax]);
150 e1->SetMinimum(1);
151 e1->SetMaximum(fEmittanceArray[fNmax]*1.1);
152 e1->SetStats(0);
153 e1->GetXaxis()->SetTitle("t (h)");
154 e1->GetYaxis()->SetTitle("rel. Emittance (t)");
155
156 TH1 *e2 = new TH1F("e2","Long. Emittance",fNmax,0,t[fNmax]);
157 e2->SetMinimum(1);
158 e2->SetMaximum(fEmittanceLArray[fNmax]*1.1);
159 e2->SetStats(0);
160 e2->GetXaxis()->SetTitle("t (h)");
161 e2->GetYaxis()->SetTitle("rel. Emittance (t)");
162
163
164 TGraph* grE = new TGraph(fNmax, t, fEmittanceArray);
165 grE->SetHistogram(e1);
166 TGraph* grEl = new TGraph(fNmax, t, fEmittanceLArray);
167 grEl->SetHistogram(e2);
168 grEl->SetLineStyle(2);
169
170 TMultiGraph* mg = new TMultiGraph();
171 mg->Add(grE);
172 mg->Add(grEl);
173
174 TCanvas *c2 = new TCanvas("c2","Emittance Increase", 200, 10, 700, 500);
175 c2->SetGrid();
176 mg->Draw("AC");
177 mg->GetXaxis()->SetTitle("t (h)");
178 mg->GetYaxis()->SetTitle("rel. Emittance(t)");
179 mg->Draw("AC");
180}
181
182
183
184