]> git.uio.no Git - u/mrichter/AliRoot.git/blob - LHC/AliLHC.cxx
New version of SPD raw-data reconstruction. The format now correponds to the actual...
[u/mrichter/AliRoot.git] / LHC / AliLHC.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 // Class for a simple description of the LHC.
20 // The LHC is described by two beams,
21 // the interaction regions and the
22 // beam loss processes.
23 // Run paramters can be set in order to simulate the time evolution
24 // of emittance, number of particles per bunch and luminosity.
25 // Author: Andreas Morsch
26 // andreas.morsch@cern.ch
27
28 #include "AliLHC.h"
29 #include "AliLhcIRegion.h"
30 #include "AliLhcProcess.h"
31 #include "AliLhcBeam.h"
32
33 ClassImp(AliLHC)
34
35
36 AliLHC::AliLHC():
37     fNRegions(0),
38     fNProcesses(0),
39     fIRegions(new TList()),
40     fProcesses(new TList()),
41     fBeams(new TObjArray(2)),
42     fRadius(0.),
43     fAverageBeta(0.),
44     fAverageDisp(0.),
45     fNt(0),
46     fNmax(0),
47     fTime(0.),
48     fTimeA(0),
49     fTimeStep(0.),
50     fTimeMax(0.),
51     fFillingTime(0.),
52     fSetUpTime(0.)
53 {
54 // Constructor
55     fBeams->AddAt(0,0);
56     fBeams->AddAt(0,1);    
57 }
58
59 AliLHC::AliLHC(const AliLHC& lhc):
60     TObject(lhc),
61     fNRegions(0),
62     fNProcesses(0),
63     fIRegions(0),
64     fProcesses(0),
65     fBeams(0),
66     fRadius(0.),
67     fAverageBeta(0.),
68     fAverageDisp(0.),
69     fNt(0),
70     fNmax(0),
71     fTime(0.),
72     fTimeA(0),
73     fTimeStep(0.),
74     fTimeMax(0.),
75     fFillingTime(0.),
76     fSetUpTime(0.)
77 {
78 // Copy constructor
79 }
80
81 AliLHC::~AliLHC()
82 {
83 // Destructor
84     delete fIRegions;
85     delete fProcesses;
86     delete fBeams;
87 }
88
89 void AliLHC::AddIRegion(AliLhcIRegion *region)
90 {
91 //
92 //  Add region to list   
93      fIRegions->Add(region);
94      fNRegions++;
95  }
96
97 void AliLHC::AddProcess(AliLhcProcess *process)
98 {
99 //
100 //  Add process to list   
101      fProcesses->Add(process);
102      fNProcesses++;
103  }
104
105 void AliLHC::SetBeams(AliLhcBeam *beam1, AliLhcBeam *beam2 )
106 {
107
108 //
109 //  Set the beams   
110
111     (*fBeams)[0] = beam1;
112     (*fBeams)[1] = beam2;
113 }
114
115   void AliLHC::Init()
116 {
117 // Initialisation
118     fNt    = 0;
119     fNmax  = Int_t(fTimeMax/fTimeStep); 
120     fTimeA = new Float_t[fNmax];
121     //
122     Beam(0)->Init();
123     Beam(1)->Init();
124   
125     TIter next(fIRegions);
126     AliLhcIRegion *region;
127     //
128     // Loop over generators and initialize
129     while((region = (AliLhcIRegion*)next())) {
130         region->Init();
131         region->SetMonitor(fNmax);
132     }
133     
134     Beam(0)->SetMonitor(fNmax);
135     Beam(1)->SetMonitor(fNmax);
136
137     TIter nextp(fProcesses);
138     AliLhcProcess *process;
139     //
140     // Loop over generators and initialize
141     while((process = (AliLhcProcess*)nextp())) {
142         process->Init();
143         process->SetMonitor(fNmax);
144     }  
145 }
146
147  void AliLHC::EvolveTime()
148 {
149 //
150 // Simulate Time Evolution
151 //
152     while (fTime <= fTimeMax) {
153         printf("\n Time: %f %f", fTime, fTimeStep);
154         //
155         //  Processes
156         //      
157         TIter next(fProcesses);
158         AliLhcProcess *process;
159         //
160         // Evolve for each process
161         while((process = (AliLhcProcess*)next())) {
162             process->Evolve(fTimeStep);
163             process->Record();
164         }  
165         //
166         // Update and Monitoring
167         //
168         TIter nextregion(fIRegions);
169         AliLhcIRegion *region;
170         //
171         while((region = (AliLhcIRegion*)nextregion())) {
172           printf("\n Region: %s, Luminosity %10.3e", 
173                  region->GetName(), region->Luminosity());
174           region->Update();
175           region->Record();
176         }
177         Beam(0)->Record();
178         fTimeA[fNt] = fTime/3600.;
179         fTime+=fTimeStep;
180         fNt++;
181     }
182 }
183
184 void AliLHC::Evaluate()
185 {
186   // Evaluation of the results
187   TIter nextregion(fIRegions);
188   AliLhcIRegion *region;
189   //
190   // Loop over generators and initialize
191   while((region = (AliLhcIRegion*)nextregion())) {
192     region->DrawPlots();
193   }
194   
195   TIter next(fProcesses);
196   AliLhcProcess *process;
197   //
198   // Evolve for each process
199   while((process = (AliLhcProcess*)next())) {
200     process->DrawPlots();
201   }
202   
203   Beam(0)->DrawPlots();
204 }
205    
206 AliLHC& AliLHC::operator=(const  AliLHC & /*rhs*/)
207 {
208 // Assignment operator
209     return *this;
210 }
211
212