]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ALIFAST/AliFast.cxx
negative indexes allowed
[u/mrichter/AliRoot.git] / ALIFAST / AliFast.cxx
1
2 //////////////////////////////////////////////////////////////////////////
3 //                                                                      //
4 // AliFast                                                              //
5 //                                                                      //
6 // Main class to control the AliFast program.                           //
7 //                                                                      //
8 // This class is a general framework for programs that needs to:        //
9 //    - Initialise some parameters                                      //
10 //    - Loop on events                                                  //
11 //    - Print results and save histograms, etc                          //
12 //                                                                      //
13 // The event processor AliFast::Make loops on a list of Makers          //
14 // where each maker performs some task on the event data and generates  //
15 // results.                                                             //
16 // New Makers can be inserted by a user without modifying this class.   //
17 // Note that the order in which the Makers are called is the order      //
18 // of insertion in the list of Makers.                                  //
19 // Each Maker is responsible for creating its branch of the Tree.       //
20 // The following table shows the list of makers currently implemented   //
21 // The default option to Save the Maker info in the Tree is mentioned.  //
22 //                                                                      //
23 //    Maker name        Save in Tree                                    //
24 //    ==========        ============                                    //
25 //    MCMaker             NO                                            //
26 //    TrackMaker          NO                                            //
27 //                                                                      //
28 // Makers must derive from the base class AliFMaker.                    //
29 // AliFMaker provides a common interface to all Makers.                 //
30 // Each Maker is responsible for defining its own parameters and        //
31 // histograms.                                                          //
32 // Each Maker has its own list of histograms.                           //
33 // Each Maker has an associated companion class corresponding to the    //
34 // type of physics object reconstructed by the Maker.                   //
35 // For example, AliFClusterMaker creates AliFCluster objects.           //
36 //              AliFTriggerMaker creates one single AliFTrigger object. //
37 // The pointer supporting the created object(s) is defined in AliFMaker //
38 //   fFruits may point to a single object (eg. AliFTrigger) or to a    //
39 //           TClonesArray of objects (eg. AliFCluster).                 //
40 //                                                                      //
41 // The function AliFast::Maketree must be called after the creation     //
42 // of the AliFast object to create a Root Tree.                         //
43 //                                                                      //
44 // An example of main program/macro to use AliFast is given below:      //
45 //========================================================================
46 //void umain(Int_t nevents=100)
47 //{
48 //   gROOT->Reset();
49 //   gSystem->Load("libalifast.so");  // dynamically link the compiled shared library
50 //
51 //   // Open the root output file
52 //   TFile file("alifast.root","recreate","AliFast root file",2);
53 //   
54 //   AliFast alifast("alifast");     // create main object to run alifast
55 //
56 //   User user;           // create an object of the User class defined in user.C
57 //
58 //   alifast.Init();      // Initialise event (maker histograms,etc)
59 //   alifast.MakeTree();  // Create the Root tree
60 //
61 //   gROOT->LoadMacro("user.C");  // compile/interpret user file
62 //
63 //   for (Int_t i=0; i<nevents; i++) {
64 //      if (i%100 == 0) printf("In loop:%d\n",i);
65 //      alifast.Make(i);       // Generate and reconstruct event
66 //      user.FillHistograms(); // User has possibility to decide if store event here!
67 //      alifast.FillTree();
68 //      alifast.Clear();       // Clear all event lists
69 //   }
70 //   alifast.Finish();
71 //
72 //   // save objects in Root file
73 //   alifast.Write();  //save main alifast object (and run parameters)
74 //}
75 //========================================================================
76 //                                                                      //
77 // This example illustrates how to:                                     //
78 //    - Load a shared library                                           //
79 //    - Open a Root file                                                //
80 //    - Initialise AliFast                                              //
81 //    - Load some user code (interpreted)                               //
82 //      This user code may redefine some Maker parameters               //
83 //    - Make a loop on events                                           //
84 //    - Save histograms and the main AliFast object and its Makers      //
85 //                                                                      //
86 //========================================================================
87 //  An example of a User class is given below:                          //
88 //========================================================================
89 //
90 //#ifndef user_H
91 //#define user_H
92 //
93 //////////////////////////////////////////////////////////////////////////
94 //                                                                      //
95 // User                                                                 //
96 //                                                                      //
97 // Example of a user class to perform user specific tasks when running  //
98 // the ALIfast program.                                                 //
99 //                                                                      //
100 // This class illustrates:                                              //
101 //   - How to set run parameters                                        //
102 //   - How to create and fill histograms                                //
103 //                                                                      //
104 //////////////////////////////////////////////////////////////////////////
105 //
106 //class TH1F;
107 //class AliFast;
108 //class AliFClusterMaker;
109 //class AliFPhotonMaker;
110 //
111 //class User {
112 //
113 //private:
114 //   TH1F             *fhist1;       //pointer to histogram
115 //   TH1F             *fhist2;       //pointer to histogram
116 //   TH1F             *fhist3;       //pointer to histogram
117 //public:
118 //               User();
119 //   void        FillHistograms();
120 //   void        SetRunParameters();
121 //
122 //#endif
123 //};
124 //
125 //_________________________________________________________________________
126 //User::User() 
127 //{
128 //   SetRunParameters();  //change default parameters
129 //
130 //         Create a few histograms
131 //   fhist1 = new TH1F("hist1","Number of tracks per event",100,0,100);
132 //   fhist2 = new TH1F("hist2","Number of clusters",100,0,100);
133 //   fhist3 = new TH1F("hist3","Number of isolated muons",20,0,20);
134 //}
135 //
136 //_________________________________________________________________________
137 //void User::FillHistograms()
138 //{
139 ////   fhist1.Fill(event->GetNtracks());
140 ////   fhist2.Fill(event->GetNclusters));
141 ////   fhist3.Fill(event->GetNIsoMuons());
142 //}
143 //
144 //_________________________________________________________________________
145 //void User::SetRunParameters()
146 //{
147 //  // change Alifast default parameters
148 //
149 //   gAliFast->SetSmearMuonOpt(0);
150 //   gAliFast->ClusterMaker()->SetGranBarrelEta(0.12);
151 //   gAliFast->PhotonMaker()->SetMinPT(6.);
152 //   gAliFast->TriggerMaker()->SetMuoEtaCoverage(2.8);
153 //
154 //}
155 //======================end of User class=================================
156 //
157 //////////////////////////////////////////////////////////////////////////
158
159 #include <TROOT.h>
160 #include <TChain.h>
161 #include <TTree.h>
162 #include <TBrowser.h>
163 #include <TClonesArray.h>
164 #include "AliRun.h"
165 #include "AliFast.h"
166 //#include "AliFMCMaker.h"
167 #include "AliFTrackMaker.h"
168 #include "AliFHistBrowser.h"
169 #include "AliFBigBang.h"
170 #include "AliFVirtualDisplay.h"
171
172
173 R__EXTERN AliRun * gAlice;
174 AliFast *gAliFast;
175
176 ClassImp(AliFast)
177
178
179 //_____________________________________________________________________________
180   //AliFast::AliFast() : TNamed("alifast","The ALICE fast simulation")
181   AliFast::AliFast() : AliRun("alifast","The ALICE fast simulation")
182 {
183
184    fTree          = 0;
185    fMakers        = 0;
186    fMode          = 0;
187    fMCMaker       = 0;
188    fTrackMaker    = 0;
189    fDisplay       = 0;
190    fDet           = new AliFDet("Detector","Make AliFast detector");
191    gAliFast     = this;
192    gAlice       = (AliRun*)this;
193 }
194
195 //_____________________________________________________________________________
196 //AliFast::AliFast(const char *name, const char *title): TNamed(name,title)
197 AliFast::AliFast(const char *name, const char *title): AliRun(name,title)
198 {
199
200    gAliFast      = this;
201    fVersion     = 001;       //AliFAST  version number and release date
202    fVersionDate = 150399;
203    fTree        = 0;
204    fMode        = 0;
205    fDisplay     = 0;
206    
207    SetDefaultParameters();
208
209    gROOT->GetListOfBrowsables()->Add(this,"AliFast");
210
211 // create the support list for the various lists of AliFast objects
212    fMakers  = new TList();
213
214 // create "standard" makers and add them to the list of makers (in AliFMaker constructor
215 // Note that the order in which makers are added to the list of makers is important
216 // makers will be processed in this order !!
217
218    //fMCMaker       = new AliFMCMaker("MCMaker","Make MC events");
219    fTrackMaker    = new AliFTrackMaker("TrackMaker","Make AliFast tracks");
220 //create detector
221    fDet           = new AliFDet("Detector","Make AliFast detector");
222 }
223
224 //_____________________________________________________________________________
225 AliFast::~AliFast()
226 {
227 //   fMakers->Delete();
228 //   delete fMakers;
229 }
230
231
232 //______________________________________________________________________________
233 void AliFast::Browse(TBrowser *b)
234 {
235
236   if( b == 0) return;
237
238   if (fTree) b->Add(fTree,fTree->GetName());
239   // fca
240   b->Add(&fHistBrowser, "Histograms");
241   b->Add(&fBigBang, "BigBang");
242   // fca
243
244   TIter next(fMakers);
245   AliFMaker *maker;
246   while ((maker = (AliFMaker*)next())) {
247      b->Add(maker,maker->GetName());
248    }
249 }
250
251 //_____________________________________________________________________________
252 void AliFast::Clear(Option_t *option)
253 {
254 //    Reset lists of event objects
255    TIter next(fMakers);
256    AliFMaker *maker;
257    while ((maker = (AliFMaker*)next())) {
258       maker->Clear(option);
259    }
260    //fca   if (fDisplay) fDisplay->Clear();
261 }
262
263 //_____________________________________________________________________________
264 void AliFast::Draw(Option_t *option)
265 {
266 //    Insert current event in graphics pad list
267
268     // Check if the Event Display object has been created
269    if (!fDisplay) {
270       Error("Draw","You must create an AliFDisplay object first");
271       return;
272    }
273
274    //fca   fDisplay->Draw(option);
275 }
276
277 //_____________________________________________________________________________
278 void  AliFast::GetTreeEvent(Int_t event)
279 {
280 //    Read event from Tree
281    if (fTree) fTree->GetEvent(event);
282    fEvent = event;  
283 }
284
285 //_____________________________________________________________________________
286 void AliFast::Init()
287 {
288 //  Initialise detector
289    AliFDet *detector=gAliFast->Detector();  
290    detector->InitDetParam(); 
291
292 //    Initialise makers
293    TIter next(fMakers);
294    AliFMaker *maker;
295    TObject *objfirst, *objlast;
296    while ((maker = (AliFMaker*)next())) {
297      // save last created histogram in current Root directory
298       objlast = gDirectory->GetList()->Last();
299
300      // Initialise maker
301       maker->Init();
302
303      // Add Maker histograms in Maker list of histograms
304       if (objlast) objfirst = gDirectory->GetList()->After(objlast);
305       else         objfirst = gDirectory->GetList()->First();
306       while (objfirst) {
307          maker->Histograms()->Add(objfirst);
308          objfirst = gDirectory->GetList()->After(objfirst);
309       }
310    }
311
312 }
313
314 //_____________________________________________________________________________
315 void AliFast::Paint(Option_t *option)
316 {
317 //    Paint AliFast objects
318
319   //fca   fDisplay->Paint(option);
320 }
321
322 //_____________________________________________________________________________
323 void AliFast::PrintInfo()
324 {
325 //     Gives information about versions etc.
326    printf("\n\n");
327    printf("**************************************************************\n");
328    printf("*             AliFast version:00%1d     last update  %6d    *\n",
329                                                        fVersion, fVersionDate);
330    printf("**************************************************************\n");
331    printf("*                                                            *\n");
332    printf("*     Simulates and reconstructs  events on particle level   *\n");
333    printf("*     Package by: Yiota Foka and Elzbieta Richter-Was        *\n");
334    printf("*           Design based on ATLFast++                        *\n");
335    printf("*         by R. Brun and E. Richter-Was                      *\n");
336    printf("**************************************************************\n");
337    printf("\n\n");
338
339 //     Print info for detector geometry
340    AliFDet *detector=gAliFast->Detector();
341    detector->PrintDetInfo(); 
342
343 //     Print info for all defined Makers
344    TIter next(fMakers);
345    AliFMaker *maker;
346    while ((maker = (AliFMaker*)next())) {
347       maker->PrintInfo();
348    }
349 }
350
351 //_____________________________________________________________________________
352 void AliFast::FillTree()
353 {
354 //  Fill the ROOT tree, looping on all active branches
355
356   // Clean generated particles (depending on option Save)
357   //MCMaker()->CleanParticles();
358
359   // Now ready to fill the Root Tree
360    if(fTree) fTree->Fill();
361 }
362
363 //_____________________________________________________________________________
364 void AliFast::InitChain(TChain *chain)
365 {
366 //  Initialize branch addresses for all makers in a TChain
367
368    if (chain == 0) return;
369
370    fTree = chain;
371
372    TIter next(fMakers);
373    AliFMaker *maker;
374    while ((maker = (AliFMaker*)next())) {
375       maker->SetChainAddress(chain);
376    }
377 }
378
379 //_____________________________________________________________________________
380 void AliFast::MakeTree(const char* name, const char*title)
381 {
382 //  Create a ROOT tree
383 //  Loop on all makers to create the Root branch (if any)
384
385    if (fTree) return;
386
387    fTree = new TTree(name,title);
388
389    TIter next(fMakers);
390    AliFMaker *maker;
391    while ((maker = (AliFMaker*)next())) {
392       maker->MakeBranch();
393    }
394 }
395
396 //_____________________________________________________________________________
397 void AliFast::SetDefaultParameters()
398 {
399
400 //    Setters for flags and switches
401    SetLuminosity();
402    SetBfield();
403    SetSmearing();
404    SetSUSYcodeLSP();
405    SetTrackFinding();
406 }
407
408 //_____________________________________________________________________________
409 void AliFast::Make(Int_t i)
410 {
411    fEvent = i;
412
413 //   Loop on all makers
414    TIter next(fMakers);
415    AliFMaker *maker;
416    while ((maker = (AliFMaker*)next())) {
417       maker->Make();
418    }
419
420 }
421
422 //_____________________________________________________________________________
423 void AliFast::FillClone()
424 {
425    // Fill Makers fruits clones
426    
427    TIter next(fMakers);
428    AliFMaker *maker;
429    while ((maker = (AliFMaker*)next())) {
430       maker->FillClone();
431    }
432 }
433
434 //_____________________________________________________________________________
435 void AliFast::Finish()
436 {
437 //    Terminate a run
438 //   place to make operations on histograms, normalization,etc.
439
440    TIter next(fMakers);
441    AliFMaker *maker;
442    while ((maker = (AliFMaker*)next())) {
443       maker->Finish();
444    }
445 }
446
447 //_____________________________________________________________________________
448 void AliFast::SortDown(Int_t n1, Float_t *a, Int_t *index, Bool_t down)
449 {
450 //  sort the n1 elements of array a.
451 //  In output the array index contains the indices of the sorted array.
452 //  if down is false sort in increasing order (default is decreasing order)
453 //   This is a translation of the CERNLIB routine sortzv (M101)
454 //   based on the quicksort algorithm
455
456    Int_t i,i1,n,i2,i3,i33,i222,iswap,n2;
457    Int_t i22 = 0;
458    Float_t ai;
459    n = n1;
460    if (n <= 0) return;
461    if (n == 1) {index[0] = 0; return;}
462    for (i=0;i<n;i++) index[i] = i+1;
463    for (i1=2;i1<=n;i1++) {
464       i3 = i1;
465       i33 = index[i3-1];
466       ai  = a[i33-1];
467       while(1) {
468          i2 = i3/2;
469          if (i2 <= 0) break;
470          i22 = index[i2-1];
471          if (ai <= a[i22-1]) break;
472          index[i3-1] = i22;
473          i3 = i2;
474       }
475       index[i3-1] = i33;
476    }
477
478    while(1) {
479       i3 = index[n-1];
480       index[n-1] = index[0];
481       ai = a[i3-1];
482       n--;
483       if(n-1 < 0) {index[0] = i3; break;}
484       i1 = 1;
485       while(2) {
486          i2 = i1+i1;
487          if (i2 <= n) i22 = index[i2-1];
488          if (i2-n > 0) {index[i1-1] = i3; break;}
489          if (i2-n < 0) {
490             i222 = index[i2];
491             if (a[i22-1] - a[i222-1] < 0) {
492                 i2++;
493                 i22 = i222;
494             }
495          }
496          if (ai - a[i22-1] > 0) {index[i1-1] = i3; break;}
497          index[i1-1] = i22;
498          i1 = i2;
499       }
500    }
501    if (!down) return;
502    n2 = n1/2;
503    for (i=0;i<n1;i++) index[i]--;
504    for (i=0;i<n2;i++) {
505       iswap         = index[i];
506       index[i]      = index[n1-i-1];
507       index[n1-i-1] = iswap;
508    }
509 }
510 //______________________________________________________________________________
511 void AliFast::Streamer(TBuffer &R__b)
512 {
513    // Stream an object of class AliFast.
514
515    if (R__b.IsReading()) {
516       UInt_t R__s, R__c;
517       if (!gAliFast) gAliFast = this;
518       gROOT->GetListOfBrowsables()->Add(this,"AliFast");
519       
520       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
521       
522       AliFast::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c);
523
524       fTree = (TTree*)gDirectory->Get("T");
525    } else {
526        AliFast::Class()->WriteBuffer(R__b,this);
527    }
528 }
529
530
531
532
533