]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALClusterizerv1.cxx
Correct the check of an already existing branch with the same name in Exec()
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALClusterizerv1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /* $Log:
19    1 October 2000. Yuri Kharlov:
20      AreNeighbours()
21      PPSD upper layer is considered if number of layers>1
22
23    18 October 2000. Yuri Kharlov:
24      AliEMCALClusterizerv1()
25      CPV clusterizing parameters added
26
27      MakeClusters()
28      After first PPSD digit remove EMC digits only once
29 */
30 //*-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
31 //////////////////////////////////////////////////////////////////////////////
32 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
33 //  unfolds the clusters having several local maxima.  
34 //  Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
35 //  EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all 
36 //  parameters including input digits branch title, thresholds etc.)
37 //  This TTask is normally called from Reconstructioner, but can as well be used in 
38 //  standalone mode.
39 // Use Case:
40 //  root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")  
41 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 //               //reads gAlice from header file "..."                      
43 //  root [1] cl->ExecuteTask()  
44 //               //finds RecPoints in all events stored in galice.root
45 //  root [2] cl->SetDigitsBranch("digits2") 
46 //               //sets another title for Digitis (input) branch
47 //  root [3] cl->SetRecPointsBranch("recp2")  
48 //               //sets another title four output branches
49 //  root [4] cl->SetTowerLocalMaxCut(0.03)  
50 //               //set clusterization parameters
51 //  root [5] cl->ExecuteTask("deb all time")  
52 //               //once more finds RecPoints options are 
53 //               // deb - print number of found rec points
54 //               // deb all - print number of found RecPoints and some their characteristics 
55 //               // time - print benchmarking results
56
57 // --- ROOT system ---
58
59 #include "TROOT.h" 
60 #include "TFile.h" 
61 #include "TFolder.h" 
62 #include "TMath.h" 
63 #include "TMinuit.h"
64 #include "TTree.h" 
65 #include "TSystem.h" 
66 #include "TBenchmark.h"
67
68 // --- Standard library ---
69
70 #include <iostream.h>
71 #include <iomanip.h>
72
73 // --- AliRoot header files ---
74
75 #include "AliEMCALClusterizerv1.h"
76 #include "AliEMCALDigit.h"
77 #include "AliEMCALDigitizer.h"
78 #include "AliEMCALTowerRecPoint.h"
79 #include "AliEMCAL.h"
80 #include "AliEMCALGetter.h"
81 #include "AliRun.h"
82
83 ClassImp(AliEMCALClusterizerv1)
84   
85 //____________________________________________________________________________
86   AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
87 {
88   // default ctor (to be used mainly by Streamer)
89   
90   fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ; 
91   
92   fPreShoClusteringThreshold = 0.0001;
93   fTowerClusteringThreshold  = 0.2;   
94   
95   fTowerLocMaxCut  = 0.03 ;
96   fPreShoLocMaxCut = 0.03 ;
97   
98   fW0     = 4.5 ;
99   fW0CPV  = 4.0 ;
100
101   fTimeGate = 1.e-8 ; 
102
103   fToUnfold = kFALSE ;
104
105   fHeaderFileName = fDigitsBranchTitle  = "" ;
106   fRecPointsInRun          = 0 ;   
107 }
108
109 //____________________________________________________________________________
110 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const char* headerFile,const char* name)
111 :AliEMCALClusterizer(headerFile, name)
112 {
113   // ctor with the indication of the file where header Tree and digits Tree are stored
114   
115
116   fNumberOfPreShoClusters = fNumberOfTowerClusters = 0 ; 
117
118  
119   
120   fPreShoClusteringThreshold  = 0.0001;
121   fTowerClusteringThreshold   = 0.2;   
122   
123   fTowerLocMaxCut  = 0.03 ;
124   fPreShoLocMaxCut = 0.03 ;
125   
126   fW0     = 4.5 ;
127   fW0CPV  = 4.0 ;
128
129   fTimeGate = 1.e-8 ; 
130   
131   fToUnfold = kFALSE ;
132   
133   fHeaderFileName     = GetTitle() ; 
134   fDigitsBranchTitle  = GetName() ;
135   
136   TString clusterizerName( GetName()) ; 
137   clusterizerName.Append(":") ; 
138   clusterizerName.Append(Version()) ; 
139   SetName(clusterizerName) ;
140   fRecPointsInRun          = 0 ; 
141
142   Init() ;
143
144 }
145 //____________________________________________________________________________
146   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
147 {
148 }
149 //____________________________________________________________________________
150 Float_t  AliEMCALClusterizerv1::Calibrate(Int_t amp, Bool_t inpresho) const
151 {
152   if ( inpresho ) // calibrate as pre shower
153      return -fADCpedestalPreSho + amp * fADCchannelPreSho ; 
154
155   else //calibrate as tower 
156     return -fADCpedestalTower + amp * fADCchannelTower ;                
157 }
158 //____________________________________________________________________________
159 void AliEMCALClusterizerv1::Exec(Option_t * option)
160 {
161   // Steering method
162
163   if( strcmp(GetName(), "")== 0 ) 
164     Init() ;
165
166   if(strstr(option,"tim"))
167     gBenchmark->Start("EMCALClusterizer"); 
168   
169   if(strstr(option,"print"))
170     Print("") ; 
171
172   gAlice->GetEvent(0) ;
173   
174  //check, if the branch with name of this" already exits?
175   TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
176   TIter next(lob) ; 
177   TBranch * branch = 0 ;  
178   Bool_t emcaltowerfound = kFALSE, emcalpreshofound = kFALSE, clusterizerfound = kFALSE ; 
179
180   TString branchname = GetName() ;
181   branchname.Remove(branchname.Index(Version())-1) ;
182   
183   while ( (branch = (TBranch*)next()) && (!emcaltowerfound || !emcalpreshofound || !clusterizerfound) ) {
184     if ( (strcmp(branch->GetName(), "EMCALTowerRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
185       emcaltowerfound = kTRUE ;
186     
187     else if ( (strcmp(branch->GetName(), "EMCALPreShoRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
188       emcalpreshofound = kTRUE ;
189    
190     else if ((strcmp(branch->GetName(), "AliEMCALClusterizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
191       clusterizerfound = kTRUE ; 
192   }
193
194   if ( emcalpreshofound || emcaltowerfound || clusterizerfound ) {
195     cerr << "WARNING: AliEMCALClusterizer::Exec -> Tower(PreSho)RecPoints and/or Clusterizer branch with name " 
196          << branchname.Data() << " already exits" << endl ;
197     return ; 
198   }       
199
200   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
201   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
202   Int_t ievent ;
203
204   for(ievent = 0; ievent < nevents; ievent++){
205
206     if(ievent == 0)
207       GetCalibrationParameters() ;
208
209     fNumberOfTowerClusters = fNumberOfPreShoClusters = 0 ;
210    
211     gime->Event(ievent,"D") ;
212         
213     MakeClusters() ;
214     
215     if(fToUnfold)
216       MakeUnfolding() ;
217
218     WriteRecPoints(ievent) ;
219
220     if(strstr(option,"deb"))  
221       PrintRecPoints(option) ;
222
223     //increment the total number of digits per run 
224     fRecPointsInRun += gime->TowerRecPoints()->GetEntriesFast() ;  
225     fRecPointsInRun += gime->PreShowerRecPoints()->GetEntriesFast() ;  
226  }
227   
228   if(strstr(option,"tim")){
229     gBenchmark->Stop("EMCALClusterizer");
230     cout << "AliEMCALClusterizer:" << endl ;
231     cout << "  took " << gBenchmark->GetCpuTime("EMCALClusterizer") << " seconds for Clusterizing " 
232          <<  gBenchmark->GetCpuTime("EMCALClusterizer")/nevents << " seconds per event " << endl ;
233     cout << endl ;
234   }
235   
236 }
237
238 //____________________________________________________________________________
239 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALTowerRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
240                                     Int_t nPar, Float_t * fitparameters) const
241
242   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
243   // The initial values for fitting procedure are set equal to the positions of local maxima.
244   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
245
246   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
247   TClonesArray * digits = gime->Digits() ; 
248   
249
250   gMinuit->mncler();                     // Reset Minuit's list of paramters
251   gMinuit->SetPrintLevel(-1) ;           // No Printout
252   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;  
253                                          // To set the address of the minimization function 
254
255   TList * toMinuit = new TList();
256   toMinuit->AddAt(emcRP,0) ;
257   toMinuit->AddAt(digits,1) ;
258   
259   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
260
261   // filling initial values for fit parameters
262   AliEMCALDigit * digit ;
263
264   Int_t ierflg  = 0; 
265   Int_t index   = 0 ;
266   Int_t nDigits = (Int_t) nPar / 3 ;
267
268   Int_t iDigit ;
269
270   AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
271
272   for(iDigit = 0; iDigit < nDigits; iDigit++){
273     digit = (AliEMCALDigit *) maxAt[iDigit]; 
274
275     Int_t relid[4] ;
276     Float_t x = 0.;
277     Float_t z = 0.;
278     geom->AbsToRelNumbering(digit->GetId(), relid) ;
279     geom->PosInAlice(relid, x, z) ;
280
281     Float_t energy = maxAtEnergy[iDigit] ;
282
283     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
284     index++ ;   
285     if(ierflg != 0){ 
286       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
287       return kFALSE;
288     }
289     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
290     index++ ;   
291     if(ierflg != 0){
292       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
293       return kFALSE;
294     }
295     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
296     index++ ;   
297     if(ierflg != 0){
298       cout << "EMCAL Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
299       return kFALSE;
300     }
301   }
302
303   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
304                       //  depends on it. 
305   Double_t p1 = 1.0 ;
306   Double_t p2 = 0.0 ;
307
308   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
309   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
310   gMinuit->SetMaxIterations(5);
311   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
312
313   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
314
315   if(ierflg == 4){  // Minimum not found   
316     cout << "EMCAL Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
317     return kFALSE ;
318   }            
319   for(index = 0; index < nPar; index++){
320     Double_t err ;
321     Double_t val ;
322     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
323     fitparameters[index] = val ;
324    }
325
326   delete toMinuit ;
327   return kTRUE;
328
329 }
330
331 //____________________________________________________________________________
332 void AliEMCALClusterizerv1::GetCalibrationParameters() 
333 {
334   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
335   TString branchname = GetName() ;
336   branchname.Remove(branchname.Index(Version())-1) ;
337
338   AliEMCALDigitizer * dig = gime->Digitizer(branchname) ;
339
340   fADCchannelTower   = dig->GetTowerchannel() ;
341   fADCpedestalTower  = dig->GetTowerpedestal();
342
343   fADCchannelPreSho  = dig->GetPreShochannel() ;
344   fADCpedestalPreSho = dig->GetPreShopedestal() ; 
345
346 }
347 //____________________________________________________________________________
348 void AliEMCALClusterizerv1::Init()
349 {
350   // Make all memory allocations which can not be done in default constructor.
351   // Attach the Clusterizer task to the list of EMCAL tasks
352   
353   if ( strcmp(GetTitle(), "") == 0 )
354     SetTitle("galice.root") ;
355
356   TString branchname = GetName() ;
357   branchname.Remove(branchname.Index(Version())-1) ;
358
359   AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), branchname, "update") ; 
360   if ( gime == 0 ) {
361     cerr << "ERROR: AliEMCALClusterizerv1::Init -> Could not obtain the Getter object !" << endl ; 
362     return ;
363   } 
364     
365   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
366   fNTowers = geom->GetNZ() *  geom->GetNPhi() ;
367
368   if(!gMinuit) 
369     gMinuit = new TMinuit(100) ;
370
371   gime->PostClusterizer(this) ;
372   // create a folder on the white board 
373   gime->PostRecPoints(branchname ) ;
374
375   gime->PostDigits(branchname) ;
376   gime->PostDigitizer(branchname) ;
377   
378 }
379
380 //____________________________________________________________________________
381 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
382 {
383   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
384   //                                       = 1 are neighbour
385   //                                       = 2 are not neighbour but do not continue searching
386   // neighbours are defined as digits having at least a common vertex 
387   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
388   //                                      which is compared to a digit (d2)  not yet in a cluster  
389
390    AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ;
391
392   Int_t rv = 0 ; 
393
394   Int_t relid1[4] ; 
395   geom->AbsToRelNumbering(d1->GetId(), relid1) ; 
396
397   Int_t relid2[4] ; 
398   geom->AbsToRelNumbering(d2->GetId(), relid2) ; 
399  
400   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same EMCAL Arm 
401     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
402     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
403     
404     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
405       if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fTimeGate))
406       rv = 1 ; 
407     }
408     else {
409       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
410         rv = 2; //  Difference in row numbers is too large to look further 
411     }
412
413   } 
414   else {
415     
416     if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )  
417       rv=2 ;
418
419   }
420
421   return rv ; 
422 }
423
424
425 //____________________________________________________________________________
426 Bool_t AliEMCALClusterizerv1::IsInTower(AliEMCALDigit * digit) const
427 {
428   // Tells if (true) or not (false) the digit is in a EMCAL-Tower 
429  
430   Bool_t rv = kFALSE ; 
431   if (!digit->IsInPreShower()) 
432     rv = kTRUE; 
433   return rv ; 
434 }
435
436 //____________________________________________________________________________
437 Bool_t AliEMCALClusterizerv1::IsInPreShower(AliEMCALDigit * digit) const
438 {
439   // Tells if (true) or not (false) the digit is in a EMCAL-PreShower
440  
441   Bool_t rv = kFALSE ; 
442   if (digit->IsInPreShower()) 
443     rv = kTRUE; 
444   return rv ; 
445 }
446
447 //____________________________________________________________________________
448 void AliEMCALClusterizerv1::WriteRecPoints(Int_t event)
449 {
450
451   // Creates new branches with given title
452   // fills and writes into TreeR.
453
454   TString branchName(GetName() ) ;
455   branchName.Remove(branchName.Index(Version())-1) ;
456
457   AliEMCALGetter *gime = AliEMCALGetter::GetInstance() ; 
458   TObjArray * towerRecPoints = gime->TowerRecPoints(branchName) ; 
459   TObjArray * preshoRecPoints = gime->PreShowerRecPoints(branchName) ; 
460   TClonesArray * digits = gime->Digits(branchName) ; 
461
462   Int_t index ;
463   //Evaluate position, dispersion and other RecPoint properties...
464   for(index = 0; index < towerRecPoints->GetEntries(); index++)
465     (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->EvalAll(fW0,digits) ;
466
467   towerRecPoints->Sort() ;
468
469   for(index = 0; index < towerRecPoints->GetEntries(); index++)
470     (dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(index)))->SetIndexInList(index) ;
471
472   towerRecPoints->Expand(towerRecPoints->GetEntriesFast()) ; 
473
474   //Now the same for CPV
475   for(index = 0; index < preshoRecPoints->GetEntries(); index++)
476     (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->EvalAll(fW0CPV,digits)  ;
477
478   preshoRecPoints->Sort() ;
479
480   for(index = 0; index < preshoRecPoints->GetEntries(); index++)
481     (dynamic_cast<AliEMCALRecPoint *>(preshoRecPoints->At(index)))->SetIndexInList(index) ;
482
483   preshoRecPoints->Expand(preshoRecPoints->GetEntriesFast()) ;
484   
485   //Make branches in TreeR for RecPoints and Clusterizer
486   char * filename = 0;
487   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
488     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
489     sprintf(filename,"%s/EMCAL.Reco.root",gAlice->GetBaseFile()) ;
490   }
491   
492   //Make new branches
493   TDirectory *cwd = gDirectory;
494   
495  
496   Int_t bufferSize = 32000 ;    
497   Int_t splitlevel = 0 ;
498
499   //First EMC
500   TBranch * emcBranch = gAlice->TreeR()->Branch("EMCALTowerRP","TObjArray",&towerRecPoints,bufferSize,splitlevel);
501   emcBranch->SetTitle(branchName);
502   if (filename) {
503     emcBranch->SetFile(filename);
504     TIter next( emcBranch->GetListOfBranches());
505     TBranch * sb ;
506     while ((sb=(TBranch*)next())) {
507       sb->SetFile(filename);
508     }   
509     
510     cwd->cd();
511   }
512     
513   //Now CPV branch
514   TBranch * cpvBranch = gAlice->TreeR()->Branch("EMCALPreShoRP","TObjArray",&preshoRecPoints,bufferSize,splitlevel);
515   cpvBranch->SetTitle(branchName);
516   if (filename) {
517     cpvBranch->SetFile(filename);
518     TIter next( cpvBranch->GetListOfBranches());
519     TBranch * sb;
520     while ((sb=(TBranch*)next())) {
521       sb->SetFile(filename);
522     }   
523     cwd->cd();
524   } 
525     
526   //And Finally  clusterizer branch
527   AliEMCALClusterizerv1 * cl = (AliEMCALClusterizerv1*)gime->Clusterizer(branchName) ;
528   TBranch * clusterizerBranch = gAlice->TreeR()->Branch("AliEMCALClusterizer","AliEMCALClusterizerv1",
529                                               &cl,bufferSize,splitlevel);
530   clusterizerBranch->SetTitle(branchName);
531   if (filename) {
532     clusterizerBranch->SetFile(filename);
533     TIter next( clusterizerBranch->GetListOfBranches());
534     TBranch * sb ;
535     while ((sb=(TBranch*)next())) {
536       sb->SetFile(filename);
537     }   
538     cwd->cd();
539   }
540   emcBranch        ->Fill() ;
541   cpvBranch        ->Fill() ;
542   clusterizerBranch->Fill() ;
543
544   gAlice->TreeR()->Write(0,kOverwrite) ;  
545   
546 }
547
548 //____________________________________________________________________________
549 void AliEMCALClusterizerv1::MakeClusters()
550 {
551   // Steering method to construct the clusters stored in a list of Reconstructed Points
552   // A cluster is defined as a list of neighbour digits
553   
554   TString branchName(GetName()) ; 
555   branchName.Remove(branchName.Index(Version())-1) ; 
556   
557   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
558   
559   TObjArray * towerRecPoints  = gime->TowerRecPoints(branchName) ; 
560   TObjArray * preshoRecPoints = gime->PreShowerRecPoints(branchName) ; 
561   towerRecPoints->Delete() ;
562   preshoRecPoints->Delete() ;
563   
564   TClonesArray * digits = gime->Digits(branchName) ; 
565   TClonesArray * digitsC =  dynamic_cast<TClonesArray*>(digits->Clone()) ;
566   
567   
568   // Clusterization starts  
569   
570   TIter nextdigit(digitsC) ; 
571   AliEMCALDigit * digit ; 
572   Bool_t notremoved = kTRUE ;
573   
574   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
575     AliEMCALRecPoint * clu = 0 ; 
576     
577     TArrayI clusterdigitslist(1500) ;   
578     Int_t index ;
579  
580     if (( IsInTower (digit)  && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fTowerClusteringThreshold  ) || 
581         ( IsInPreShower (digit) && Calibrate(digit->GetAmp(),digit->IsInPreShower()) > fPreShoClusteringThreshold  ) ) {
582       
583       Int_t iDigitInCluster = 0 ; 
584       
585       if  ( IsInTower(digit) ) {   
586         // start a new Tower RecPoint
587         if(fNumberOfTowerClusters >= towerRecPoints->GetSize()) 
588           towerRecPoints->Expand(2*fNumberOfTowerClusters+1) ;
589         
590         towerRecPoints->AddAt(new  AliEMCALTowerRecPoint(""), fNumberOfTowerClusters) ;
591         clu = dynamic_cast<AliEMCALTowerRecPoint *>(towerRecPoints->At(fNumberOfTowerClusters)) ; 
592         fNumberOfTowerClusters++ ; 
593         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower())) ; 
594         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;  
595         iDigitInCluster++ ; 
596         digitsC->Remove(digit) ; 
597         
598       } else { 
599         
600         // start a new Pre Shower cluster
601         if(fNumberOfPreShoClusters >= preshoRecPoints->GetSize()) 
602           preshoRecPoints->Expand(2*fNumberOfPreShoClusters+1);
603         
604         preshoRecPoints->AddAt(new AliEMCALTowerRecPoint(""), fNumberOfPreShoClusters) ;
605         
606         clu =  dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(fNumberOfPreShoClusters))  ;  
607         fNumberOfPreShoClusters++ ; 
608         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->IsInPreShower() ) );     
609         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList()  ; 
610         iDigitInCluster++ ; 
611         digitsC->Remove(digit) ; 
612         nextdigit.Reset() ;
613         
614         // Here we remove remaining Tower digits, which cannot make a cluster
615         
616         if( notremoved ) { 
617           while( ( digit = dynamic_cast<AliEMCALDigit *>(nextdigit()) ) ) {
618             if( IsInTower(digit) )
619               digitsC->Remove(digit) ;
620             else 
621               break ; 
622           }
623           notremoved = kFALSE ;
624         }
625         
626       } // else        
627       
628       nextdigit.Reset() ;
629       
630       AliEMCALDigit * digitN ; 
631       index = 0 ;
632       while (index < iDigitInCluster){ // scan over digits already in cluster 
633         digit =  (AliEMCALDigit*)digits->At(clusterdigitslist[index])  ;      
634         index++ ; 
635         while ( (digitN = (AliEMCALDigit *)nextdigit()) ) { // scan over the reduced list of digits 
636           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!!
637          switch (ineb ) {
638           case 0 :   // not a neighbour
639             break ;
640           case 1 :   // are neighbours 
641             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->IsInPreShower() ) ) ;
642             clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ; 
643             iDigitInCluster++ ; 
644             digitsC->Remove(digitN) ;
645             break ;
646           case 2 :   // too far from each other
647             goto endofloop;   
648           } // switch
649           
650         } // while digitN
651         
652       endofloop: ;
653         nextdigit.Reset() ; 
654         
655       } // loop over cluster     
656       
657     } // energy theshold  
658     
659     
660   } // while digit
661   
662   delete digitsC ;
663   
664 }
665
666 //____________________________________________________________________________
667 void AliEMCALClusterizerv1::MakeUnfolding()
668 {
669   Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
670   
671 //   // Unfolds clusters using the shape of an ElectroMagnetic shower
672 //   // Performs unfolding of all EMC/CPV clusters
673
674 //   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
675   
676 //   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
677 //   TObjArray * emcRecPoints = gime->TowerRecPoints() ; 
678 //   TObjArray * cpvRecPoints = gime->PreShoRecPoints() ; 
679 //   TClonesArray * digits = gime->Digits() ; 
680   
681 //   // Unfold first EMC clusters 
682 //   if(fNumberOfTowerClusters > 0){
683
684 //     Int_t nModulesToUnfold = geom->GetNModules() ; 
685
686 //     Int_t numberofNotUnfolded = fNumberOfTowerClusters ; 
687 //     Int_t index ;   
688 //     for(index = 0 ; index < numberofNotUnfolded ; index++){
689       
690 //       AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint *) emcRecPoints->At(index) ;
691 //       if(emcRecPoint->GetEMCALMod()> nModulesToUnfold)
692 //      break ;
693       
694 //       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
695 //       Int_t * maxAt = new Int_t[nMultipl] ;
696 //       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
697 //       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fTowerLocMaxCut,digits) ;
698       
699 //       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
700 //      UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
701 //      emcRecPoints->Remove(emcRecPoint); 
702 //      emcRecPoints->Compress() ;
703 //      index-- ;
704 //      fNumberOfTowerClusters -- ;
705 //      numberofNotUnfolded-- ;
706 //       }
707       
708 //       delete[] maxAt ; 
709 //       delete[] maxAtEnergy ; 
710 //     }
711 //   } 
712 //   // Unfolding of EMC clusters finished
713
714
715 //   // Unfold now CPV clusters
716 //   if(fNumberOfPreShoClusters > 0){
717     
718 //     Int_t nModulesToUnfold = geom->GetNModules() ;
719
720 //     Int_t numberofPreShoNotUnfolded = fNumberOfPreShoClusters ;     
721 //     Int_t index ;   
722 //     for(index = 0 ; index < numberofPreShoNotUnfolded ; index++){
723       
724 //       AliEMCALRecPoint * recPoint = (AliEMCALRecPoint *) cpvRecPoints->At(index) ;
725
726 //       if(recPoint->GetEMCALMod()> nModulesToUnfold)
727 //      break ;
728       
729 //       AliEMCALTowerRecPoint * emcRecPoint = (AliEMCALTowerRecPoint*) recPoint ; 
730       
731 //       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
732 //       Int_t * maxAt = new Int_t[nMultipl] ;
733 //       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
734 //       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fPreShoLocMaxCut,digits) ;
735       
736 //       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
737 //      UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
738 //      cpvRecPoints->Remove(emcRecPoint); 
739 //      cpvRecPoints->Compress() ;
740 //      index-- ;
741 //      numberofPreShoNotUnfolded-- ;
742 //      fNumberOfPreShoClusters-- ;
743 //       }
744       
745 //       delete[] maxAt ; 
746 //       delete[] maxAtEnergy ; 
747 //     } 
748 //   }
749 //   //Unfolding of PreSho clusters finished
750   
751 }
752
753 //____________________________________________________________________________
754 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t r)
755
756   // Shape of the shower (see EMCAL TDR)
757   // If you change this function, change also the gradient evaluation in ChiSquare()
758
759   Double_t r4    = r*r*r*r ;
760   Double_t r295  = TMath::Power(r, 2.95) ;
761   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
762   return shape ;
763 }
764
765 //____________________________________________________________________________
766 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALTowerRecPoint * iniTower, 
767                                                  Int_t nMax, 
768                                                  int * maxAt, 
769                                                  Float_t * maxAtEnergy)
770 {
771   // Performs the unfolding of a cluster with nMax overlapping showers 
772   
773   Fatal("AliEMCALClusterizerv1::UnfoldCluster", "--> Unfolding not implemented") ;
774
775  //  AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
776 //   const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
777 //   const TClonesArray * digits = gime->Digits() ; 
778 //   TObjArray * emcRecPoints = gime->TowerRecPoints() ; 
779 //   TObjArray * cpvRecPoints = gime->PreShoRecPoints() ; 
780
781 //   Int_t nPar = 3 * nMax ;
782 //   Float_t * fitparameters = new Float_t[nPar] ;
783
784 //   Bool_t rv = FindFit(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;
785 //   if( !rv ) {
786 //     // Fit failed, return and remove cluster
787 //     delete[] fitparameters ; 
788 //     return ;
789 //   }
790
791 //   // create ufolded rec points and fill them with new energy lists
792 //   // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
793 //   // and later correct this number in acordance with actual energy deposition
794
795 //   Int_t nDigits = iniTower->GetMultiplicity() ;  
796 //   Float_t * efit = new Float_t[nDigits] ;
797 //   Float_t xDigit=0.,zDigit=0.,distance=0. ;
798 //   Float_t xpar=0.,zpar=0.,epar=0.  ;
799 //   Int_t relid[4] ;
800 //   AliEMCALDigit * digit = 0 ;
801 //   Int_t * emcDigits = iniTower->GetDigitsList() ;
802
803 //   Int_t iparam ;
804 //   Int_t iDigit ;
805 //   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
806 //     digit = (AliEMCALDigit*) digits->At(emcDigits[iDigit] ) ;   
807 //     geom->AbsToRelNumbering(digit->GetId(), relid) ;
808 //     geom->RelPosInModule(relid, xDigit, zDigit) ;
809 //     efit[iDigit] = 0;
810
811 //     iparam = 0 ;    
812 //     while(iparam < nPar ){
813 //       xpar = fitparameters[iparam] ;
814 //       zpar = fitparameters[iparam+1] ;
815 //       epar = fitparameters[iparam+2] ;
816 //       iparam += 3 ;
817 //       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
818 //       distance =  TMath::Sqrt(distance) ;
819 //       efit[iDigit] += epar * ShowerShape(distance) ;
820 //     }
821 //   }
822   
823
824 //   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
825 //   // so that energy deposited in each cell is distributed betwin new clusters proportionally
826 //   // to its contribution to efit
827
828 //   Float_t * emcEnergies = iniTower->GetEnergiesList() ;
829 //   Float_t ratio ;
830
831 //   iparam = 0 ;
832 //   while(iparam < nPar ){
833 //     xpar = fitparameters[iparam] ;
834 //     zpar = fitparameters[iparam+1] ;
835 //     epar = fitparameters[iparam+2] ;
836 //     iparam += 3 ;    
837     
838 //     AliEMCALTowerRecPoint * emcRP = 0 ;  
839
840 //     if(iniTower->IsTower()){ //create new entries in fTowerRecPoints...
841       
842 //       if(fNumberOfTowerClusters >= emcRecPoints->GetSize())
843 //      emcRecPoints->Expand(2*fNumberOfTowerClusters) ;
844       
845 //       (*emcRecPoints)[fNumberOfTowerClusters] = new AliEMCALTowerRecPoint("") ;
846 //       emcRP = (AliEMCALTowerRecPoint *) emcRecPoints->At(fNumberOfTowerClusters);
847 //       fNumberOfTowerClusters++ ;
848 //     }
849 //     else{//create new entries in fPreShoRecPoints
850 //       if(fNumberOfPreShoClusters >= cpvRecPoints->GetSize())
851 //      cpvRecPoints->Expand(2*fNumberOfPreShoClusters) ;
852       
853 //       (*cpvRecPoints)[fNumberOfPreShoClusters] = new AliEMCALPreShoRecPoint("") ;
854 //       emcRP = (AliEMCALTowerRecPoint *) cpvRecPoints->At(fNumberOfPreShoClusters);
855 //       fNumberOfPreShoClusters++ ;
856 //     }
857     
858 //     Float_t eDigit ;
859 //     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
860 //       digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ; 
861 //       geom->AbsToRelNumbering(digit->GetId(), relid) ;
862 //       geom->RelPosInModule(relid, xDigit, zDigit) ;
863 //       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
864 //       distance =  TMath::Sqrt(distance) ;
865 //       ratio = epar * ShowerShape(distance) / efit[iDigit] ; 
866 //       eDigit = emcEnergies[iDigit] * ratio ;
867 //       emcRP->AddDigit( *digit, eDigit ) ;
868 //     }        
869 //   }
870  
871 //   delete[] fitparameters ; 
872 //   delete[] efit ; 
873   
874 }
875
876 //_____________________________________________________________________________
877 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
878 {
879   // Calculates the Chi square for the cluster unfolding minimization
880   // Number of parameters, Gradient, Chi squared, parameters, what to do
881
882   abort() ; 
883  //  Fatal("AliEMCALClusterizerv1::UnfoldingChiSquare","-->Unfolding not implemented") ;
884
885 //   TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
886
887 //   AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint*) toMinuit->At(0)  ;
888 //   TClonesArray * digits = (TClonesArray*)toMinuit->At(1)  ;
889
890
891   
892 //   //  AliEMCALTowerRecPoint * emcRP = (AliEMCALTowerRecPoint *) gMinuit->GetObjectFit() ; // TowerRecPoint to fit
893
894 //   Int_t * emcDigits     = emcRP->GetDigitsList() ;
895
896 //   Int_t nOdigits = emcRP->GetDigitsMultiplicity() ; 
897
898 //   Float_t * emcEnergies = emcRP->GetEnergiesList() ;
899
900 //   const AliEMCALGeometry * geom = AliEMCALGetter::GetInstance()->EMCALGeometry() ; 
901 //   fret = 0. ;     
902 //   Int_t iparam ;
903
904 //   if(iflag == 2)
905 //     for(iparam = 0 ; iparam < nPar ; iparam++)    
906 //       Grad[iparam] = 0 ; // Will evaluate gradient
907   
908 //   Double_t efit ;    
909
910 //   AliEMCALDigit * digit ;
911 //   Int_t iDigit ;
912
913 //   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
914
915 //     digit = (AliEMCALDigit*) digits->At( emcDigits[iDigit] ) ; 
916
917 //     Int_t relid[4] ;
918 //     Float_t xDigit ;
919 //     Float_t zDigit ;
920
921 //     geom->AbsToRelNumbering(digit->GetId(), relid) ;
922
923 //     geom->RelPosInModule(relid, xDigit, zDigit) ;
924
925 //      if(iflag == 2){  // calculate gradient
926 //        Int_t iParam = 0 ;
927 //        efit = 0 ;
928 //        while(iParam < nPar ){
929 //       Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
930 //       iParam++ ; 
931 //       distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ; 
932 //       distance = TMath::Sqrt( distance ) ; 
933 //       iParam++ ;      
934 //       efit += x[iParam] * ShowerShape(distance) ;
935 //       iParam++ ;
936 //        }
937 //        Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E) 
938 //        iParam = 0 ;
939 //        while(iParam < nPar ){
940 //       Double_t xpar = x[iParam] ;
941 //       Double_t zpar = x[iParam+1] ;
942 //       Double_t epar = x[iParam+2] ;
943 //       Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
944 //       Double_t shape = sum * ShowerShape(dr) ;
945 //       Double_t r4 = dr*dr*dr*dr ;
946 //       Double_t r295 = TMath::Power(dr,2.95) ;
947 //       Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
948 //                                       0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
949          
950 //       Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x    
951 //       iParam++ ; 
952 //       Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z         
953 //       iParam++ ; 
954 //       Grad[iParam] += shape ;                                  // Derivative over energy             
955 //       iParam++ ; 
956 //        }
957 //      }
958 //      efit = 0;
959 //      iparam = 0 ;
960
961 //      while(iparam < nPar ){
962 //        Double_t xpar = x[iparam] ;
963 //        Double_t zpar = x[iparam+1] ;
964 //        Double_t epar = x[iparam+2] ;
965 //        iparam += 3 ;
966 //        Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
967 //        distance =  TMath::Sqrt(distance) ;
968 //        efit += epar * ShowerShape(distance) ;
969 //      }
970
971 //      fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ; 
972 //      // Here we assume, that sigma = sqrt(E)
973 //   }
974
975 }
976
977 //____________________________________________________________________________
978 void AliEMCALClusterizerv1::Print(Option_t * option)const
979 {
980   // Print clusterizer parameters
981
982   if( strcmp(GetName(), "") !=0 ){
983     
984     // Print parameters
985  
986     TString taskName(GetName()) ; 
987     taskName.ReplaceAll(Version(), "") ;
988
989     cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl 
990          << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl 
991          << "                           Branch: " << fDigitsBranchTitle.Data() << endl 
992          << endl 
993          << "                       EMC Clustering threshold = " << fTowerClusteringThreshold << endl
994          << "                       EMC Local Maximum cut    = " << fTowerLocMaxCut << endl
995          << "                       EMC Logarothmic weight   = " << fW0 << endl
996          << endl
997          << "                       CPV Clustering threshold = " << fPreShoClusteringThreshold << endl
998          << "                       CPV Local Maximum cut    = " << fPreShoLocMaxCut << endl
999        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
1000          << endl ;
1001     if(fToUnfold)
1002       cout << " Unfolding on " << endl ;
1003     else
1004       cout << " Unfolding off " << endl ;
1005     
1006     cout << "------------------------------------------------------------------" <<endl ;
1007   }
1008   else
1009     cout << " AliEMCALClusterizerv1 not initialized " << endl ;
1010 }
1011 //____________________________________________________________________________
1012 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
1013 {
1014   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
1015
1016   TObjArray * towerRecPoints = AliEMCALGetter::GetInstance()->TowerRecPoints() ; 
1017   TObjArray * preshoRecPoints = AliEMCALGetter::GetInstance()->PreShowerRecPoints() ; 
1018
1019   cout << "AliEMCALClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1020   cout << "       Found "<< towerRecPoints->GetEntriesFast() << " TOWER Rec Points and " 
1021            << preshoRecPoints->GetEntriesFast() << " PRE SHOWER RecPoints" << endl ;
1022
1023   fRecPointsInRun +=  towerRecPoints->GetEntriesFast() ; 
1024   fRecPointsInRun +=  preshoRecPoints->GetEntriesFast() ; 
1025
1026   if(strstr(option,"all")) {
1027
1028     cout << "Tower clusters " << endl ;
1029     cout << " Index  Ene(MeV)   Multi  Module     phi     r  theta    Lambda 1   Lambda 2  # of prim  Primaries list "      <<  endl;      
1030     
1031     Int_t index ;
1032     for (index = 0 ; index < towerRecPoints->GetEntries() ; index++) {
1033       AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint * >(towerRecPoints->At(index)) ; 
1034       TVector3  globalpos;  
1035       rp->GetGlobalPosition(globalpos);
1036       Float_t lambda[2]; 
1037       rp->GetElipsAxis(lambda);
1038       Int_t * primaries; 
1039       Int_t nprimaries;
1040       primaries = rp->GetPrimaries(nprimaries);
1041
1042       cout << setw(4) << rp->GetIndexInList() << "   " 
1043            << setw(7) << setprecision(3) << rp->GetEnergy() << "      "
1044            << setw(3) << rp->GetMultiplicity() << "      " 
1045            << setw(1) << rp->GetEMCALArm() << "     " 
1046            << setw(5) << setprecision(4) << globalpos.X() << "  " 
1047            << setw(5) << setprecision(4) << globalpos.Y() << "  " 
1048            << setw(5) << setprecision(4) << globalpos.Z() << "     "
1049            << setw(4) << setprecision(2) << lambda[0]  << "  "
1050            << setw(4) << setprecision(2) << lambda[1]  << "  "
1051            << setw(2) << nprimaries << "  " ;
1052      
1053       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1054         cout << setw(4) <<   primaries[iprimary] << "  "  ;
1055       cout << endl ;     
1056     }
1057
1058     //Now plot Pre shower recPoints
1059
1060     cout << "-----------------------------------------------------------------------"<<endl ;
1061
1062     cout << "PreShower clusters " << endl ;
1063     cout << " Index  Ene(MeV)   Multi  Module     phi     r  theta    Lambda 1   Lambda 2  # of prim  Primaries list "      <<  endl;      
1064     
1065     for (index = 0 ; index < preshoRecPoints->GetEntries() ; index++) {
1066       AliEMCALTowerRecPoint * rp = dynamic_cast<AliEMCALTowerRecPoint *>(preshoRecPoints->At(index)) ; 
1067       TVector3  globalpos;  
1068       rp->GetGlobalPosition(globalpos);
1069       Float_t lambda[2]; 
1070       rp->GetElipsAxis(lambda);
1071       Int_t * primaries; 
1072       Int_t nprimaries;
1073       primaries = rp->GetPrimaries(nprimaries);
1074
1075       cout << setw(4) << rp->GetIndexInList() << "   " 
1076            << setw(7) << setprecision(3) << rp->GetEnergy() << "      "
1077            << setw(3) << rp->GetMultiplicity() << "      " 
1078            << setw(1) << rp->GetEMCALArm() << "     " 
1079            << setw(5) << setprecision(4) << globalpos.X() << "  " 
1080            << setw(5) << setprecision(4) << globalpos.Y() << "  " 
1081            << setw(5) << setprecision(4) << globalpos.Z() << "     "
1082            << setw(4) << setprecision(2) << lambda[0]  << "  "
1083            << setw(4) << setprecision(2) << lambda[1]  << "  "
1084            << setw(2) << nprimaries << "  " ;
1085      
1086       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1087         cout << setw(4) <<   primaries[iprimary] << "  "  ;
1088       cout << endl ;     
1089     }
1090
1091     cout << "-----------------------------------------------------------------------"<<endl ;
1092   }
1093 }
1094