]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/physics/AliHLTCaloHistoInvMass.cxx
-Moved PHOS Physics histogram producers to HLT/global/physics
[u/mrichter/AliRoot.git] / HLT / global / physics / AliHLTCaloHistoInvMass.cxx
1 //-*- Mode: C++ -*-
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        * 
4  * All rights reserved.                                                   *
5  *                                                                        *
6  * Primary Authors: Albin Gaignette, Svein Lindal                         *
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   AliHLTCaloHistoInvMass
19  * @author Albin Gaignette rewritten Svein Lindal <slindal@fys.uio.no>
20  * @date 
21  * @brief  
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 "AliHLTCaloHistoInvMass.h"
31 #include "AliESDCaloCluster.h"
32 #include "TMath.h"
33
34 // #include "AliHLTCaloClusterDataStruct.h"
35 // #include "AliHLTCaloClusterReader.h"
36 #include "TObjArray.h"
37 #include "AliESDEvent.h"
38 #include "AliESDCaloCluster.h"
39 #include "TRefArray.h"
40 #include "TH1F.h"
41 #include "TString.h"
42
43 AliHLTCaloHistoInvMass::AliHLTCaloHistoInvMass(TString det) :
44   fHistTwoClusterInvMass(NULL),
45   fHistArrayPtr(NULL)
46 {
47   // See header file for documentation
48
49   fHistArrayPtr = new TObjArray;
50   
51   fHistTwoClusterInvMass = new TH1F(Form("%s fHistTwoClusterInvMass", det.Data()), Form("%s Invariant mass of two clusters PHOS", det.Data()), 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
59 AliHLTCaloHistoInvMass::~AliHLTCaloHistoInvMass()
60 {
61   if(fHistTwoClusterInvMass)
62     delete fHistTwoClusterInvMass;
63   fHistTwoClusterInvMass = NULL;
64  
65   if(fHistArrayPtr)
66     delete fHistArrayPtr;
67   fHistArrayPtr = NULL;
68
69 }
70
71 TObjArray* AliHLTCaloHistoInvMass::GetHistograms()
72 {  
73   // See header file for documentation
74   return fHistArrayPtr;
75 }
76
77 Int_t AliHLTCaloHistoInvMass::FillHistograms(Int_t nc, TRefArray * clustersArray) {
78   //See header file for documentation
79   
80   Float_t cPos[nc][3];
81   Float_t cEnergy[nc];
82
83   for(int ic = 0; ic < nc; ic++) {
84     AliESDCaloCluster * cluster = static_cast<AliESDCaloCluster*>(clustersArray->At(ic));
85     cluster->GetPosition(cPos[ic]);
86     cEnergy[ic] = cluster->E();
87   }
88
89   for(Int_t ipho = 0; ipho<(nc-1); ipho++) { 
90     for(Int_t jpho = ipho+1; jpho<nc; jpho++) { 
91       // Calculate the theta angle between two photons
92       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));
93       
94       // Calculate the mass m of the pion candidate
95       Double_t m =(TMath::Sqrt(2 * cEnergy[ipho]* cEnergy[jpho]*(1-TMath::Cos(theta))));
96       
97       fHistTwoClusterInvMass->Fill(m);
98     }
99   }
100
101   return 0;
102 }
103
104 // Int_t AliHLTCaloHistoInvMass::DoEvent(AliHLTCaloClusterHeaderStruct* cHeader) {   
105   
106 //   fClusterReader->SetMemory(cHeader);
107   
108 //   int ncls = cHeader->fNClusters;
109 //   Float_t* cPos[ncls];
110 //   Float_t cEnergy[ncls];
111   
112 //   AliHLTCaloClusterDataStruct* cluster;
113 //   Int_t icls = 0;
114 //   while ( ( cluster = fClusterReader->NextCluster() ) ) {
115     
116 //     cPos[icls] = cluster->fGlobalPos;
117 //     cEnergy[icls] = cluster->fEnergy; 
118     
119 //     icls++;
120 //   }  
121   
122 //   for(Int_t ipho = 0; ipho<(ncls-1); ipho++) { 
123 //     for(Int_t jpho = ipho+1 ; jpho<ncls ; jpho++) { 
124 //       // Calcul of the theta angle between two photons
125 //       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));
126       
127 //       // Calcul of the mass m of the pion 
128 //       Double_t m =(TMath::Sqrt(2 * cEnergy[ipho]* cEnergy[jpho]*(1-TMath::Cos(theta))));
129       
130
131 //       //BALLE
132 //       // fHistTwoClusterInvMass->Fill(m);
133 //     }
134 //   }
135   
136 //   return 0;
137 // }
138   
139  
140  
141