]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis/AliFMDAnalysisTaskDndeta.cxx
e9998a4b7fef1583bc2d9b084ebe3d88888a3ba1
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis / AliFMDAnalysisTaskDndeta.cxx
1  
2 #include <TROOT.h>
3 #include <TSystem.h>
4 #include <TInterpreter.h>
5 #include <TChain.h>
6 #include <TFile.h>
7 #include <TList.h>
8 #include <iostream>
9 #include "TH1F.h"
10 #include "TH2F.h"
11 #include "AliFMDAnalysisTaskDndeta.h"
12 #include "AliAnalysisManager.h"
13 #include "AliESDFMD.h"
14 #include "AliESDEvent.h"
15 #include "AliAODEvent.h"
16 #include "AliAODHandler.h"
17 #include "AliMCEventHandler.h"
18 #include "AliStack.h"
19 #include "AliLog.h"
20 #include "AliESDVertex.h"
21 #include "TMath.h"
22 #include "AliFMDAnaParameters.h"
23 //#include "AliFMDGeometry.h"
24 #include "AliGenEventHeader.h"
25 #include "AliGenPythiaEventHeader.h"
26 #include "AliHeader.h"
27 //#include "TDatabasePDG.h"
28 //#include "TParticlePDG.h"
29 #include "AliFMDStripIndex.h"
30 #include "AliESDInputHandler.h"
31 #include "AliGenDPMjetEventHeader.h"
32 #include "AliLog.h"
33 ClassImp(AliFMDAnalysisTaskDndeta)
34
35
36 AliFMDAnalysisTaskDndeta::AliFMDAnalysisTaskDndeta()
37 : fDebug(0),
38   fOutputList(0),
39   fInputList(0),
40   fVertexString(0x0),
41   fNevents(),
42   fNNSDevents(),
43   fNMCevents(),
44   fNMCNSDevents(),
45   fStandalone(kTRUE),
46   fLastTrackByStrip(0),
47   fVtxEff(1),
48   fVtxEffNSD(1)
49 {
50   // Default constructor
51   DefineInput (0, TList::Class());
52   DefineOutput(0, TList::Class());
53 }
54 //_____________________________________________________________________
55 AliFMDAnalysisTaskDndeta::AliFMDAnalysisTaskDndeta(const char* name, Bool_t SE):
56     AliAnalysisTask(name, "Density"),
57     fDebug(0),
58     fOutputList(0),
59     fInputList(0),
60     fVertexString(0x0),
61     fNevents(),
62     fNNSDevents(),
63     fNMCevents(),
64     fNMCNSDevents(),
65     fStandalone(kTRUE),
66     fLastTrackByStrip(0),
67     fVtxEff(1),
68     fVtxEffNSD(1)
69 {
70   fStandalone = SE;
71   if(fStandalone) {
72     DefineInput (0, TList::Class());
73     DefineInput(1, TObjString::Class());
74     DefineOutput(0, TList::Class());
75     
76   }
77 }
78 //_____________________________________________________________________
79 void AliFMDAnalysisTaskDndeta::CreateOutputObjects()
80 {
81   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
82   
83   if(!fOutputList)
84     fOutputList = new TList();
85   fOutputList->SetName("BackgroundCorrected");
86   
87   
88   TH2F* hMult = 0;
89   TH2F* hMultNSD = 0;
90   TH1F* hHits = 0;
91   TH2F* hMultTrVtx = 0;
92   TH1F* hPrimVertexBin = 0;
93   TH1F* hPrimVertexBinNSD = 0;
94   
95   TH2F* hBgTmp   = pars->GetBackgroundCorrection(1, 'I', 0);
96   TH1F* hPrimary = new TH1F("hMultvsEta","hMultvsEta",
97                             hBgTmp->GetNbinsX(),
98                             hBgTmp->GetXaxis()->GetXmin(),
99                             hBgTmp->GetXaxis()->GetXmax());
100   hPrimary->Sumw2();
101   fOutputList->Add(hPrimary);
102   TH1F* hPrimaryNSD = new TH1F("hMultvsEtaNSD","hMultvsEtaNSD",
103                             hBgTmp->GetNbinsX(),
104                             hBgTmp->GetXaxis()->GetXmin(),
105                             hBgTmp->GetXaxis()->GetXmax());
106   hPrimaryNSD->Sumw2();
107   fOutputList->Add(hPrimaryNSD);
108   Int_t nVtxbins = pars->GetNvtxBins();
109   TH2F* hBg = 0;
110   for(Int_t i = 0; i< nVtxbins; i++) {
111      
112     for(Int_t det =1; det<=3;det++)
113       {
114         Int_t nRings = (det==1 ? 1 : 2);
115         for(Int_t ring = 0;ring<nRings;ring++)
116           {
117             Char_t ringChar = (ring == 0 ? 'I' : 'O');
118             Int_t  nSec     = (ring == 0 ? 20 : 40);
119             
120             
121             
122             hBg = pars->GetBackgroundCorrection(det, ringChar, i);
123             hMult  = new TH2F(Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,i),Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,i),
124                               hBg->GetNbinsX(),
125                               hBg->GetXaxis()->GetXmin(),
126                               hBg->GetXaxis()->GetXmax(),
127                               nSec, 0, 2*TMath::Pi());
128             hMultTrVtx  = new TH2F(Form("dNdetaTrVtx_FMD%d%c_vtxbin%d",det,ringChar,i),Form("dNdetaTrVtx_FMD%d%c_vtxbin%d",det,ringChar,i),
129                                    hBg->GetNbinsX(),
130                                    hBg->GetXaxis()->GetXmin(),
131                                    hBg->GetXaxis()->GetXmax(),
132                                    nSec, 0, 2*TMath::Pi());
133             hMultNSD  = new TH2F(Form("dNdetaNSD_FMD%d%c_vtxbin%d",det,ringChar,i),Form("dNdetaNSD_FMD%d%c_vtxbin%d",det,ringChar,i),
134                                    hBg->GetNbinsX(),
135                                    hBg->GetXaxis()->GetXmin(),
136                                    hBg->GetXaxis()->GetXmax(),
137                                    nSec, 0, 2*TMath::Pi());
138             hHits  = new TH1F(Form("hMCHits_FMD%d%c_vtxbin%d",det,ringChar,i),Form("hMCHits_FMD%d%c_vtxbin%d",det,ringChar,i),
139                               hBg->GetNbinsX(),
140                               hBg->GetXaxis()->GetXmin(),
141                               hBg->GetXaxis()->GetXmax());
142             hHits->Sumw2();
143             fOutputList->Add(hHits);
144             
145             hMult->Sumw2();
146             fOutputList->Add(hMult);
147
148             hMultTrVtx->Sumw2();
149             fOutputList->Add(hMultTrVtx);
150             
151             hMultNSD->Sumw2();
152             fOutputList->Add(hMultNSD);
153             
154           }
155       } 
156   }
157   
158   for(Int_t i = 0; i< nVtxbins; i++) {
159    
160     hPrimVertexBin = new TH1F(Form("primmult_vtxbin%d",i),
161                               Form("primmult_vtxbin%d",i),
162                               hBgTmp->GetNbinsX(),
163                               hBgTmp->GetXaxis()->GetXmin(),
164                               hBgTmp->GetXaxis()->GetXmax());
165     hPrimVertexBin->Sumw2();
166     fOutputList->Add(hPrimVertexBin);
167     
168     hPrimVertexBinNSD = new TH1F(Form("primmult_NSD_vtxbin%d",i),
169                                  Form("primmult_NSD_vtxbin%d",i),
170                                  hBgTmp->GetNbinsX(),
171                                  hBgTmp->GetXaxis()->GetXmin(),
172                                  hBgTmp->GetXaxis()->GetXmax());
173     hPrimVertexBinNSD->Sumw2();
174     fOutputList->Add(hPrimVertexBinNSD);
175     
176
177     //SPD part
178     TH2F* hSPDMult  = new TH2F(Form("dNdeta_SPD_vtxbin%d",i),Form("dNdeta_SPD_vtxbin%d",i),
179                          hBgTmp->GetNbinsX(),
180                          hBgTmp->GetXaxis()->GetXmin(),
181                          hBgTmp->GetXaxis()->GetXmax(),
182                          20, 0, 2*TMath::Pi());
183     hSPDMult->Sumw2();
184     fOutputList->Add(hSPDMult);
185     TH2F* hSPDMultTrVtx  = new TH2F(Form("dNdetaTrVtx_SPD_vtxbin%d",i),Form("dNdetaTrVtx_SPD_vtxbin%d",i),
186                          hBgTmp->GetNbinsX(),
187                          hBgTmp->GetXaxis()->GetXmin(),
188                          hBgTmp->GetXaxis()->GetXmax(),
189                          20, 0, 2*TMath::Pi());
190     hSPDMultTrVtx->Sumw2();
191     fOutputList->Add(hSPDMultTrVtx);
192
193     TH2F* hSPDMultNSD  = new TH2F(Form("dNdetaNSD_SPD_vtxbin%d",i),Form("dNdetaNSD_SPD_vtxbin%d",i),
194                          hBgTmp->GetNbinsX(),
195                          hBgTmp->GetXaxis()->GetXmin(),
196                          hBgTmp->GetXaxis()->GetXmax(),
197                          20, 0, 2*TMath::Pi());
198     hSPDMultNSD->Sumw2();
199     fOutputList->Add(hSPDMultNSD);
200
201   }
202   
203   //AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
204   
205   //  TH2F* dNdetadphiHistogramTotal = new TH2F("dNdetadphiHistogramTotal","dNdetadphiHistogram;#eta;#Phi",pars->GetNetaBins(),-6,6,20,0,2*TMath::Pi());
206   //dNdetadphiHistogramTotal->SetErrorOption("g");
207   //fOutputList->Add(dNdetadphiHistogramTotal);
208   
209   
210   
211   fNevents.SetBins(nVtxbins,0,nVtxbins);
212   fNevents.SetName("nEvents");
213   fNNSDevents.SetBins(nVtxbins,0,nVtxbins);
214   fNNSDevents.SetName("nNSDEvents");
215   
216   fNMCevents.SetBins(nVtxbins,0,nVtxbins);
217   fNMCevents.SetName("nMCEvents");
218   fNMCNSDevents.SetBins(nVtxbins,0,nVtxbins);
219   fNMCNSDevents.SetName("nMCNSDEvents");
220   fOutputList->Add(&fNevents);
221   fOutputList->Add(&fNNSDevents);
222   fOutputList->Add(&fNMCevents);
223   fOutputList->Add(&fNMCNSDevents);
224 }
225 //_____________________________________________________________________
226 void AliFMDAnalysisTaskDndeta::ConnectInputData(Option_t */*option*/)
227 {
228   if(fStandalone) {
229     fInputList   = (TList*)GetInputData(0);
230   }
231 }
232 //_____________________________________________________________________
233 void AliFMDAnalysisTaskDndeta::Exec(Option_t */*option*/)
234 {
235   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
236   fVertexString = (TObjString*)fInputList->At(0);
237   
238   Int_t vtxbin   = fVertexString->GetString().Atoi();
239   
240   fNevents.Fill(vtxbin);
241   
242   //AliESDInputHandler* eventHandler = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
243   //AliESDEvent* esd = eventHandler->GetEvent();
244   Bool_t nsd = pars->IsEventTriggered(AliFMDAnaParameters::kNSD);
245   if(nsd) fNNSDevents.Fill(vtxbin);
246   
247   for(UShort_t det=1;det<=3;det++) {
248     Int_t nRings = (det==1 ? 1 : 2);
249     for (UShort_t ir = 0; ir < nRings; ir++) {
250       Char_t ringChar = (ir == 0 ? 'I' : 'O');
251             
252       TH2F* hMultTotal      = (TH2F*)fOutputList->FindObject(Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
253       TH2F* hMultTotalTrVtx = (TH2F*)fOutputList->FindObject(Form("dNdetaTrVtx_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
254       TH2F* hMultTotalNSD   = (TH2F*)fOutputList->FindObject(Form("dNdetaNSD_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
255       
256       
257       TH2F* hMultInput      = (TH2F*)fInputList->FindObject(Form("mult_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
258       TH2F* hMultInputTrVtx = (TH2F*)fInputList->FindObject(Form("multTrVtx_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
259       TH2F* hMultInputNSD   = (TH2F*)fInputList->FindObject(Form("multNSD_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
260       
261       hMultTotal->Add(hMultInput);
262       hMultTotalTrVtx->Add(hMultInputTrVtx);
263       if(nsd)
264         hMultTotalNSD->Add(hMultInputNSD);
265       
266     }
267   }
268   
269   //SPD code
270   TH2F* hMultSPDTotal      = (TH2F*)fOutputList->FindObject(Form("dNdeta_SPD_vtxbin%d",vtxbin));
271   TH2F* hMultSPDTotalTrVtx = (TH2F*)fOutputList->FindObject(Form("dNdetaTrVtx_SPD_vtxbin%d",vtxbin));
272   TH2F* hMultSPDTotalNSD   = (TH2F*)fOutputList->FindObject(Form("dNdetaNSD_SPD_vtxbin%d",vtxbin));
273   TH2F* hMultSPDInput      = (TH2F*)fInputList->FindObject(Form("mult_SPD_vtxbin%d",vtxbin));
274   TH2F* hMultSPDInputTrVtx = (TH2F*)fInputList->FindObject(Form("multTrVtx_SPD_vtxbin%d",vtxbin));
275   TH2F* hMultSPDInputNSD   = (TH2F*)fInputList->FindObject(Form("multNSD_SPD_vtxbin%d",vtxbin));
276   
277   hMultSPDTotal->Add(hMultSPDInput);
278   hMultSPDTotalTrVtx->Add(hMultSPDInputTrVtx);
279   if(nsd)
280     hMultSPDTotalNSD->Add(hMultSPDInputNSD);
281   
282   if(pars->GetProcessPrimary())
283     ProcessPrimary();
284   
285   //TH2F* dNdetadphiHistogram = (TH2F*)fOutputList->FindObject("dNdetadphiHistogramSPDTrVtx");
286   
287   //  TH2F* dNdetadphiHistogramTotal = (TH2F*)fOutputList->FindObject("dNdetadphiHistogramTotal");
288   
289   //  if(vtxbin == 4)
290   //  dNdetadphiHistogramTotal->Add(dNdetadphiHistogram);
291   
292   
293   
294   
295   if(fStandalone) {
296     PostData(0, fOutputList); 
297   }
298   
299 }
300 //_____________________________________________________________________
301 void AliFMDAnalysisTaskDndeta::Terminate(Option_t */*option*/) {
302   
303   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
304   
305   Int_t nVtxbins = pars->GetNvtxBins();
306   
307   for(UShort_t det=1;det<=3;det++) {
308     Int_t nRings = (det==1 ? 1 : 2);
309     for (UShort_t ir = 0; ir < nRings; ir++) {
310       Char_t ringChar = (ir == 0 ? 'I' : 'O');
311       for(Int_t i =0; i<nVtxbins; i++) {
312         
313         TH2F* hMultTotal = (TH2F*)fOutputList->FindObject(Form("dNdeta_FMD%d%c_vtxbin%d",det,ringChar,i));
314         fOutputList->Add(hMultTotal->Clone(Form("%s_orig", hMultTotal->GetName())));
315         if(fVtxEff) { 
316           hMultTotal->Scale(fVtxEff);
317         }
318         
319         TH2F* hMultTotalNSD = (TH2F*)fOutputList->FindObject(Form("dNdetaNSD_FMD%d%c_vtxbin%d",det,ringChar,i));
320         fOutputList->Add(hMultTotalNSD->Clone(Form("%s_orig", hMultTotalNSD->GetName())));
321         if(fVtxEffNSD) { 
322           hMultTotalNSD->Scale(fVtxEffNSD);
323         }
324         
325         //TH2F* hMultTrVtx = (TH2F*)hMultTotal->Clone(Form("dNdeta_FMD%d%c_TrVtx_vtxbin%d",det,ringChar,i));
326         TH2F* hMultTrVtx = (TH2F*)fOutputList->FindObject(Form("dNdetaTrVtx_FMD%d%c_vtxbin%d",det,ringChar,i));
327         
328         TH1D* hMultProj   = hMultTotal->ProjectionX(Form("dNdeta_FMD%d%c_vtxbin%d_proj",det,ringChar,i),1,hMultTotal->GetNbinsY());
329         TH1D* hMultProjTrVtx   = hMultTrVtx->ProjectionX(Form("dNdeta_FMD%d%c_TrVtx_vtxbin%d_proj",det,ringChar,i),1,hMultTrVtx->GetNbinsY());
330         TH1D* hMultProjNSD   = hMultTotalNSD->ProjectionX(Form("dNdetaNSD_FMD%d%c_vtxbin%d_proj",det,ringChar,i),1,hMultTotalNSD->GetNbinsY());
331         
332         //fOutputList->Add(hMultTrVtx);
333         fOutputList->Add(hMultProj);
334         fOutputList->Add(hMultProjTrVtx);
335         fOutputList->Add(hMultProjNSD);
336       }
337     }
338   }
339   
340   for(Int_t i =0; i<nVtxbins; i++) {
341     
342     TH2F* hSPDMult = (TH2F*)fOutputList->FindObject(Form("dNdeta_SPD_vtxbin%d",i));
343     TH2F* hSPDMultTrVtx = (TH2F*)fOutputList->FindObject(Form("dNdetaTrVtx_SPD_vtxbin%d",i));
344     TH2F* hSPDMultNSD = (TH2F*)fOutputList->FindObject(Form("dNdetaNSD_SPD_vtxbin%d",i));
345     
346     if(fVtxEff)
347       hSPDMult->Scale(fVtxEff);
348         
349     if(fVtxEffNSD)
350       hSPDMultNSD->Scale(fVtxEffNSD);
351     
352     TH1D* hMultProj   = hSPDMult->ProjectionX(Form("dNdeta_SPD_vtxbin%d_proj",i),1,hSPDMult->GetNbinsY());
353     TH1D* hMultProjTrVtx   = hSPDMultTrVtx->ProjectionX(Form("dNdetaTrVtx_SPD_vtxbin%d_proj",i),1,hSPDMultTrVtx->GetNbinsY());
354     TH1D* hMultProjNSD   = hSPDMultNSD->ProjectionX(Form("dNdetaNSD_SPD_vtxbin%d_proj",i),1,hSPDMultNSD->GetNbinsY());
355     fOutputList->Add(hMultProj);
356     fOutputList->Add(hMultProjTrVtx);
357     fOutputList->Add(hMultProjNSD);
358   
359   }
360   
361   std::cout<<"FMD analysis accepted "<<fNevents.GetEntries()<<" INEL events and "<<fNNSDevents.GetEntries()<<" NSD events"<<std::endl;
362 }
363 //_____________________________________________________________________
364 void AliFMDAnalysisTaskDndeta::ProcessPrimary() {
365   
366   AliMCEventHandler* eventHandler = 
367     dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()
368                                       ->GetMCtruthEventHandler());
369   if (!eventHandler) return;
370
371   AliMCEvent* mcEvent = eventHandler->MCEvent();
372   if(!mcEvent) return;
373   
374   fLastTrackByStrip.Reset(-1);
375   
376   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
377   
378   AliMCParticle* particle = 0;
379   AliStack* stack = mcEvent->Stack();
380   
381   TH1F* hPrimary = (TH1F*)fOutputList->FindObject("hMultvsEta");
382   TH1F* hPrimaryNSD = (TH1F*)fOutputList->FindObject("hMultvsEtaNSD");
383   AliHeader* header            = mcEvent->Header();
384   AliGenEventHeader* genHeader = header->GenEventHeader();
385   
386   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
387   AliGenDPMjetEventHeader* dpmHeader = dynamic_cast<AliGenDPMjetEventHeader*>(header->GenEventHeader());
388   Bool_t nsd = kTRUE;
389   
390   if (!pythiaGenHeader && !dpmHeader) {
391     std::cout<<" no pythia or dpm header!"<<std::endl;
392   }
393   else {
394     if(pythiaGenHeader) {
395       Int_t pythiaType = pythiaGenHeader->ProcessType();
396       
397       if(pythiaType==92||pythiaType==93)
398         nsd = kFALSE;
399       
400     }
401     if(dpmHeader) {
402       Int_t processType = dpmHeader->ProcessType();
403       if(processType == 5 || processType == 6)  
404         nsd = kFALSE;
405     
406     }
407   }
408   
409   TArrayF vertex;
410   genHeader->PrimaryVertex(vertex);
411   if(TMath::Abs(vertex.At(2)) > pars->GetVtxCutZ())
412     return;
413   Double_t delta           = 2*pars->GetVtxCutZ()/pars->GetNvtxBins();
414   Double_t vertexBinDouble = (vertex.At(2) + pars->GetVtxCutZ()) / delta;
415   Int_t    vertexBin       = (Int_t)vertexBinDouble;
416     
417   Bool_t firstTrack = kTRUE;
418   Bool_t firstTrackNSD = kTRUE;
419   
420   TH1F* hPrimVtxBinNSD = (TH1F*)fOutputList->FindObject(Form("primmult_NSD_vtxbin%d",vertexBin));
421   TH1F* hPrimVtxBin = (TH1F*)fOutputList->FindObject(Form("primmult_vtxbin%d",vertexBin));
422   
423   // we loop over the primaries only unless we need the hits (diagnostics running slowly)
424   Int_t nTracks = stack->GetNprimary();
425   if(pars->GetProcessHits())
426     nTracks = stack->GetNtrack();
427   
428   for(Int_t i = 0 ;i<nTracks;i++) {
429     particle = (AliMCParticle*) mcEvent->GetTrack(i);
430     if(!particle)
431       continue;
432    
433     if(stack->IsPhysicalPrimary(i) && particle->Charge() != 0) {
434       hPrimary->Fill(particle->Eta());
435       hPrimVtxBin->Fill(particle->Eta());
436       if(firstTrack) {
437         fNMCevents.Fill(vertexBin);
438         firstTrack = kFALSE;
439       }
440       
441       if(nsd) {
442         hPrimaryNSD->Fill(particle->Eta());
443         hPrimVtxBinNSD->Fill(particle->Eta());
444         if(firstTrackNSD) {
445           fNMCNSDevents.Fill(vertexBin);
446           firstTrackNSD = kFALSE;
447         }
448       }
449       
450       
451     
452     }
453     if(pars->GetProcessHits()) {
454            
455       for(Int_t j=0; j<particle->GetNumberOfTrackReferences();j++) {
456         
457         AliTrackReference* ref = particle->GetTrackReference(j);
458         UShort_t det,sec,strip;
459         Char_t   ring;
460         if(ref->DetectorId() != AliTrackReference::kFMD)
461           continue;
462         AliFMDStripIndex::Unpack(ref->UserId(),det,ring,sec,strip);
463         Float_t thisStripTrack = fLastTrackByStrip(det,ring,sec,strip);
464         if(particle->Charge() != 0 && i != thisStripTrack ) {
465           //Double_t x,y,z;
466           
467           Float_t   eta   = pars->GetEtaFromStrip(det,ring,sec,strip,vertex.At(2));//-1*TMath::Log(TMath::Tan(0.5*theta));
468           TH1F* hHits = (TH1F*)fOutputList->FindObject(Form("hMCHits_FMD%d%c_vtxbin%d",det,ring,vertexBin));
469           
470         
471           hHits->Fill(eta);
472           
473           Float_t nstrips = (ring =='O' ? 256 : 512);
474           
475           fLastTrackByStrip(det,ring,sec,strip) = (Float_t)i;
476         
477           if(strip >0)
478             fLastTrackByStrip(det,ring,sec,strip-1) = (Float_t)i;
479           if(strip < (nstrips - 1))
480             fLastTrackByStrip(det,ring,sec,strip+1) = (Float_t)i;
481           
482         }
483       
484         
485       }
486       
487       
488     }
489
490       
491       
492   }
493 }
494 //_____________________________________________________________________
495 //
496 //
497 // EOF
498 // EOF