]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALUnfolding.cxx
bug fix - Adam
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALUnfolding.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 \r
16 //_________________________________________________________________________\r
17 //  Base class for the cluster unfolding algorithm \r
18 //*-- Author: Adam Matyja (SUBATECH)\r
19 //  Based on unfolding in clusterizerv1 done by Cynthia Hadjidakis\r
20 //-- Unfolding for eta~0: Cynthia Hadjidakis - still in AliEMCALCLusterizerv1\r
21 //-- Unfolding extension for whole EMCAL: Adam Matyja (SUBATECH & INP PAN)\r
22 //\r
23 //  unfolds the clusters having several local maxima. \r
24 //////////////////////////////////////////////////////////////////////////////\r
25 \r
26 // --- ROOT system ---\r
27 #include "TClonesArray.h"\r
28 #include <TMath.h> \r
29 #include <TMinuit.h>\r
30 \r
31 // --- Standard library ---\r
32 #include <cassert>\r
33 \r
34 // --- AliRoot header files ---\r
35 #include "AliEMCALUnfolding.h"\r
36 #include "AliEMCALGeometry.h"\r
37 #include "AliRunLoader.h"\r
38 #include "AliRun.h"\r
39 #include "AliEMCAL.h"\r
40 #include "AliEMCALRecParam.h"\r
41 #include "AliEMCALRecPoint.h"\r
42 #include "AliEMCALDigit.h"\r
43 #include "AliEMCALReconstructor.h"\r
44 \r
45 #include "AliLog.h"\r
46 #include "AliCDBManager.h"\r
47 class AliCDBStorage;\r
48 #include "AliCDBEntry.h"\r
49 \r
50 Double_t AliEMCALUnfolding::fgSSPars[8]={0.9262,3.365,1.548,0.1625,-0.4195,0.,0.,2.332};\r
51 Double_t AliEMCALUnfolding::fgPar5[3]={12.31,-0.007381,-0.06936};\r
52 Double_t AliEMCALUnfolding::fgPar6[3]={0.05452,0.0001228,0.001361};\r
53 \r
54 ClassImp(AliEMCALUnfolding)\r
55   \r
56 //____________________________________________________________________________\r
57 AliEMCALUnfolding::AliEMCALUnfolding():\r
58   fNumberOfECAClusters(0),\r
59   fECALocMaxCut(0),\r
60   fThreshold(0.01),//10 MeV\r
61   fGeom(NULL),\r
62   fRecPoints(NULL),\r
63   fDigitsArr(NULL)\r
64 {\r
65   // ctor with the indication of the file where header Tree and digits Tree are stored\r
66  \r
67   Init() ;\r
68 }\r
69 \r
70 //____________________________________________________________________________\r
71 AliEMCALUnfolding::AliEMCALUnfolding(AliEMCALGeometry* geometry):\r
72   fNumberOfECAClusters(0),\r
73   fECALocMaxCut(0),\r
74   fThreshold(0.01),//10 MeV\r
75   fGeom(geometry),\r
76   fRecPoints(NULL),\r
77   fDigitsArr(NULL)\r
78 {\r
79   // ctor with the indication of the file where header Tree and digits Tree are stored\r
80   // use this contructor to avoid usage of Init() which uses runloader\r
81   // change needed by HLT - MP\r
82   if (!fGeom)\r
83   {\r
84     AliFatal("AliEMCALUnfolding: Geometry not initialized.");\r
85   }\r
86   \r
87 }\r
88 \r
89 //____________________________________________________________________________\r
90 AliEMCALUnfolding::AliEMCALUnfolding(AliEMCALGeometry* geometry,Float_t ECALocMaxCut,Double_t *SSPars,Double_t *Par5,Double_t *Par6):\r
91   fNumberOfECAClusters(0),\r
92   fECALocMaxCut(ECALocMaxCut),\r
93   fThreshold(0.01),//10 MeV\r
94   fGeom(geometry),\r
95   fRecPoints(NULL),\r
96   fDigitsArr(NULL)\r
97 {\r
98   // ctor with the indication of the file where header Tree and digits Tree are stored\r
99   // use this contructor to avoid usage of Init() which uses runloader\r
100   // change needed by HLT - MP\r
101   if (!fGeom)\r
102   {\r
103     AliFatal("AliEMCALUnfolding: Geometry not initialized.");\r
104   }\r
105   Int_t i=0;\r
106   for (i = 0; i < 8; i++) fgSSPars[i] = SSPars[i];\r
107   for (i = 0; i < 3; i++) {\r
108     fgPar5[i] = Par5[i];\r
109     fgPar6[i] = Par6[i];\r
110   }\r
111 \r
112 }\r
113 \r
114 //____________________________________________________________________________\r
115 void AliEMCALUnfolding::Init()\r
116 {\r
117   // Make all memory allocations which can not be done in default constructor.\r
118   // Attach the Clusterizer task to the list of EMCAL tasks\r
119 \r
120   AliRunLoader *rl = AliRunLoader::Instance();\r
121   if (rl && rl->GetAliRun()){\r
122     AliEMCAL* emcal = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"));\r
123     if(emcal)fGeom = emcal->GetGeometry();\r
124   }\r
125   \r
126   if(!fGeom)\r
127     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());\r
128   \r
129   AliDebug(1,Form("geom %p",fGeom));\r
130   \r
131   if(!gMinuit) \r
132     gMinuit = new TMinuit(100) ;\r
133   \r
134 }\r
135 \r
136 //____________________________________________________________________________\r
137   AliEMCALUnfolding::~AliEMCALUnfolding()\r
138 {\r
139   // dtor\r
140 }\r
141 \r
142 //____________________________________________________________________________\r
143 void AliEMCALUnfolding::SetInput(Int_t numberOfECAClusters,TObjArray *recPoints,TClonesArray *digitsArr)\r
144 {\r
145   //\r
146   //Set input for unfolding purposes\r
147   SetNumberOfECAClusters(numberOfECAClusters);\r
148   SetRecPoints(recPoints);\r
149   SetDigitsArr(digitsArr);\r
150 }\r
151 \r
152 //____________________________________________________________________________\r
153 void AliEMCALUnfolding::MakeUnfolding()\r
154 {\r
155   // Unfolds clusters using the shape of an ElectroMagnetic shower\r
156   // Performs unfolding of all clusters\r
157   \r
158   if(fNumberOfECAClusters > 0){\r
159     if (fGeom==0)\r
160       AliFatal("Did not get geometry from EMCALLoader") ;\r
161     //Int_t nModulesToUnfold = fGeom->GetNCells();\r
162     \r
163     Int_t numberofNotUnfolded = fNumberOfECAClusters ;\r
164     Int_t index ;\r
165     for(index = 0 ; index < numberofNotUnfolded ; index++){\r
166       AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(index) ) ;\r
167       if(recPoint){\r
168         \r
169         Int_t nMultipl = recPoint->GetMultiplicity() ;\r
170         AliEMCALDigit ** maxAt = new AliEMCALDigit*[nMultipl] ;\r
171         Float_t * maxAtEnergy = new Float_t[nMultipl] ;\r
172         Int_t nMax = recPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fECALocMaxCut,fDigitsArr) ;\r
173         \r
174         if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0\r
175           if(UnfoldClusterV2(recPoint, nMax, maxAt, maxAtEnergy) ){\r
176             fRecPoints->Remove(recPoint);\r
177             fRecPoints->Compress() ;//is it really needed\r
178             index-- ;\r
179             fNumberOfECAClusters-- ;\r
180             numberofNotUnfolded-- ;\r
181           }\r
182         }\r
183         else{\r
184           recPoint->SetNExMax(1) ; //Only one local maximum\r
185         }\r
186         \r
187         delete[] maxAt ;\r
188         delete[] maxAtEnergy ;\r
189       } else AliError("RecPoint NULL");\r
190     } // rec point loop\r
191   }\r
192   // End of Unfolding of clusters\r
193 }\r
194 \r
195 //____________________________________________________________________________\r
196 Bool_t AliEMCALUnfolding::UnfoldClusterV2(AliEMCALRecPoint * iniTower, \r
197                                           Int_t nMax, \r
198                                           AliEMCALDigit ** maxAt, \r
199                                           Float_t * maxAtEnergy)\r
200 {\r
201   // Extended to whole EMCAL \r
202   \r
203   //**************************** part 1 *******************************************\r
204   // Performs the unfolding of a cluster with nMax overlapping showers \r
205   \r
206   //printf("Original cluster E %f\n",iniTower->GetEnergy());\r
207   \r
208   Int_t nPar = 3 * nMax ;\r
209   Float_t * fitparameters = new Float_t[nPar] ;\r
210   \r
211   if (fGeom==0)\r
212     AliFatal("Did not get geometry from EMCALLoader") ;\r
213   \r
214   Bool_t rv = FindFitV2(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;\r
215   if( !rv ) \r
216   {\r
217     // Fit failed, return (and remove cluster? - why? I leave the cluster)\r
218     iniTower->SetNExMax(-1) ;\r
219     delete[] fitparameters ;\r
220     return kFALSE;\r
221   }\r
222   \r
223   //**************************** part 2 *******************************************\r
224   // create unfolded rec points and fill them with new energy lists\r
225   // First calculate energy deposited in each sell in accordance with\r
226   // fit (without fluctuations): efit[]\r
227   // and later correct this number in acordance with actual energy\r
228   // deposition\r
229   \r
230   Int_t nDigits = iniTower->GetMultiplicity() ;\r
231   Float_t * efit = new Float_t[nDigits] ;//new fitted energy in cells\r
232   Float_t xpar=0.,zpar=0.,epar=0.  ;//center of gravity in cell units\r
233   \r
234   AliEMCALDigit * digit = 0 ;\r
235   Int_t * digitsList = iniTower->GetDigitsList() ;\r
236   \r
237   Int_t iSupMod =  0 ;\r
238   Int_t iTower  =  0 ;\r
239   Int_t iIphi   =  0 ;\r
240   Int_t iIeta   =  0 ;\r
241   Int_t iphi    =  0 ;//x direction\r
242   Int_t ieta    =  0 ;//z direstion\r
243   \r
244   Int_t iparam = 0 ;\r
245   Int_t iDigit = 0 ;\r
246   \r
247   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++)\r
248   {\r
249     digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At(digitsList[iDigit] ) ) ;\r
250     if(digit)\r
251     {\r
252       fGeom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
253       fGeom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
254                                          iIphi, iIeta,iphi,ieta);\r
255       EvalParsPhiDependence(digit->GetId(),fGeom);\r
256       \r
257       efit[iDigit] = 0.;\r
258       iparam = 0;\r
259       while(iparam < nPar )\r
260       {\r
261         xpar = fitparameters[iparam] ;\r
262         zpar = fitparameters[iparam+1] ;\r
263         epar = fitparameters[iparam+2] ;\r
264 \r
265         efit[iDigit] += epar * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) ;\r
266 \r
267         iparam += 3 ;\r
268       }\r
269 \r
270     } else AliDebug(1,"Digit NULL part 2!");\r
271     \r
272   }//digit loop\r
273   \r
274   //**************************** part 3 *******************************************\r
275   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations\r
276   // so that energy deposited in each cell is distributed between new clusters proportionally\r
277   // to its contribution to efit\r
278   \r
279   Float_t * energiesList = iniTower->GetEnergiesList() ;\r
280   Float_t ratio = 0 ;\r
281   Float_t eDigit = 0. ;\r
282   Int_t nSplittedClusters=(Int_t)nPar/3;\r
283   \r
284   Float_t * correctedEnergyList = new Float_t[nDigits*nSplittedClusters];\r
285   //above - temporary table with energies after unfolding.\r
286   //the order is following: \r
287   //first cluster <first cell - last cell>, \r
288   //second cluster <first cell - last cell>, etc.\r
289   \r
290   //**************************** sub-part 3.1 *************************************\r
291   //here we check if energy of the cell in the cluster after unfolding is above threshold. \r
292   //If not the energy from a given cell in the cluster is divided in correct proportions \r
293   //in accordance to the other clusters and added to them and set to 0.\r
294   \r
295   iparam = 0 ;\r
296   while(iparam < nPar )\r
297   {\r
298     xpar = fitparameters[iparam] ;\r
299     zpar = fitparameters[iparam+1] ;\r
300     epar = fitparameters[iparam+2] ;\r
301     \r
302     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++)\r
303     {\r
304       digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;\r
305       if(digit)\r
306       {\r
307         fGeom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
308         fGeom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
309                                            iIphi, iIeta,iphi,ieta);\r
310         \r
311         EvalParsPhiDependence(digit->GetId(),fGeom);\r
312        \r
313  \r
314         if(efit[iDigit]==0) \r
315         {//just for sure\r
316           correctedEnergyList[iparam/3*nDigits+iDigit] = 0;//correction here\r
317           continue;\r
318         }\r
319         \r
320         ratio = epar * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) / efit[iDigit] ;\r
321         eDigit = energiesList[iDigit] * ratio ;\r
322         \r
323         //add energy to temporary matrix\r
324         correctedEnergyList[iparam/3*nDigits+iDigit] = eDigit;\r
325         \r
326       } else AliDebug(1,"NULL digit part 3");\r
327     }//digit loop \r
328     iparam += 3 ;\r
329   }//while\r
330   \r
331   //**************************** sub-part 3.2 *************************************\r
332   //here we correct energy for each cell and cluster\r
333   Float_t maximumEne=0;\r
334   Int_t maximumIndex=0;\r
335   Bool_t isAnyBelowThreshold=kFALSE;\r
336   //  Float_t Threshold=0.01;\r
337   Float_t * energyFraction = new Float_t[nSplittedClusters];\r
338   Int_t iparam2 = 0 ;\r
339   for(iDigit = 0 ; iDigit < nDigits ; iDigit++)\r
340   {\r
341     isAnyBelowThreshold=kFALSE;\r
342     maximumEne=0;\r
343     for(iparam = 0 ; iparam < nPar ; iparam+=3)\r
344     {\r
345 \r
346       if(correctedEnergyList[iparam/3*nDigits+iDigit] < fThreshold ) isAnyBelowThreshold = kTRUE;\r
347       if(correctedEnergyList[iparam/3*nDigits+iDigit] > maximumEne) \r
348       {\r
349         maximumEne = correctedEnergyList[iparam/3*nDigits+iDigit];\r
350         maximumIndex = iparam;\r
351       }\r
352     }//end of loop over clusters after unfolding\r
353     \r
354     if(!isAnyBelowThreshold) continue; //no cluster-cell below threshold \r
355     \r
356     //    printf("Correct E cell %f < %f, org Digit index %d, e = %f\n",correctedEnergyList[iparam/3*nDigits+iDigit],fThreshold,iDigit, energiesList[iDigit]);\r
357     //if( energiesList[iDigit] < correctedEnergyList[iparam/3*nDigits+iDigit]) printf("\t What? \n");\r
358 \r
359     if(maximumEne < fThreshold) \r
360     {//add all cluster cells and put energy into max index, other set to 0\r
361       maximumEne=0.;\r
362       for(iparam = 0 ; iparam < nPar ; iparam+=3)\r
363       {\r
364         maximumEne+=correctedEnergyList[iparam/3*nDigits+iDigit];\r
365         correctedEnergyList[iparam/3*nDigits+iDigit]=0;\r
366       }\r
367       correctedEnergyList[maximumIndex/3*nDigits+iDigit]=maximumEne;\r
368       continue;\r
369     }//end if\r
370     \r
371     //divide energy of cell below threshold in the correct proportion and add to other cells\r
372     maximumEne=0;//not used any more so use it for the energy sum\r
373     for(iparam = 0 ; iparam < nPar ; iparam+=3)\r
374     {//calculate energy sum\r
375       if(correctedEnergyList[iparam/3*nDigits+iDigit] < fThreshold) energyFraction[iparam/3]=0;\r
376       else \r
377       {\r
378         energyFraction[iparam/3]=1;\r
379         maximumEne+=correctedEnergyList[iparam/3*nDigits+iDigit];\r
380       }\r
381     }//end of loop over clusters after unfolding\r
382     if(maximumEne>0)\r
383     {\r
384       for(iparam = 0 ; iparam < nPar ; iparam+=3){//calculate fraction\r
385         energyFraction[iparam/3] = energyFraction[iparam/3] * correctedEnergyList[iparam/3*nDigits+iDigit] / maximumEne;\r
386       }\r
387       \r
388       for(iparam = 0 ; iparam < nPar ; iparam+=3)\r
389       {//add energy from cells below threshold to others\r
390         if(energyFraction[iparam/3]>0) continue;\r
391         else\r
392         {\r
393           for(iparam2 = 0 ; iparam2 < nPar ; iparam2+=3)\r
394           {\r
395             correctedEnergyList[iparam2/3*nDigits+iDigit] += (energyFraction[iparam2/3] * \r
396                                                       correctedEnergyList[iparam/3*nDigits+iDigit]) ;\r
397           }//inner loop\r
398           correctedEnergyList[iparam/3*nDigits+iDigit] = 0;\r
399         }\r
400       }\r
401     }\r
402     else\r
403     {\r
404       //digit energy to be set to 0\r
405       for(iparam = 0 ; iparam < nPar ; iparam+=3)\r
406       {\r
407         correctedEnergyList[iparam/3*nDigits+iDigit] = 0;\r
408       }\r
409     }//new adam correction for is energy>0\r
410     \r
411   }//end of loop over digits\r
412   delete[] energyFraction;\r
413   \r
414   //**************************** sub-part 3.3 *************************************\r
415   //here we add digits to recpoints with corrected energy\r
416   iparam = 0 ;\r
417   while(iparam < nPar )\r
418   {\r
419     AliEMCALRecPoint * recPoint = 0 ;\r
420     \r
421     if(fNumberOfECAClusters >= fRecPoints->GetSize())\r
422       fRecPoints->Expand(2*fNumberOfECAClusters) ;\r
423     \r
424     //add recpoint\r
425     (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("") ;\r
426     recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(fNumberOfECAClusters) ) ;\r
427     \r
428     if(recPoint)\r
429     {\r
430       fNumberOfECAClusters++ ;\r
431       recPoint->SetNExMax(nSplittedClusters) ;\r
432       \r
433       for(iDigit = 0 ; iDigit < nDigits ; iDigit ++)\r
434       {\r
435         digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;\r
436         \r
437         \r
438         if(digit && correctedEnergyList[iparam/3*nDigits+iDigit]>0. )\r
439         {\r
440           if(correctedEnergyList[iparam/3*nDigits+iDigit]<fThreshold) printf("Final E cell %f < %f\n",correctedEnergyList[iparam/3*nDigits+iDigit],fThreshold);\r
441           recPoint->AddDigit( *digit, correctedEnergyList[iparam/3*nDigits+iDigit], kFALSE ) ; //FIXME, need to study the shared case\r
442         } else \r
443         {\r
444           AliDebug(1,Form("NULL digit part3.3 or NULL energy=%f",correctedEnergyList[iparam/3*nDigits+iDigit]));\r
445         }\r
446       }//digit loop \r
447     } else AliError("NULL RecPoint");\r
448     \r
449     //protection from recpoint with no digits\r
450     if(recPoint->GetMultiplicity()==0)\r
451     {\r
452       delete (*fRecPoints)[fNumberOfECAClusters];\r
453       //cout<<"size fRecPoints before "<<fRecPoints->GetSize()<<endl;\r
454       fRecPoints->RemoveAt(fNumberOfECAClusters);\r
455       //cout<<"size fRecPoints after "<<fRecPoints->GetSize()<<endl;\r
456       fNumberOfECAClusters--;\r
457       nSplittedClusters--;\r
458       \r
459     }\r
460     \r
461     iparam += 3 ;\r
462   }//while\r
463   \r
464   delete[] fitparameters ;\r
465   delete[] efit ;\r
466   delete[] correctedEnergyList ;\r
467   \r
468   return kTRUE;\r
469 }\r
470 \r
471 \r
472 //____________________________________________________________________________\r
473 Bool_t AliEMCALUnfolding::UnfoldClusterV2old(AliEMCALRecPoint * iniTower, \r
474                                           Int_t nMax, \r
475                                           AliEMCALDigit ** maxAt, \r
476                                           Float_t * maxAtEnergy)\r
477 {\r
478   // Extended to whole EMCAL \r
479   // Performs the unfolding of a cluster with nMax overlapping showers \r
480   \r
481   Int_t nPar = 3 * nMax ;\r
482   Float_t * fitparameters = new Float_t[nPar] ;\r
483   \r
484   if (fGeom==0)\r
485     AliFatal("Did not get geometry from EMCALLoader") ;\r
486   \r
487   Bool_t rv = FindFitV2(iniTower, maxAt, maxAtEnergy, nPar, fitparameters) ;\r
488   if( !rv ) {\r
489     // Fit failed, return (and remove cluster? - why? I leave the cluster)\r
490     iniTower->SetNExMax(-1) ;\r
491     delete[] fitparameters ;\r
492     return kFALSE;\r
493   }\r
494   \r
495   // create unfolded rec points and fill them with new energy lists\r
496   // First calculate energy deposited in each sell in accordance with\r
497   // fit (without fluctuations): efit[]\r
498   // and later correct this number in acordance with actual energy\r
499   // deposition\r
500   \r
501   Int_t nDigits = iniTower->GetMultiplicity() ;\r
502   Float_t * efit = new Float_t[nDigits] ;//new fitted energy in cells\r
503   Float_t xpar=0.,zpar=0.,epar=0.  ;//center of gravity in cell units\r
504   \r
505   AliEMCALDigit * digit = 0 ;\r
506   Int_t * digitsList = iniTower->GetDigitsList() ;\r
507   \r
508   Int_t iSupMod =  0 ;\r
509   Int_t iTower  =  0 ;\r
510   Int_t iIphi   =  0 ;\r
511   Int_t iIeta   =  0 ;\r
512   Int_t iphi    =  0 ;//x direction\r
513   Int_t ieta    =  0 ;//z direstion\r
514   \r
515   Int_t iparam = 0 ;\r
516   Int_t iDigit = 0 ;\r
517   \r
518   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){\r
519     digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At(digitsList[iDigit] ) ) ;\r
520     if(digit){\r
521       fGeom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
522       fGeom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
523                                          iIphi, iIeta,iphi,ieta);\r
524       EvalParsPhiDependence(digit->GetId(),fGeom);\r
525       \r
526       efit[iDigit] = 0.;\r
527       iparam = 0;\r
528       while(iparam < nPar ){\r
529         xpar = fitparameters[iparam] ;\r
530         zpar = fitparameters[iparam+1] ;\r
531         epar = fitparameters[iparam+2] ;\r
532         iparam += 3 ;\r
533         \r
534         efit[iDigit] += epar * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) ;\r
535       }\r
536     } else AliError("Digit NULL!");\r
537     \r
538   }//digit loop\r
539   \r
540   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations\r
541   // so that energy deposited in each cell is distributed between new clusters proportionally\r
542   // to its contribution to efit\r
543   \r
544   Float_t * energiesList = iniTower->GetEnergiesList() ;\r
545   Float_t ratio = 0 ;\r
546   \r
547   iparam = 0 ;\r
548   while(iparam < nPar ){\r
549     xpar = fitparameters[iparam] ;\r
550     zpar = fitparameters[iparam+1] ;\r
551     epar = fitparameters[iparam+2] ;\r
552     iparam += 3 ;\r
553     \r
554     AliEMCALRecPoint * recPoint = 0 ;\r
555     \r
556     if(fNumberOfECAClusters >= fRecPoints->GetSize())\r
557       fRecPoints->Expand(2*fNumberOfECAClusters) ;\r
558     \r
559     //add recpoint\r
560     (*fRecPoints)[fNumberOfECAClusters] = new AliEMCALRecPoint("") ;\r
561     recPoint = dynamic_cast<AliEMCALRecPoint *>( fRecPoints->At(fNumberOfECAClusters) ) ;\r
562     \r
563     if(recPoint){\r
564       \r
565       fNumberOfECAClusters++ ;\r
566       recPoint->SetNExMax((Int_t)nPar/3) ;\r
567       \r
568       Float_t eDigit = 0. ;\r
569       for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){\r
570         digit = dynamic_cast<AliEMCALDigit*>( fDigitsArr->At( digitsList[iDigit] ) ) ;\r
571         if(digit){\r
572           fGeom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
573           fGeom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
574                                              iIphi, iIeta,iphi,ieta);\r
575           EvalParsPhiDependence(digit->GetId(),fGeom);\r
576           if(efit[iDigit]==0) continue;//just for sure\r
577           ratio = epar * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) / efit[iDigit] ;\r
578           eDigit = energiesList[iDigit] * ratio ;\r
579           recPoint->AddDigit( *digit, eDigit, kFALSE ) ; //FIXME, need to study the shared case\r
580         } else AliError("NULL digit");\r
581       }//digit loop \r
582     } else AliError("NULL RecPoint");\r
583   }//while\r
584   \r
585   delete[] fitparameters ;\r
586   delete[] efit ;\r
587   \r
588   return kTRUE;\r
589 }\r
590 \r
591 \r
592 //____________________________________________________________________________\r
593 Bool_t AliEMCALUnfolding::FindFitV2(AliEMCALRecPoint * recPoint, AliEMCALDigit ** maxAt, \r
594                                         const Float_t* maxAtEnergy,\r
595                                         Int_t nPar, Float_t * fitparameters) const\r
596 {\r
597   // Calls TMinuit to fit the energy distribution of a cluster with several maxima\r
598   // The initial values for fitting procedure are set equal to the\r
599   // positions of local maxima.       \r
600   // Cluster will be fitted as a superposition of nPar/3\r
601   // electromagnetic showers\r
602 \r
603   if (fGeom==0) AliFatal("Did not get geometry from EMCALLoader");\r
604         \r
605   if(!gMinuit)\r
606     gMinuit = new TMinuit(100) ;//max 100 parameters\r
607 \r
608   gMinuit->mncler();                     // Reset Minuit's list of paramters\r
609   gMinuit->SetPrintLevel(-1) ;           // No Printout\r
610   gMinuit->SetFCN(AliEMCALUnfolding::UnfoldingChiSquareV2) ;\r
611   // To set the address of the minimization function\r
612   TList * toMinuit = new TList();\r
613   toMinuit->AddAt(recPoint,0) ;\r
614   toMinuit->AddAt(fDigitsArr,1) ;\r
615   toMinuit->AddAt(fGeom,2) ;\r
616 \r
617   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare\r
618 \r
619   // filling initial values for fit parameters\r
620   AliEMCALDigit * digit ;\r
621 \r
622   Int_t ierflg  = 0;\r
623   Int_t index   = 0 ;\r
624   Int_t nDigits = (Int_t) nPar / 3 ;\r
625 \r
626   Int_t iDigit ;\r
627 \r
628   Int_t iSupMod =  0 ;\r
629   Int_t iTower  =  0 ;\r
630   Int_t iIphi   =  0 ;\r
631   Int_t iIeta   =  0 ;\r
632   Int_t iphi    =  0 ;//x direction\r
633   Int_t ieta    =  0 ;//z direstion\r
634 \r
635   for(iDigit = 0; iDigit < nDigits; iDigit++){\r
636     digit = maxAt[iDigit];\r
637     if(digit==0) AliError("energy of digit = 0!");\r
638     fGeom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
639     fGeom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
640                                        iIphi, iIeta,iphi,ieta);\r
641 \r
642     Float_t energy = maxAtEnergy[iDigit] ;\r
643 \r
644     //gMinuit->mnparm(index, "x",  iphi, 0.1, 0, 0, ierflg) ;//original\r
645     gMinuit->mnparm(index, "x",  iphi, 0.05, 0, 0, ierflg) ;\r
646     index++ ;\r
647     if(ierflg != 0){\r
648       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %d", iphi ) ;\r
649       toMinuit->Clear();\r
650       delete toMinuit ;\r
651       return kFALSE;\r
652     }\r
653     //gMinuit->mnparm(index, "z",  ieta, 0.1, 0, 0, ierflg) ;//original\r
654     gMinuit->mnparm(index, "z",  ieta, 0.05, 0, 0, ierflg) ;\r
655     index++ ;\r
656     if(ierflg != 0){\r
657       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %d", ieta) ;\r
658       toMinuit->Clear();\r
659       delete toMinuit ;\r
660       return kFALSE;\r
661     }\r
662     //gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;//original\r
663     gMinuit->mnparm(index, "Energy",  energy , 0.001*energy, 0., 5.*energy, ierflg) ;//was 0.05\r
664     index++ ;\r
665     if(ierflg != 0){\r
666       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;\r
667       toMinuit->Clear();\r
668       delete toMinuit ;\r
669       return kFALSE;\r
670     }\r
671   }\r
672 \r
673   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; \r
674                       // The number of function call slightly depends on it.\r
675   //  Double_t p1 = 1.0 ;// par to gradient \r
676   Double_t p2 = 0.0 ;\r
677   //  Double_t p3 = 3.0 ;\r
678   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls\r
679   //  gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient\r
680   gMinuit->SetMaxIterations(5);//was 5\r
681   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings\r
682   //gMinuit->mnexcm("SET PRI", &p3 , 3, ierflg) ;  // printouts\r
683 \r
684   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize\r
685   //gMinuit->mnexcm("MINI", &p0, 0, ierflg) ;    // minimize\r
686   if(ierflg == 4){  // Minimum not found\r
687     AliDebug(1,"EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;\r
688     toMinuit->Clear();\r
689     delete toMinuit ;\r
690     return kFALSE ;\r
691   }\r
692   for(index = 0; index < nPar; index++){\r
693     Double_t err = 0. ;\r
694     Double_t val = 0. ;\r
695     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index\r
696     fitparameters[index] = val ;\r
697   }\r
698 \r
699   toMinuit->Clear();\r
700   delete toMinuit ;\r
701   return kTRUE;\r
702 \r
703 }\r
704 \r
705 //____________________________________________________________________________\r
706 Double_t  AliEMCALUnfolding::ShowerShapeV2(Double_t x, Double_t y)\r
707\r
708   // extended to whole EMCAL \r
709   // Shape of the shower\r
710   // If you change this function, change also the gradient evaluation in ChiSquare()\r
711 \r
712   Double_t r = fgSSPars[7]*TMath::Sqrt(x*x+y*y);\r
713   Double_t rp1  = TMath::Power(r, fgSSPars[1]) ;\r
714   Double_t rp5  = TMath::Power(r, fgSSPars[5]) ;\r
715   Double_t shape = fgSSPars[0]*TMath::Exp( -rp1 * (1. / (fgSSPars[2] + fgSSPars[3] * rp1) + fgSSPars[4] / (1 + fgSSPars[6] * rp5) ) ) ;\r
716   return shape ;\r
717 }\r
718 \r
719 //____________________________________________________________________________\r
720 void AliEMCALUnfolding::UnfoldingChiSquareV2(Int_t & nPar, Double_t * Grad,\r
721                                              Double_t & fret,\r
722                                              Double_t * x, Int_t iflag)\r
723 {\r
724   // Calculates the Chi square for the cluster unfolding minimization\r
725   // Number of parameters, Gradient, Chi squared, parameters, what to do\r
726   \r
727   TList * toMinuit = dynamic_cast<TList*>( gMinuit->GetObjectFit() ) ;\r
728   if(toMinuit){\r
729     AliEMCALRecPoint * recPoint = dynamic_cast<AliEMCALRecPoint*>( toMinuit->At(0) )  ;\r
730     TClonesArray * digits = dynamic_cast<TClonesArray*>( toMinuit->At(1) )  ;\r
731     // A bit buggy way to get an access to the geometry\r
732     // To be revised!\r
733     AliEMCALGeometry *geom = dynamic_cast<AliEMCALGeometry *>(toMinuit->At(2));\r
734     \r
735     if(recPoint && digits && geom){\r
736       \r
737       Int_t * digitsList     = recPoint->GetDigitsList() ;\r
738       \r
739       Int_t nOdigits = recPoint->GetDigitsMultiplicity() ;\r
740       \r
741       Float_t * energiesList = recPoint->GetEnergiesList() ;\r
742       \r
743       fret = 0. ;\r
744       Int_t iparam = 0 ;\r
745       \r
746       if(iflag == 2)\r
747         for(iparam = 0 ; iparam < nPar ; iparam++)\r
748           Grad[iparam] = 0 ; // Will evaluate gradient\r
749       \r
750       Double_t efit = 0. ;\r
751       \r
752       AliEMCALDigit * digit ;\r
753       Int_t iDigit ;\r
754       \r
755       Int_t iSupMod =  0 ;\r
756       Int_t iTower  =  0 ;\r
757       Int_t iIphi   =  0 ;\r
758       Int_t iIeta   =  0 ;\r
759       Int_t iphi    =  0 ;//x direction\r
760       Int_t ieta    =  0 ;//z direstion\r
761       \r
762       \r
763       for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {\r
764         if(energiesList[iDigit]==0) continue;\r
765         \r
766         digit = dynamic_cast<AliEMCALDigit*>( digits->At( digitsList[iDigit] ) );\r
767         \r
768         if(digit){\r
769         geom->GetCellIndex(digit->GetId(),iSupMod,iTower,iIphi,iIeta); \r
770         geom->GetCellPhiEtaIndexInSModule(iSupMod,iTower,\r
771                                           iIphi, iIeta,iphi,ieta);\r
772         EvalParsPhiDependence(digit->GetId(),geom);\r
773         \r
774         if(iflag == 2){  // calculate gradient\r
775           Int_t iParam = 0 ;\r
776           efit = 0. ;\r
777           while(iParam < nPar ){\r
778             Double_t dx = ((Float_t)iphi - x[iParam]) ;\r
779             iParam++ ;\r
780             Double_t dz = ((Float_t)ieta - x[iParam]) ;\r
781             iParam++ ;\r
782             efit += x[iParam] * ShowerShapeV2(dx,dz) ;\r
783             iParam++ ;\r
784           }\r
785           \r
786           Double_t sum = 2. * (efit - energiesList[iDigit]) / energiesList[iDigit] ; // Here we assume, that sigma = sqrt(E)\r
787           iParam = 0 ;\r
788           while(iParam < nPar ){\r
789             Double_t xpar = x[iParam] ;\r
790             Double_t zpar = x[iParam+1] ;\r
791             Double_t epar = x[iParam+2] ;\r
792             \r
793             Double_t dr = fgSSPars[7]*TMath::Sqrt( ((Float_t)iphi - xpar) * ((Float_t)iphi - xpar) + ((Float_t)ieta - zpar) * ((Float_t)ieta - zpar) );\r
794             Double_t shape = sum * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) ;\r
795             Double_t rp1  = TMath::Power(dr, fgSSPars[1]) ;\r
796             Double_t rp5  = TMath::Power(dr, fgSSPars[5]) ;\r
797             \r
798             Double_t deriv = -2 * TMath::Power(dr,fgSSPars[1]-2.) * fgSSPars[7] * fgSSPars[7] * \r
799             (fgSSPars[1] * ( 1/(fgSSPars[2]+fgSSPars[3]*rp1) + fgSSPars[4]/(1+fgSSPars[6]*rp5) ) - \r
800              (fgSSPars[1]*fgSSPars[3]*rp1/( (fgSSPars[2]+fgSSPars[3]*rp1)*(fgSSPars[2]+fgSSPars[3]*rp1) ) + \r
801               fgSSPars[4]*fgSSPars[5]*fgSSPars[6]*rp5/( (1+fgSSPars[6]*rp5)*(1+fgSSPars[6]*rp5) ) ) );\r
802             \r
803             //Double_t deriv =-1.33 * TMath::Power(dr,0.33)*dr * ( 1.57 / ( (1.57 + 0.0860 * r133) * (1.57 + 0.0860 * r133) )\r
804             //                                                   - 0.55 / (1 + 0.000563 * r669) / ( (1 + 0.000563 * r669) * (1 + 0.000563 * r669) ) ) ;\r
805             \r
806             Grad[iParam] += epar * shape * deriv * ((Float_t)iphi - xpar) ;  // Derivative over x\r
807             iParam++ ;\r
808             Grad[iParam] += epar * shape * deriv * ((Float_t)ieta - zpar) ;  // Derivative over z\r
809             iParam++ ;\r
810             Grad[iParam] += shape ;                                  // Derivative over energy\r
811             iParam++ ;\r
812           }\r
813         }\r
814         efit = 0;\r
815         iparam = 0 ;\r
816         \r
817         while(iparam < nPar ){\r
818           Double_t xpar = x[iparam] ;\r
819           Double_t zpar = x[iparam+1] ;\r
820           Double_t epar = x[iparam+2] ;\r
821           iparam += 3 ;\r
822           efit += epar * ShowerShapeV2((Float_t)iphi - xpar,(Float_t)ieta - zpar) ;\r
823         }\r
824         \r
825         fret += (efit-energiesList[iDigit])*(efit-energiesList[iDigit])/energiesList[iDigit] ;\r
826         // Here we assume, that sigma = sqrt(E) \r
827         } else printf("AliEMCALUnfoding::UnfoldingChiSquareV2 - NULL digit!, nPar %d \n", nPar); // put nPar here to cheat coverity and rule checker\r
828       } // digit loop\r
829     } // recpoint, digits and geom not NULL\r
830   }// List is not NULL\r
831   \r
832 }\r
833 \r
834 \r
835 //____________________________________________________________________________\r
836 void AliEMCALUnfolding::SetShowerShapeParams(Double_t *pars){\r
837   for(UInt_t i=0;i<7;++i)\r
838     fgSSPars[i]=pars[i];\r
839   if(pars[2]==0. && pars[3]==0.) fgSSPars[2]=1.;//to avoid dividing by 0\r
840 }\r
841 \r
842 //____________________________________________________________________________\r
843 void AliEMCALUnfolding::SetPar5(Double_t *pars){\r
844   for(UInt_t i=0;i<3;++i)\r
845     fgPar5[i]=pars[i];\r
846 }\r
847 \r
848 //____________________________________________________________________________\r
849 void AliEMCALUnfolding::SetPar6(Double_t *pars){\r
850   for(UInt_t i=0;i<3;++i)\r
851     fgPar6[i]=pars[i];\r
852 }\r
853 \r
854 //____________________________________________________________________________\r
855 void AliEMCALUnfolding::EvalPar5(Double_t phi){\r
856   //\r
857   //Evaluate the 5th parameter of the shower shape function\r
858   //phi in degrees range (-10,10)\r
859   //\r
860   //fSSPars[5] = 12.31 - phi*0.007381 - phi*phi*0.06936;\r
861   fgSSPars[5] = fgPar5[0] + phi * fgPar5[1] + phi*phi * fgPar5[2];\r
862 }\r
863 \r
864 //____________________________________________________________________________\r
865 void AliEMCALUnfolding::EvalPar6(Double_t phi){\r
866   //\r
867   //Evaluate the 6th parameter of the shower shape function\r
868   //phi in degrees range (-10,10)\r
869   //\r
870   //fSSPars[6] = 0.05452 + phi*0.0001228 + phi*phi*0.001361;\r
871   fgSSPars[6] = fgPar6[0] + phi * fgPar6[1] + phi*phi * fgPar6[2];\r
872 }\r
873 \r
874 //____________________________________________________________________________\r
875 void AliEMCALUnfolding::EvalParsPhiDependence(Int_t absId, const AliEMCALGeometry *geom){\r
876   //\r
877   // calculate params p5 and p6 depending on the phi angle in global coordinate\r
878   // for the cell with given absId index\r
879   //\r
880   Double_t etaGlob = 0.;//eta in global c.s. - unused\r
881   Double_t phiGlob = 0.;//phi in global c.s. in radians\r
882   geom->EtaPhiFromIndex(absId, etaGlob, phiGlob);\r
883   phiGlob*=180./TMath::Pi();\r
884   phiGlob-=90.;\r
885   phiGlob-= (Double_t)((Int_t)geom->GetSuperModuleNumber(absId)/2 * 20);\r
886 \r
887   EvalPar5(phiGlob);\r
888   EvalPar6(phiGlob);\r
889 }\r
890 \r