]> git.uio.no Git - u/mrichter/AliRoot.git/blob - LHC/AliLHC.cxx
renamed CorrectionMatrix class
[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 AliLHC::AliLHC()
36 {
37 // Constructor
38     fIRegions   = new TList;
39     fProcesses  = new TList;
40     fNRegions   = 0;
41     fNProcesses = 0;
42     fBeams = new TObjArray(2);
43     //PH    (*fBeams)[0] = 0;
44     //PH    (*fBeams)[1] = 0;    
45     fBeams->AddAt(0,0);
46     fBeams->AddAt(0,1);    
47     fTime = 0;
48     fTimeMax = 0;
49     fTimeA = 0;
50 }
51
52 AliLHC::AliLHC(const AliLHC& lhc)
53     : TObject(lhc)
54 {
55 // copy constructor
56 }
57
58 AliLHC::~AliLHC()
59 {
60 // Destructor
61     delete fIRegions;
62     delete fProcesses;
63     delete fBeams;
64 }
65
66 void AliLHC::AddIRegion(AliLhcIRegion *region)
67 {
68 //
69 //  Add region to list   
70      fIRegions->Add(region);
71      fNRegions++;
72  }
73
74 void AliLHC::AddProcess(AliLhcProcess *process)
75 {
76 //
77 //  Add process to list   
78      fProcesses->Add(process);
79      fNProcesses++;
80  }
81
82 void AliLHC::SetBeams(AliLhcBeam *beam1, AliLhcBeam *beam2 )
83 {
84
85 //
86 //  Set the beams   
87
88     (*fBeams)[0] = beam1;
89     (*fBeams)[1] = beam2;
90 }
91
92   void AliLHC::Init()
93 {
94 // Initialisation
95     fNt    = 0;
96     fNmax  = Int_t(fTimeMax/fTimeStep); 
97     fTimeA = new Float_t[fNmax];
98     //
99     Beam(0)->Init();
100     Beam(1)->Init();
101   
102     TIter next(fIRegions);
103     AliLhcIRegion *region;
104     //
105     // Loop over generators and initialize
106     while((region = (AliLhcIRegion*)next())) {
107         region->Init();
108         region->SetMonitor(fNmax);
109     }
110     
111     Beam(0)->SetMonitor(fNmax);
112     Beam(1)->SetMonitor(fNmax);
113
114     TIter nextp(fProcesses);
115     AliLhcProcess *process;
116     //
117     // Loop over generators and initialize
118     while((process = (AliLhcProcess*)nextp())) {
119         process->Init();
120         process->SetMonitor(fNmax);
121     }  
122 }
123
124  void AliLHC::EvolveTime()
125 {
126 //
127 // Simulate Time Evolution
128 //
129     while (fTime <= fTimeMax) {
130         printf("\n Time: %f %f", fTime, fTimeStep);
131         //
132         //  Processes
133         //      
134         TIter next(fProcesses);
135         AliLhcProcess *process;
136         //
137         // Evolve for each process
138         while((process = (AliLhcProcess*)next())) {
139             process->Evolve(fTimeStep);
140             process->Record();
141         }  
142         //
143         // Update and Monitoring
144         //
145         TIter nextregion(fIRegions);
146         AliLhcIRegion *region;
147         //
148         while((region = (AliLhcIRegion*)nextregion())) {
149           printf("\n Region: %s, Luminosity %10.3e", 
150                  region->GetName(), region->Luminosity());
151           region->Update();
152           region->Record();
153         }
154         Beam(0)->Record();
155         fTimeA[fNt] = fTime/3600.;
156         fTime+=fTimeStep;
157         fNt++;
158     }
159 }
160
161 void AliLHC::Evaluate()
162 {
163   // Evaluation of the results
164   TIter nextregion(fIRegions);
165   AliLhcIRegion *region;
166   //
167   // Loop over generators and initialize
168   while((region = (AliLhcIRegion*)nextregion())) {
169     region->DrawPlots();
170   }
171   
172   TIter next(fProcesses);
173   AliLhcProcess *process;
174   //
175   // Evolve for each process
176   while((process = (AliLhcProcess*)next())) {
177     process->DrawPlots();
178   }
179   
180   Beam(0)->DrawPlots();
181 }
182    
183 AliLHC& AliLHC::operator=(const  AliLHC & /*rhs*/)
184 {
185 // Assignment operator
186     return *this;
187 }
188
189