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