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