]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSClusterizerv1.cxx
can now only write mixed digits into one of the file that contains summable digits...
[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 file name 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 input file
46 // root [3] cl->SetRecPointsBranch("recp2")  
47 //               //sets another aouput file
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)
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     file = new TFile(fHeaderFileName.Data(),"update") ;
143     gAlice = (AliRun *) file->Get("gAlice") ;
144   }
145   
146   AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
147   fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
148   
149   fDigits = new TClonesArray("AliPHOSDigit",10) ;
150   fDigitizer = new AliPHOSDigitizer() ;
151   fEmcRecPoints = new TObjArray(200) ;
152   fCpvRecPoints = new TObjArray(200) ;
153   
154   if(!gMinuit) gMinuit = new TMinuit(100) ;
155   
156   // add Task to //root/Tasks folder
157   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
158   roottasks->Add(this) ; 
159
160
161   fIsInitialized = kTRUE ;
162
163 }
164 //____________________________________________________________________________
165 void AliPHOSClusterizerv1::Exec(Option_t * option){
166   // Steerign function
167
168   if(!fIsInitialized) Init() ;
169
170   if(strstr(option,"tim"))
171     gBenchmark->Start("PHOSClusterizer");  
172   
173   Int_t nEvents = (Int_t) gAlice->TreeE()->GetEntries() ;
174   
175   for(fEvent = 0;fEvent< nEvents; fEvent++){
176     if(!ReadDigits())  //reads digits for event fEvent
177       return;
178     MakeClusters() ;
179     
180     if(fToUnfold) MakeUnfolding() ;
181     WriteRecPoints() ;
182     if(strstr(option,"deb"))
183       PrintRecPoints(option) ;
184   }
185   
186   if(strstr(option,"tim")){
187     gBenchmark->Stop("PHOSClusterizer");
188     cout << "AliPHOSClusterizer:" << endl ;
189     cout << "  took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing " 
190          <<  gBenchmark->GetCpuTime("PHOSClusterizer")/nEvents << " seconds per event " << endl ;
191     cout << endl ;
192   }
193   
194   
195 }
196
197 //____________________________________________________________________________
198 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
199                                     Int_t nPar, Float_t * fitparameters)
200
201   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
202
203   gMinuit->mncler();                     // Reset Minuit's list of paramters
204   gMinuit->SetPrintLevel(-1) ;           // No Printout
205   gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;  
206                                          // To set the address of the minimization function 
207
208   TList * toMinuit = new TList();
209   toMinuit->AddAt(emcRP,0) ;
210   toMinuit->AddAt(fDigits,1) ;
211   
212   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
213
214   // filling initial values for fit parameters
215   AliPHOSDigit * digit ;
216
217   Int_t ierflg  = 0; 
218   Int_t index   = 0 ;
219   Int_t nDigits = (Int_t) nPar / 3 ;
220
221   Int_t iDigit ;
222
223
224   for(iDigit = 0; iDigit < nDigits; iDigit++){
225     digit = (AliPHOSDigit *) maxAt[iDigit]; 
226
227     Int_t relid[4] ;
228     Float_t x ;
229     Float_t z ;
230     fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
231     fGeom->RelPosInModule(relid, x, z) ;
232
233     Float_t energy = maxAtEnergy[iDigit] ;
234
235     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
236     index++ ;   
237     if(ierflg != 0){ 
238       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
239       return kFALSE;
240     }
241     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
242     index++ ;   
243     if(ierflg != 0){
244       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
245       return kFALSE;
246     }
247     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
248     index++ ;   
249     if(ierflg != 0){
250       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
251       return kFALSE;
252     }
253   }
254
255   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
256                       //  depends on it. 
257   Double_t p1 = 1.0 ;
258   Double_t p2 = 0.0 ;
259
260   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
261   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
262   gMinuit->SetMaxIterations(5);
263   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
264
265   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
266
267   if(ierflg == 4){  // Minimum not found   
268     cout << "PHOS Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
269     return kFALSE ;
270   }            
271   for(index = 0; index < nPar; index++){
272     Double_t err ;
273     Double_t val ;
274     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
275     fitparameters[index] = val ;
276    }
277
278   delete toMinuit ;
279   return kTRUE;
280
281 }
282
283 //____________________________________________________________________________
284 void AliPHOSClusterizerv1::Init(){
285
286   if(!fIsInitialized){
287     if(fHeaderFileName.IsNull())
288       fHeaderFileName = "galice.root" ;
289     
290
291     TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
292
293     if(file == 0){
294       file = new TFile(fHeaderFileName.Data(),"update") ;
295       gAlice = (AliRun *) file->Get("gAlice") ;
296     }
297
298     AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
299     fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
300
301     fDigits = new TClonesArray("AliPHOSDigit",10) ;
302     fDigitizer = new AliPHOSDigitizer() ;
303     fEmcRecPoints = new TObjArray(200) ;
304     fCpvRecPoints = new TObjArray(200) ;
305      
306     if(!gMinuit) gMinuit = new TMinuit(100) ;
307
308     // add Task to //root/Tasks folder
309     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
310     roottasks->Add(this) ; 
311
312     fIsInitialized = kTRUE ;
313   }
314 }
315 //____________________________________________________________________________
316 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
317 {
318   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
319   //                                       = 1 are neighbour
320   //                                       = 2 are not neighbour but do not continue searching
321   // neighbours are defined as digits having at least common vertex 
322   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
323   //                                      which is compared to a digit (d2)  not yet in a cluster  
324
325   Int_t rv = 0 ; 
326
327   Int_t relid1[4] ; 
328   fGeom->AbsToRelNumbering(d1->GetId(), relid1) ; 
329
330   Int_t relid2[4] ; 
331   fGeom->AbsToRelNumbering(d2->GetId(), relid2) ; 
332  
333   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module and the same PPSD Module 
334     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
335     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
336     
337     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
338       rv = 1 ; 
339     }
340     else {
341       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
342         rv = 2; //  Difference in row numbers is too large to look further 
343     }
344
345   } 
346   else {
347     
348     if( (relid1[0] < relid2[0]) || (relid1[1] < relid2[1]) )  
349       rv=2 ;
350
351   }
352
353   //Do NOT clusterize upper PPSD  
354   if( IsInPpsd(d1) && IsInPpsd(d2) &&
355      relid1[1] > 0                 &&
356      relid1[1] < fGeom->GetNumberOfPadsPhi()*fGeom->GetNumberOfPadsPhi() ) rv = 2 ;
357
358   return rv ; 
359 }
360
361
362 //____________________________________________________________________________
363 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
364 {
365   // Tells if (true) or not (false) the digit is in a PHOS-EMC module
366  
367   Bool_t rv = kFALSE ; 
368
369   Int_t relid[4] ; 
370   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
371
372   if ( relid[1] == 0  ) rv = kTRUE; 
373
374   return rv ; 
375 }
376
377 //____________________________________________________________________________
378 Bool_t AliPHOSClusterizerv1::IsInPpsd(AliPHOSDigit * digit) const
379 {
380   // Tells if (true) or not (false) the digit is in a PHOS-PPSD module
381  
382   Bool_t rv = kFALSE ; 
383
384   Int_t relid[4] ; 
385   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
386
387   if ( relid[1] > 0 && relid[0] > fGeom->GetNCPVModules() ) rv = kTRUE; 
388
389   return rv ; 
390 }
391
392 //____________________________________________________________________________
393 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
394 {
395   // Tells if (true) or not (false) the digit is in a PHOS-CPV module
396  
397   Bool_t rv = kFALSE ; 
398
399   Int_t relid[4] ; 
400   fGeom->AbsToRelNumbering(digit->GetId(), relid) ; 
401
402   if ( relid[1] > 0 && relid[0] <= fGeom->GetNCPVModules() ) rv = kTRUE; 
403
404   return rv ; 
405 }
406 //____________________________________________________________________________
407 Bool_t AliPHOSClusterizerv1::ReadDigits(){
408
409   fNumberOfEmcClusters  = 0 ;
410   fNumberOfCpvClusters  = 0 ;
411
412   // Get Digits Tree header from file
413   char treeName[20]; 
414   sprintf(treeName,"TreeD%d",fEvent);
415   gAlice->GetEvent(fEvent) ;
416   gAlice->SetEvent(fEvent) ;
417
418   TTree * treeD = gAlice->TreeD()  ; // (TTree*)file->Get(treeName);
419
420   if(treeD==0){
421     cout << "Error in AliPHOSClusterizerv1 : no "<<treeName << endl  ;
422     cout << "    Do nothing " << endl ;
423     return kFALSE ;
424   }
425
426   TBranch * digitsBranch = 0;
427   TBranch * digitizerBranch = 0;
428
429   TObjArray * branches = treeD->GetListOfBranches() ;
430   Int_t ibranch;
431   Bool_t phosNotFound = kTRUE ;
432   Bool_t digitizerNotFound = kTRUE ;
433   
434   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
435
436     if(phosNotFound){
437       digitsBranch=(TBranch *) branches->At(ibranch) ;
438       if( fDigitsBranchTitle.CompareTo(digitsBranch->GetTitle())==0 )
439         if( strcmp(digitsBranch->GetName(),"PHOS") == 0)
440           phosNotFound = kFALSE ;
441     }
442     
443     if(digitizerNotFound){
444       digitizerBranch = (TBranch *) branches->At(ibranch) ;
445       if( fDigitsBranchTitle.CompareTo(digitizerBranch->GetTitle()) == 0)
446         if( strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) 
447           digitizerNotFound = kFALSE ;
448     }
449     
450   }
451   
452   if(digitizerNotFound || phosNotFound){
453     cout << "ERROR in AliPHOSClusterizerv1: " << endl ;
454     cout << "      Can't find Branch with digits or Digitizer "<< endl ; ;
455     cout << "      Do nothing" <<endl  ;
456     return kFALSE ;
457   }
458   
459   digitsBranch->SetAddress(&fDigits) ;
460   digitizerBranch->SetAddress(&fDigitizer) ;
461   
462   treeD->GetEvent(0) ;
463   
464   fPedestal = fDigitizer->GetPedestal() ;
465   fSlope    = fDigitizer->GetSlope() ;
466   return kTRUE ;
467 }
468
469 //____________________________________________________________________________
470 void AliPHOSClusterizerv1::WriteRecPoints(){
471   
472   Int_t index ;
473   //Evaluate poisition, dispersion and other RecPoint properties...
474   for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
475     ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->EvalAll(fW0,fDigits) ;
476
477   fEmcRecPoints->Sort() ;
478
479   for(index = 0; index < fEmcRecPoints->GetEntries(); index++)
480     ((AliPHOSEmcRecPoint *)fEmcRecPoints->At(index))->SetIndexInList(index) ;
481
482   //Now the same for CPV
483   for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
484     ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->EvalAll(fW0CPV,fDigits)  ;
485
486   fCpvRecPoints->Sort() ;
487
488   for(index = 0; index < fCpvRecPoints->GetEntries(); index++)
489     ((AliPHOSRecPoint *)fCpvRecPoints->At(index))->SetIndexInList(index) ;
490
491   if(gAlice->TreeR()==0)
492     gAlice->MakeTree("R") ;
493   
494
495   //Check, if branches already exist
496   TBranch * emcBranch = 0;
497   TBranch * cpvBranch = 0;
498   TBranch * clusterizerBranch = 0;
499   
500   TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
501   Int_t ibranch;
502   Bool_t emcNotFound = kTRUE ;
503   Bool_t cpvNotFound = kTRUE ;  
504   Bool_t clusterizerNotFound = kTRUE ;
505   
506   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
507     
508     if(emcNotFound){
509       emcBranch=(TBranch *) branches->At(ibranch) ;
510       if( fRecPointsBranchTitle.CompareTo(emcBranch->GetTitle())==0 ){
511         if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) {
512           emcNotFound = kFALSE ;
513         }
514       }
515     }
516     if(cpvNotFound){
517       cpvBranch=(TBranch *) branches->At(ibranch) ;
518       if( fRecPointsBranchTitle.CompareTo(cpvBranch->GetTitle())==0 ){
519         if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) {
520           cpvNotFound = kFALSE ;
521         }
522       }
523     }
524
525     if(clusterizerNotFound){
526       clusterizerBranch = (TBranch *) branches->At(ibranch) ;
527       if( fRecPointsBranchTitle.CompareTo(clusterizerBranch->GetTitle()) == 0){
528         if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0) {
529           clusterizerNotFound = kFALSE ;
530         }
531       }
532     }
533     
534   }
535   
536   if(!(clusterizerNotFound && emcNotFound && cpvNotFound)){
537     cout << "AliPHOSClusterizer error" << endl;
538     cout << "       Branches PHOSEmcRP, PHOSCpvRP and AliPHOSClusterizer " << endl ;
539     cout << "       with title '" << fRecPointsBranchTitle.Data() <<"' already exist" << endl ;
540     cout << "       can not overwrite " << endl ;
541     return ;
542   }
543
544   //Make branches in TreeR for RecPoints and Clusterizer
545   char * filename = 0;
546   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
547     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
548     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
549   }
550   
551   //Make new branches
552   TDirectory *cwd = gDirectory;
553   
554   //First EMC
555   Int_t bufferSize = 32000 ;    
556   Int_t splitlevel = 0 ;
557   emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&fEmcRecPoints,bufferSize,splitlevel);
558   emcBranch->SetTitle(fRecPointsBranchTitle.Data());
559   if (filename) {
560     emcBranch->SetFile(filename);
561     TIter next( emcBranch->GetListOfBranches());
562     while ((emcBranch=(TBranch*)next())) {
563       emcBranch->SetFile(filename);
564     }   
565     cwd->cd();
566   }
567     
568   //Now CPV branch
569   cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&fCpvRecPoints,bufferSize,splitlevel);
570   cpvBranch->SetTitle(fRecPointsBranchTitle.Data());
571   if (filename) {
572     cpvBranch->SetFile(filename);
573     TIter next( cpvBranch->GetListOfBranches());
574     while ((cpvBranch=(TBranch*)next())) {
575       cpvBranch->SetFile(filename);
576     }   
577     cwd->cd();
578   } 
579     
580   //And Finally  clusterizer branch
581   AliPHOSClusterizerv1 * cl = this ;
582   clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
583                                               &cl,bufferSize,splitlevel);
584   clusterizerBranch->SetTitle(fRecPointsBranchTitle.Data());
585   if (filename) {
586     clusterizerBranch->SetFile(filename);
587     TIter next( clusterizerBranch->GetListOfBranches());
588     while ((clusterizerBranch=(TBranch*)next())) {
589       clusterizerBranch->SetFile(filename);
590     }   
591     cwd->cd();
592   }
593   
594   gAlice->TreeR()->Fill() ;
595   
596   gAlice->TreeR()->Write(0,kOverwrite) ;  
597   
598 }
599
600 //____________________________________________________________________________
601 void AliPHOSClusterizerv1::MakeClusters()
602 {
603   // Steering method to construct the clusters stored in a list of Reconstructed Points
604   // A cluster is defined as a list of neighbour digits
605   fEmcRecPoints->Clear() ;
606   fCpvRecPoints->Clear() ;
607   
608
609   // Clusterization starts  
610   TClonesArray * digits =  (TClonesArray*)fDigits->Clone() ;
611
612   TIter nextdigit(digits) ; 
613   AliPHOSDigit * digit ; 
614   Bool_t notremoved = kTRUE ;
615
616   while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digits
617     AliPHOSRecPoint * clu = 0 ; 
618
619     TArrayI clusterdigitslist(1000) ;   
620     Int_t index ;
621
622     if (( IsInEmc (digit) && Calibrate(digit->GetAmp()) > fEmcClusteringThreshold  ) || 
623         ( IsInPpsd(digit) && Calibrate(digit->GetAmp()) > fPpsdClusteringThreshold ) ||
624         ( IsInCpv (digit) && Calibrate(digit->GetAmp()) > fCpvClusteringThreshold  ) ) {
625       
626       Int_t iDigitInCluster = 0 ; 
627
628       if  ( IsInEmc(digit) ) {   
629         // start a new EMC RecPoint
630         if(fNumberOfEmcClusters >= fEmcRecPoints->GetSize()) fEmcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
631         fEmcRecPoints->AddAt(new  AliPHOSEmcRecPoint(), fNumberOfEmcClusters) ;
632         clu = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(fNumberOfEmcClusters) ; 
633         fNumberOfEmcClusters++ ; 
634         clu->AddDigit(*digit, Calibrate(digit->GetAmp())) ; 
635         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;  
636         iDigitInCluster++ ; 
637         digits->Remove(digit) ; 
638
639       } else { 
640         
641         // start a new PPSD/CPV cluster
642         if(fNumberOfCpvClusters >= fCpvRecPoints->GetSize()) fCpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
643         if(IsInPpsd(digit)) 
644           fCpvRecPoints->AddAt(new AliPHOSPpsdRecPoint(),fNumberOfCpvClusters) ;
645         else
646           fCpvRecPoints->AddAt(new AliPHOSCpvRecPoint(), fNumberOfCpvClusters) ;
647         clu =  (AliPHOSPpsdRecPoint *) fCpvRecPoints->At(fNumberOfCpvClusters)  ;  
648         fNumberOfCpvClusters++ ; 
649
650         clu->AddDigit(*digit, Calibrate(digit->GetAmp()) ) ;    
651         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList()  ; 
652         iDigitInCluster++ ; 
653         digits->Remove(digit) ; 
654         nextdigit.Reset() ;
655         
656         // Here we remove resting EMC digits, which cannot make cluster
657         
658         if( notremoved ) { 
659           while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
660             if( IsInEmc(digit) ) 
661               digits->Remove(digit) ;
662             else 
663               break ;
664           }
665           notremoved = kFALSE ;
666         }
667         
668       } // else        
669       
670       nextdigit.Reset() ;
671       
672       AliPHOSDigit * digitN ; 
673       index = 0 ;
674       while (index < iDigitInCluster){ // scan over digits already in cluster 
675         digit =  (AliPHOSDigit*)fDigits->At(clusterdigitslist[index])  ;      
676         index++ ; 
677         while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits 
678           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!!
679           switch (ineb ) {
680           case 0 :   // not a neighbour
681             break ;
682           case 1 :   // are neighbours 
683             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp() ) ) ;
684             clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ; 
685             iDigitInCluster++ ; 
686             digits->Remove(digitN) ;
687             break ;
688           case 2 :   // too far from each other
689             goto endofloop;   
690           } // switch
691           
692         } // while digitN
693         
694       endofloop: ;
695         nextdigit.Reset() ; 
696         
697       } // loop over cluster     
698
699     } // energy theshold  
700
701     
702   } // while digit
703
704   delete digits ;
705
706 }
707
708 //____________________________________________________________________________
709 void AliPHOSClusterizerv1::MakeUnfolding(){
710   //Unfolds clusters using the shape of ElectroMagnetic shower
711   // Performs unfolding of all EMC/CPV but NOT ppsd clusters
712
713   //Unfold first EMC clusters 
714   if(fNumberOfEmcClusters > 0){
715
716     Int_t nModulesToUnfold = fGeom->GetNModules() ; 
717
718     Int_t numberofNotUnfolded = fNumberOfEmcClusters ; 
719     Int_t index ;   
720     for(index = 0 ; index < numberofNotUnfolded ; index++){
721       
722       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(index) ;
723       if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
724         break ;
725       
726       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
727       Int_t * maxAt = new Int_t[nMultipl] ;
728       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
729       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,fDigits) ;
730       
731       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
732         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
733         fEmcRecPoints->Remove(emcRecPoint); 
734         fEmcRecPoints->Compress() ;
735         index-- ;
736         fNumberOfEmcClusters -- ;
737         numberofNotUnfolded-- ;
738       }
739       
740       delete[] maxAt ; 
741       delete[] maxAtEnergy ; 
742     }
743   } 
744   //Unfolding of EMC clusters finished
745
746
747   //Unfold now CPV clusters
748   if(fNumberOfCpvClusters > 0){
749     
750     Int_t nModulesToUnfold = fGeom->GetNCPVModules() ;
751
752     Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;     
753     Int_t index ;   
754     for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
755       
756       AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) fCpvRecPoints->At(index) ;
757
758       if(recPoint->GetPHOSMod()> nModulesToUnfold)
759         break ;
760       
761       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ; 
762       
763       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
764       Int_t * maxAt = new Int_t[nMultipl] ;
765       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
766       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,fDigits) ;
767       
768       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
769         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
770         fCpvRecPoints->Remove(emcRecPoint); 
771         fCpvRecPoints->Compress() ;
772         index-- ;
773         numberofCpvNotUnfolded-- ;
774         fNumberOfCpvClusters-- ;
775       }
776       
777       delete[] maxAt ; 
778       delete[] maxAtEnergy ; 
779     } 
780   }
781   //Unfolding of Cpv clusters finished
782   
783 }
784
785 //____________________________________________________________________________
786 void AliPHOSClusterizerv1::SetDigitsBranch(const char * title){
787   
788     fDigitsBranchTitle = title  ; 
789
790 }
791 //____________________________________________________________________________
792 void AliPHOSClusterizerv1::SetRecPointsBranch(const char * title){
793   
794     fRecPointsBranchTitle = title;
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 gradien 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 th 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   if(fIsInitialized){
1017     
1018     // Print parameters
1019     
1020     cout << "---------------"<< GetName() << " " << GetTitle()<< "-----------" << endl 
1021          << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl 
1022          << "                           Branch: " << fDigitsBranchTitle.Data() << endl 
1023          << endl 
1024          << "                       EMC Clustering threshold = " << fEmcClusteringThreshold << endl
1025          << "                       EMC Local Maximum cut    = " << fEmcLocMaxCut << endl
1026          << "                       EMC Logarothmic weight   = " << fW0 << endl
1027          << endl
1028          << "                       CPV Clustering threshold = " << fCpvClusteringThreshold << endl
1029          << "                       CPV Local Maximum cut    = " << fCpvLocMaxCut << endl
1030        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
1031          << endl
1032          << "                      PPSD  Clustering threshold = " << fPpsdClusteringThreshold << endl;
1033     if(fToUnfold)
1034       cout << " Unfolding on " << endl ;
1035     else
1036       cout << " Unfolding off " << endl ;
1037     
1038     cout << "------------------------------------------------------------------" <<endl ;
1039   }
1040   else
1041     cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1042 }
1043 //____________________________________________________________________________
1044 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option){
1045   //Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1046
1047   cout << "AliPHOSClusterizerv1: " << endl ;
1048   cout << "       Found "<< fEmcRecPoints->GetEntriesFast() << " EMC Rec Points and " 
1049            << fCpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1050
1051   if(strstr(option,"all")) {
1052     cout << "EMC clusters " << endl ;
1053     cout << "  Index    " 
1054          << "  Ene(MeV) " 
1055          << "  Multi    "
1056          << "  Module   "  
1057          << "    X      "
1058          << "    Y      "
1059          << "    Z      "
1060          << " Lambda 1  "
1061          << " Lambda 2  "
1062          << " MaxEnergy "
1063          << " # of prim "
1064          << " Primaries list "      <<  endl;      
1065     
1066     Int_t index ;
1067     for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1068       AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )fEmcRecPoints->At(index) ; 
1069
1070       cout << setw(6) << rp->GetIndexInList() << "     ";
1071       cout << setw(6) << rp->GetEnergy()      << "     ";
1072       cout << setw(6) << rp->GetMultiplicity()<< "     ";
1073       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1074
1075       TVector3  locpos;  
1076       rp->GetLocalPosition(locpos);
1077       cout << setw(8) <<  locpos.X()          << "     ";
1078       cout << setw(8) <<  locpos.Y()          << "     ";
1079       cout << setw(8) <<  locpos.Z()          << "     ";
1080
1081       Float_t lambda[2]; 
1082       rp->GetElipsAxis(lambda);
1083       cout << setw(10)<<  lambda[0]           << "     ";
1084       cout << setw(10)<<  lambda[1]           << "     ";
1085       
1086       
1087       Int_t * primaries; 
1088       Int_t nprimaries;
1089       primaries = rp->GetPrimaries(nprimaries);
1090       cout << setw(8) <<    primaries         << "     ";
1091
1092       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1093         cout << setw(4)  <<  primaries[iprimary] << " ";
1094       cout << endl;      
1095     }
1096     cout << endl ;
1097
1098     //Now plot PCV/PPSD recPoints
1099     cout << "EMC clusters " << endl ;
1100     cout << "  Index    " 
1101          << "  Multi    "
1102          << "  Module   "  
1103          << "  Layer    "  
1104          << "    X      "
1105          << "    Y      "
1106          << "    Z      "
1107          << " # of prim "
1108          << " Primaries list "      <<  endl;      
1109     
1110     for (index = 0 ; index < fEmcRecPoints->GetEntries() ; index++) {
1111       AliPHOSRecPoint * rp = (AliPHOSRecPoint * )fCpvRecPoints->At(index) ; 
1112       cout << setw(6) << rp->GetIndexInList() << "     ";
1113       cout << setw(6) << rp->GetPHOSMod()     << "     ";
1114
1115       if( (strcmp(rp->ClassName() , "AliPHOSPpsdRecPoint" )) == 0){
1116         AliPHOSPpsdRecPoint * ppsd = (AliPHOSPpsdRecPoint*) rp ;
1117         if(ppsd->GetUp())
1118           cout <<"        CPV     ";
1119         else
1120           cout <<"       PPSD     ";
1121       }
1122       else
1123         cout <<"        CPV     ";
1124       
1125       TVector3  locpos;  
1126       rp->GetLocalPosition(locpos);
1127       cout << setw(8) <<  locpos.X()          << "     ";
1128       cout << setw(8) <<  locpos.Y()          << "     ";
1129       cout << setw(8) <<  locpos.Z()          << "     ";
1130       
1131       Int_t * primaries; 
1132       Int_t nprimaries;
1133       primaries = rp->GetPrimaries(nprimaries);
1134       cout << setw(8) <<    primaries         << "     ";
1135
1136       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1137         cout << setw(4)  <<  primaries[iprimary] << " ";
1138       cout << endl;      
1139     }
1140
1141
1142     cout << "-------------------------------------------------"<<endl ;
1143   }
1144 }
1145