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