]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/PHOS/ana/IMT/AliHLTPHOSIsolatedMipTrigger.cxx
removing obsolete classes from build system
[u/mrichter/AliRoot.git] / HLT / PHOS / ana / IMT / AliHLTPHOSIsolatedMipTrigger.cxx
1 // $Id$
2
3 /**************************************************************************
4  * This file is property of and copyright by the Experimental Nuclear     *
5  * Physics Group, Dep. of Physics                                         *
6  * University of Oslo, Norway, 2007                                       *
7  *                                                                        *
8  * Author: Per Thomas Hille <perthi@fys.uio.no> for the ALICE HLT Project.*
9  * Contributors are mentioned in the code where appropriate.              *
10  * Please report bugs to perthi@fys.uio.no                                *
11  *                                                                        *
12  * Permission to use, copy, modify and distribute this software and its   *
13  * documentation strictly for non-commercial purposes is hereby granted   *
14  * without fee, provided that the above copyright notice appears in all   *
15  * copies and that both the copyright notice and this permission notice   *
16  * appear in the supporting documentation. The authors make no claims     *
17  * about the suitability of this software for any purpose. It is          *
18  * provided "as is" without express or implied warranty.                  *
19  **************************************************************************/
20 #include "AliHLTPHOSIsolatedMipTrigger.h"
21 #include "stdio.h"
22 #include "AliHLTPHOSDigit.h"
23 #include "TTree.h"
24 #include "TChain.h"
25 #include "TClonesArray.h"
26 #include "TH1F.h"
27 #include "TFile.h"
28 #include <iostream>
29 #include "AliHLTPHOSMipCluster.h"
30 #include "AliHLTPHOSMipClusterManager.h"
31
32
33 #define Q_MIN 27
34 #define Q_MAX 33
35
36 #define XBIN_LOW  0
37 #define XBIN_UP   1023
38 #define N_BINS    1023
39
40 AliHLTPHOSIsolatedMipTrigger::AliHLTPHOSIsolatedMipTrigger() :  fMipLowCut(3),  
41                                                                 fMipHighCut(150),
42                                                                 fModuleId(2)
43 {
44   fClusterManager = new AliHLTPHOSMipClusterManager();
45   Init();
46 }
47
48
49 AliHLTPHOSIsolatedMipTrigger::~AliHLTPHOSIsolatedMipTrigger()
50 {
51
52 }
53
54 /*
55 void
56 AliHLTPHOSIsolatedMipTrigger::SetFilePath(char *path)
57 {
58   sprintf(fCurrentFilePath, "%s", path);
59   printf("\n AliHLTPHOSIsolatedMipTrigger::SetFilePath, filepath was set to %s \n", fCurrentFilePath);
60 }
61 */
62
63 void
64 AliHLTPHOSIsolatedMipTrigger::SetMipLowCut(Float_t lCut)
65 {
66   fMipLowCut = lCut;
67 }
68
69 void
70 AliHLTPHOSIsolatedMipTrigger::SetMipHighCut(Float_t hCut)
71 {
72   fMipHighCut = hCut;
73 }
74
75
76 int
77 AliHLTPHOSIsolatedMipTrigger::ScanFileList()
78 {
79   int iRet = 0;
80   FILE *fp = fopen("filelist.txt", "r");
81   
82   if(fp != 0)
83     {
84       fscanf(fp, "%d\n", &fNFiles);
85       cout << "there are " << fNFiles << "  files " <<endl;
86  
87       for(int i=0; i< fNFiles; i++)
88         {
89           fscanf(fp, "%s\n"  ,fFileList[i]);
90         }
91       
92       for(int i=0; i< fNFiles; i++)
93         {
94           cout << fFileList[i] << endl;
95         }
96     }
97   else
98     {
99       cout << "ERROR opening file, filelist.txt" << endl;
100       iRet = -1;
101     } 
102   
103   return iRet;
104 }
105
106
107
108 int
109 AliHLTPHOSIsolatedMipTrigger::Analyze()
110 {
111   //  TClonesArray *digArray = new TClonesArray("AliHLTPHOSDigit" , 500); 
112
113   ScanFileList();
114   cout << "AliHLTPHOSIsolatedMipTrigger::Analyze fNFiles = " << fNFiles  <<endl;
115
116   for(int file=0; file < fNFiles; file ++)
117     {
118  
119      printf("\n AliHLTPHOSIsolatedMipTrigger::Analyze, analyzing file %s\n", fFileList[file]);   
120       
121       TFile *infile = infile->Open(fFileList[file],"read");
122       //      TFile *infile = infile->Open("/tmp/data2/digits151007/run8340_digitTree_5.root","read");
123
124       AliHLTPHOSDigit *tmpDigit = 0;  
125
126  
127       if(infile == 0)
128         {
129           cout << "ERROR could not open file" << endl;
130         }
131       else
132         {
133           
134           infile->cd();
135           cout << "File successfully opened, creating" << endl;
136           TChain *cain= new TChain("digitTree");
137           TClonesArray *digArray = new TClonesArray("AliHLTPHOSDigit" , 100);
138
139           cain->Add(fFileList[file]);
140
141           int nThrees =  cain->GetEntries() ;
142           cout << "there are "  <<  nThrees  <<  " entries in this file" << endl;
143           cain->SetBranchAddress("Digit", &digArray);
144
145           for(int i=0; i< nThrees; i++)
146             {
147               fClusterManager->ResetMipManager(); 
148               cain->GetEntry(i);
149               int nArrays = digArray->GetEntriesFast();
150
151               if(i%100 == 0)
152                 {
153                   cout << "processing event " << i << "  of file  " << fFileList[file] <<endl;
154                 }
155               
156          
157               for(int j=0; j < nArrays; j++)
158                 {
159
160                   tmpDigit = (AliHLTPHOSDigit*)digArray->At(j);
161                       
162                   if( IsMipCandidate(tmpDigit) == true)
163                     {
164                       fClusterManager->AddDigit(tmpDigit);
165                     }
166                 }
167               
168               FillClusterHistograms();
169             }
170  
171           delete digArray;
172           delete cain;
173           infile->Close(); 
174           delete infile;
175         }
176       
177       WriteHistograms();
178
179     }
180  
181   return 0;
182 }
183
184
185 bool
186 AliHLTPHOSIsolatedMipTrigger::IsMipCandidate(AliHLTPHOSDigit  *digit)
187 {
188   bool ret = false;
189   Float_t amplitude = digit->GetAmplitude();
190   Int_t *data =  digit->GetRawData();
191   Float_t sum = 0;
192   Float_t q;
193
194   if( amplitude > fMipLowCut && amplitude < fMipHighCut)
195     {
196       int start = digit->GetPreSamples();
197       int end = digit->GetTotalSamples() - start;
198  
199       for(int  i = start; i < end; i++)
200         {
201           sum+=data[i];
202         }
203       
204       q = ((Float_t)sum)/amplitude;
205
206       if( (q > Q_MIN) && (q < Q_MAX))
207         {
208           ret = true;
209         }
210     }
211   return ret;
212
213 }
214
215
216
217 void
218 AliHLTPHOSIsolatedMipTrigger::Init()
219 {
220   
221   char tmpHistoName[256];
222
223   for(int z = 0; z < N_ZROWS_MOD; z ++)
224     {
225       for(int x = 0; x < N_XCOLUMNS_MOD; x ++)
226         {
227           sprintf(tmpHistoName, "ClusterEnergies3x3_%d_%d_%d_%d",(int)fModuleId,  z, x, 1);
228           fClusterHistograms[z][x] = new TH1F( tmpHistoName, tmpHistoName, N_BINS, XBIN_LOW, XBIN_UP);
229         }
230     } 
231 }
232
233
234
235 void 
236 AliHLTPHOSIsolatedMipTrigger::FillClusterHistograms()
237 {
238   int cnt = 0;
239   AliHLTPHOSMipCluster* tmp;
240
241   for(int z = 0; z < N_ZROWS_MOD; z ++)
242     {
243       for(int x = 0; x < N_XCOLUMNS_MOD; x ++)
244         {
245           tmp = fClusterManager->GetCluster(z, x);
246           if(tmp != 0)
247             {
248               fClusterHistograms[z][x]->Fill(tmp->Get3x3Sum());
249               cnt ++;
250             }
251         }
252     } 
253   
254   //  cout << "filled  " << cnt <<"  clusters"  <<endl;
255 }
256
257 void
258 AliHLTPHOSIsolatedMipTrigger::WriteHistograms()
259 {
260   char tmpFileName[256];
261   sprintf(tmpFileName, "ClusterEnergies.root");
262
263   TFile *clusterFile =  new TFile(tmpFileName, "recreate");
264
265   for(int z = 0; z < N_ZROWS_MOD; z ++)
266     {
267       for(int x = 0; x < N_XCOLUMNS_MOD; x ++)
268         {
269           if(fClusterHistograms[z][x]->GetEntries() > 10 )
270             {
271               fClusterHistograms[z][x]->Write();
272             }
273         }
274     }  
275
276   clusterFile->Close();
277
278   delete clusterFile;
279 }