]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliInvmass.h
Introduction of the reference to Copyright and cvs Id
[u/mrichter/AliRoot.git] / RALICE / AliInvmass.h
1 #ifndef ALIINVMASS_H
2 #define ALIINVMASS_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 ////////////////////////////////////////////////////////////////////////////////
9 // Class AliInvmass
10 // Construction of invariant mass and combinatorial background.
11 //
12 // Example :
13 // ---------
14 //
15 // TObjArray* photons=new TObjArray(); // Array with photon tracks for pi0 rec.
16 //
17 // // Code to create some photon tracks from pi0 decays
18 // Int_t ntracks=200;
19 // for (Int_t i=0; i<ntracks; i++)
20 // {
21 //  photons->Add(new Alitrack);
22 //  ...
23 //  ...
24 //  ...
25 // }  
26 //
27 // // Perform the invariant mass and comb. bkg. reconstruction
28 //
29 // TObjArray* allm=q.Invmass(photons,photons); // All reconstructed invariant masses
30 //
31 // TH1F* hall=new TH1F("hall","hall",200,0,2); // Histo with M_inv of all combinations
32 //
33 // Int_t nall=0;
34 // if (allm) nall=allm->GetEntries();
35 //
36 // AliTrack* t;
37 // Float_t minv;
38 // for (Int_t j=0; j<nall; j++)
39 // {
40 //  t=(AliTrack*)allm->At(j);
41 //  if (t)
42 //  {
43 //   minv=t->GetMass();
44 //   hall->Fill(minv);
45 //  }
46 // }
47 //
48 // TObjArray* bkgm=q.CombBkg(photons,photons); // Reconstructed comb. background
49 //
50 // TH1F* hbkg=new TH1F("hbkg","hbkg",200,0,2); // Histo with M_inv. of comb. background
51 //
52 // Int_t nbkg=0;
53 // if (bkgm) nbkg=bkgm->GetEntries();
54 //
55 // for (Int_t j=0; j<nbkg; j++)
56 // {
57 //  t=(AliTrack*)bkgm->At(j);
58 //  if (t)
59 //  {
60 //   minv=t->GetMass();
61 //   hbkg->Fill(minv);
62 //  }
63 // }
64 //
65 // TH1F* hsig=new TH1F("sig","sig",200,0,2);   // Histo with the bkg. subtracted signal
66 // hsig->Sumw2();
67 // hsig->Add(hall,hbkg,1,-1);
68 //
69 //
70 // Note : By default the storage of the reconstructed information is performed
71 //        in separate TObjArrays for the signal and comb. background resp.
72 //        In order to limit the memory usage, AliInvmass::SetStorageMode(1) may be
73 //        used to activate only a single TObjArray to store the reconstructed information.
74 //        Consequently, the following statements 
75 //
76 //         TObjArray* allm=q.Invmass(photons,photons);
77 //         TObjArray* bkgm=q.CombBkg(photons,photons);
78 //
79 //        will result in the fact that after he invokation of CombBkg
80 //        the information of "allm" is lost due to the fact that the storage is
81 //        is re-used for "bkgm" in case the "single storage" option has been selected.
82 //        Usage of the, in that case invalid, pointer "allm" may cause your
83 //        program to crash.
84 //
85 //        * Thus : In case of single storage usage, all invokations of the returned
86 //                 array pointer have to be completed before invoking any memberfunction
87 //                 of the same AliInvmass object again.
88 //        
89 //        
90 //
91 //--- NvE 12-apr-1999 UU-SAP Utrecht
92 ////////////////////////////////////////////////////////////////////////////////
93
94 #include <iostream.h>
95 #include <math.h>
96  
97 #include "TObject.h"
98 #include "TObjArray.h"
99
100 #include "AliRandom.h"
101 #include "AliTrack.h"
102
103 class AliInvmass : public TObject
104 {
105  public:
106   AliInvmass();                                    // Default constructor
107   ~AliInvmass();                                   // Destructor
108   void SetStorageMode(Int_t m);                    // Set storage mode (1=single, 2=multiple)
109   void SetThetaSwitch(Int_t i=1);                  // Enable (1/0) new theta for comb. bkg. reco.
110   void SetPhiSwitch(Int_t i=1);                    // Enable (1/0) new phi for comb. bkg. reco.
111   Int_t GetStorageMode();                          // Provide storage mode
112   Int_t GetThetaSwitch();                          // Provide theta switch flag
113   Int_t GetPhiSwitch();                            // Provide phi switch flag
114   TObjArray* Invmass(TObjArray* a1,TObjArray* a2); // Two-particle inv. mass reco.
115   TObjArray* CombBkg(TObjArray* a1,TObjArray* a2); // Two-particle comb. background reco.
116
117  protected:
118   Double_t fPi;     // Value of pi
119   Int_t fMode;      // Storage mode for signal and bkg. results (2=separate arrays)
120   Int_t fBkg;       // Flag to denote comb. background processing
121   AliRandom fRndm;  // The random number generator for the comb. bkg. reconstruction
122   Int_t fNewtheta;  // Flag to denote enabling of switching theta for comb. bkg. reco.
123   Int_t fNewphi;    // Flag to denote enabling of switching phi for comb. bkg. reco.
124   TObjArray* fMinv; // Array with reconstructed invariant mass 'tracks'
125   TObjArray* fMbkg; // Array with reconstructed comb. background 'tracks'
126
127  private:
128   void Combine(TObjArray* a1,TObjArray* a2); // Make two-particle combinations
129
130  ClassDef(AliInvmass,1) // Class definition to enable ROOT I/O
131 };
132 #endif