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