]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/PHOS/AliHLTPHOSHistoProdInvMass.cxx
A fast HLT version of the SPD clusterfinder implemented.
[u/mrichter/AliRoot.git] / HLT / PHOS / AliHLTPHOSHistoProdInvMass.cxx
1
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        * 
4  * All rights reserved.                                                   *
5  *                                                                        *
6  * Primary Authors: Albin Gaignette
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          * 
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /** 
18  * @file   AliHLTPHOSHistoProdInvMass
19  * @author Albin Gaignette
20  * @date 
21  * @brief  Histogram producer for PHOS HLT 
22  */
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
29
30 #include "AliHLTPHOSHistoProdInvMass.h"
31 //#include "AliESDCaloCluster.h"
32 #include "TMath.h"
33
34 #include "AliHLTCaloClusterDataStruct.h"
35 #include "AliHLTCaloClusterReader.h"
36 #include "TObjArray.h"
37 //#include "TClonesArray.h"
38 //#include <iostream>
39 #include "TH1F.h"
40 //#include "TH2F.h"
41
42 AliHLTPHOSHistoProdInvMass::AliHLTPHOSHistoProdInvMass() :
43   fClusterReader(NULL),
44   fHistTwoClusterInvMass(0),
45   fHistArrayPtr(0)
46 {
47   // See header file for documentation
48   fHistArrayPtr = new TObjArray;
49   fClusterReader = new AliHLTCaloClusterReader();
50
51   fHistTwoClusterInvMass = new TH1F("fHistTwoClusterInvMass", "Invariant mass of two clusters", 200, 0, 1);
52   fHistTwoClusterInvMass->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
53   fHistTwoClusterInvMass->GetYaxis()->SetTitle("Number of counts");
54   fHistTwoClusterInvMass->SetMarkerStyle(21);
55   fHistArrayPtr->AddLast(fHistTwoClusterInvMass);
56 }
57
58 AliHLTPHOSHistoProdInvMass::~AliHLTPHOSHistoProdInvMass()
59 {
60   if(fHistTwoClusterInvMass)
61     {
62       delete fHistTwoClusterInvMass;
63       fHistTwoClusterInvMass = 0;
64     }
65 }
66
67 TObjArray* AliHLTPHOSHistoProdInvMass::GetHistograms()
68 {  
69   // See header file for documentation
70
71   return fHistArrayPtr;
72 }
73
74 Int_t AliHLTPHOSHistoProdInvMass::DoEvent(AliHLTCaloClusterHeaderStruct* cHeader) {   
75   
76   fClusterReader->SetMemory(cHeader);
77   
78   int ncls = cHeader->fNClusters;
79   Float_t* cPos[ncls];
80   Float_t cEnergy[ncls];
81   
82   AliHLTCaloClusterDataStruct* cluster;
83   Int_t icls = 0;
84   while ( ( cluster = fClusterReader->NextCluster() ) ) {
85     
86     cPos[icls] = cluster->fGlobalPos;
87     cEnergy[icls] = cluster->fEnergy; 
88     
89     icls++;
90   }  
91   
92   for(Int_t ipho = 0; ipho<(ncls-1); ipho++) { 
93     for(Int_t jpho = ipho+1 ; jpho<ncls ; jpho++) { 
94       // Calcul of the theta angle between two photons
95       Double_t theta = (2* asin(0.5*TMath::Sqrt((cPos[ipho][0]-cPos[jpho][0])*(cPos[ipho][0]-cPos[jpho][0]) +(cPos[ipho][1]-cPos[jpho][1])*(cPos[ipho][1]-cPos[jpho][1]))/460));
96       
97       // Calcul of the mass m of the pion 
98       Double_t m =(TMath::Sqrt(2 * cEnergy[ipho]* cEnergy[jpho]*(1-TMath::Cos(theta))));
99       
100       fHistTwoClusterInvMass->Fill(m);
101     }
102   }
103   
104   return 0;
105 }
106   
107  
108  
109