]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliPHOSTenderSupply.cxx
New version of Add task tender allowing to customize the components
[u/mrichter/AliRoot.git] / ANALYSIS / TenderSupplies / AliPHOSTenderSupply.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  PHOS tender, recalibrate PHOS clusters                                   //
20 //  and do track matching                                                    //
21 //  Author : Dmitri Peressounko (RRC KI)                                     //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include "TROOT.h"
26 #include "TH2.h"
27 #include "TFile.h"
28
29 #include <AliLog.h>
30 #include <AliESDEvent.h>
31 #include <AliAnalysisManager.h>
32 #include <AliTender.h>
33 #include <AliCDBManager.h>
34 #include "AliMagF.h"
35 #include "TGeoGlobalMagField.h"
36
37 #include "AliESDCaloCluster.h"
38 #include "AliPHOSTenderSupply.h"
39 #include "AliPHOSCalibData.h"
40 #include "AliPHOSGeometry.h"
41 #include "AliPHOSEsdCluster.h"
42 #include "AliOADBContainer.h"
43
44 ClassImp(AliPHOSTenderSupply)
45
46 AliPHOSTenderSupply::AliPHOSTenderSupply() :
47   AliTenderSupply()
48   ,fOCDBpass("local://OCDB")
49   ,fNonlinearityVersion("Default")
50   ,fPHOSGeo(0x0)
51   ,fRunNumber(0)
52   ,fRecoPass(-1)  //to be defined
53   ,fUsePrivateBadMap(0)
54   ,fUsePrivateCalib(0)
55   ,fPHOSCalibData(0x0)
56 {
57         //
58         // default ctor
59         //
60    for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
61    for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
62 }
63
64 //_____________________________________________________
65 AliPHOSTenderSupply::AliPHOSTenderSupply(const char *name, const AliTender *tender) :
66   AliTenderSupply(name,tender)
67   ,fOCDBpass("alien:///alice/cern.ch/user/p/prsnko/PHOSrecalibrations/")
68   ,fNonlinearityVersion("Default")
69   ,fPHOSGeo(0x0)
70   ,fRunNumber(0)
71   ,fRecoPass(-1)  //to be defined
72   ,fUsePrivateBadMap(0)
73   ,fUsePrivateCalib(0)
74   ,fPHOSCalibData(0x0)
75 {
76         //
77         // named ctor
78         //
79    for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
80    for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
81 }
82
83 //_____________________________________________________
84 AliPHOSTenderSupply::~AliPHOSTenderSupply()
85 {
86   //Destructor
87   if(fPHOSCalibData)
88     delete fPHOSCalibData;
89   fPHOSCalibData=0x0 ;
90 }
91
92 //_____________________________________________________
93 void AliPHOSTenderSupply::InitTender()
94 {
95   //
96   // Initialise PHOS tender
97   //
98   AliESDEvent *event=fTender->GetEvent();
99   if (!event) return;
100               
101   fRunNumber=event->GetRunNumber();           
102   if(fRecoPass<0){ //not defined yet
103     // read if from filename.
104     AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
105     TTree * t = mgr->GetTree();
106     if(t){  
107       TFile * f = t->GetCurrentFile() ;
108       if(f){  
109         TString fname(f->GetName());
110         if(fname.Contains("pass1"))
111            fRecoPass=1;
112         else 
113           if(fname.Contains("pass2"))
114            fRecoPass=2;
115           else 
116             if(fname.Contains("pass3")) 
117               fRecoPass=3;
118       }
119     }
120     if(fRecoPass<0){
121       AliError("Can not find pass number from file name, set it manually");
122     }
123   }
124    
125   //Init geometry 
126   if(!fPHOSGeo){
127     AliOADBContainer geomContainer("phosGeo");
128     geomContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSGeometry.root","PHOSRotationMatrixes");
129     TObjArray *matrixes = (TObjArray*)geomContainer.GetObject(fRunNumber,"PHOSRotationMatrixes");
130     fPHOSGeo =  AliPHOSGeometry::GetInstance("IHEP") ;
131     for(Int_t mod=0; mod<5; mod++) {
132       if(!matrixes->At(mod)) continue;
133       fPHOSGeo->SetMisalMatrix(((TGeoHMatrix*)matrixes->At(mod)),mod) ;
134       printf(".........Adding Matrix(%d), geo=%p\n",mod,fPHOSGeo) ;
135     }
136   }
137   
138   //Init Bad channels map
139   if(!fUsePrivateBadMap){
140    AliOADBContainer badmapContainer(Form("phosBadMap"));
141     badmapContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSBadMaps.root","phosBadMap");
142     TObjArray *maps = (TObjArray*)badmapContainer.GetObject(fRunNumber,"phosBadMap");
143     if(!maps){
144       AliError(Form("Can not read Bad map for run %d. \n You may choose to use your map with ForceUsingBadMap()\n",fRunNumber)) ;    
145     }
146     else{
147       AliInfo(Form("Setting PHOS bad map with name %s \n",maps->GetName())) ;
148       for(Int_t mod=0; mod<5;mod++){
149         if(fPHOSBadMap[mod]) 
150           delete fPHOSBadMap[mod] ;
151         TH2I * h = (TH2I*)maps->At(mod) ;      
152         if(h)
153           fPHOSBadMap[mod]=new TH2I(*h) ;
154       }
155     }    
156   } 
157
158   if(!fUsePrivateCalib){
159     //Init recalibration
160     //Check the pass1-pass2-pass3 reconstruction
161     AliOADBContainer calibContainer("phosRecalibration");
162     calibContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSCalibrations.root","phosRecalibration");
163     TObjArray *recalib = (TObjArray*)calibContainer.GetObject(fRunNumber,"PHOSRecalibration");
164     if(!recalib){
165       AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",fRunNumber)) ;
166     }
167     else{
168       fPHOSCalibData = (AliPHOSCalibData*)recalib->At(fRecoPass-1) ;
169       if(!fPHOSCalibData) {
170         AliFatal(Form("Can not find calibration for run %d, pass %d \n",fRunNumber, fRecoPass)) ;
171       }
172     }
173   }
174   
175 }
176
177 //_____________________________________________________
178 void AliPHOSTenderSupply::ProcessEvent()
179 {
180   //Choose PHOS clusters and recalibrate them
181   //that it recalculate energy, position and distance 
182   //to closest track extrapolation      
183
184   AliESDEvent *event=fTender->GetEvent();
185   if (!event) return;
186               
187   fRunNumber=event->GetRunNumber();           
188
189   if(!fPHOSCalibData || fTender->RunChanged()){
190     InitTender();
191     
192   }
193
194
195   const AliESDVertex *esdVertex = event->GetPrimaryVertex();
196   AliESDCaloCells * cells = event->GetPHOSCells() ;
197   TVector3 vertex(esdVertex->GetX(),esdVertex->GetY(),esdVertex->GetZ());
198   if(vertex.Mag()>99.) //vertex not defined?
199     vertex.SetXYZ(0.,0.,0.) ;
200
201   //For re-calibration
202   const Double_t logWeight=4.5 ;  
203
204   Int_t multClust = event->GetNumberOfCaloClusters();
205   for (Int_t i=0; i<multClust; i++) {
206     AliESDCaloCluster *clu = event->GetCaloCluster(i);
207     if ( !clu->IsPHOS()) continue;
208     
209     Float_t  position[3];
210     clu->GetPosition(position);
211     TVector3 global(position) ;
212     Int_t relId[4] ;
213     fPHOSGeo->GlobalPos2RelId(global,relId) ;
214     Int_t mod  = relId[0] ;
215     Int_t cellX = relId[2];
216     Int_t cellZ = relId[3] ;
217     if ( !IsGoodChannel(mod,cellX,cellZ) ) {
218       clu->SetE(0.) ;
219       continue ;
220     }  
221    //Apply re-Calibreation
222     AliPHOSEsdCluster cluPHOS1(*clu);
223     cluPHOS1.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
224     cluPHOS1.EvalAll(logWeight,vertex);         // recalculate the cluster parameters
225     cluPHOS1.SetE(CorrectNonlinearity(cluPHOS1.E()));// Users's nonlinearity
226     
227     Float_t  xyz[3];
228     cluPHOS1.GetPosition(xyz);
229     clu->SetPosition(xyz);                       //rec.point position in MARS
230     clu->SetE(cluPHOS1.E());                           //total or core particle energy
231     clu->SetDispersion(cluPHOS1.GetDispersion());  //cluster dispersion
232     //    ec->SetPID(rp->GetPID()) ;            //array of particle identification
233     clu->SetM02(cluPHOS1.GetM02()) ;               //second moment M2x
234     clu->SetM20(cluPHOS1.GetM20()) ;               //second moment M2z
235     Double_t dx=0.,dz=0. ;
236     fPHOSGeo->GlobalPos2RelId(global,relId) ;
237     TVector3 locPos;
238     fPHOSGeo->Global2Local(locPos,global,mod) ;
239
240     Double_t pttrack=0.;
241     Int_t charge=0;
242     FindTrackMatching(mod,&locPos,dx,dz,pttrack,charge) ;
243     Double_t r=TestCPV(dx, dz, pttrack,charge) ;
244     clu->SetTrackDistance(dx,dz); 
245     
246     clu->SetEmcCpvDistance(r);    
247     clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02()));                     //not yet implemented
248     clu->SetTOF(cluPHOS1.GetTOF());       
249
250   }
251
252
253 }
254 //___________________________________________________________________________________________________
255 void AliPHOSTenderSupply::FindTrackMatching(Int_t mod,TVector3 *locpos,
256                                             Double_t &dx, Double_t &dz,
257                                             Double_t &pt,Int_t &charge){
258   //Find track with closest extrapolation to cluster
259   
260   AliESDEvent *event= fTender->GetEvent();
261   Double_t  magF = event->GetMagneticField();
262   Double_t magSign = 1.0;
263   if(magF<0)magSign = -1.0;
264   
265   if (!TGeoGlobalMagField::Instance()->GetField()) {
266     AliMagF* field = new AliMagF("Maps","Maps", magSign, magSign, AliMagF::k5kG);
267     TGeoGlobalMagField::Instance()->SetField(field);
268   }
269
270   // *** Start the matching
271   Int_t nt=event->GetNumberOfTracks();
272   //Calculate actual distance to PHOS module
273   TVector3 globaPos ;
274   fPHOSGeo->Local2Global(mod, 0.,0., globaPos) ;
275   const Double_t rPHOS = globaPos.Pt() ; //Distance to center of  PHOS module
276   const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
277   const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
278   const Double_t kAlpha0=330./180.*TMath::Pi() ; //First PHOS module angular direction
279   const Double_t kAlpha= 20./180.*TMath::Pi() ; //PHOS module angular size
280   Double_t minDistance = 1.e6;
281
282
283   Double_t gposTrack[3] ; 
284
285   Double_t bz = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField();
286   bz = TMath::Sign(0.5*kAlmost0Field,bz) + bz;
287
288   Double_t b[3]; 
289   for (Int_t i=0; i<nt; i++) {
290     AliESDtrack *esdTrack=event->GetTrack(i);
291
292     // Skip the tracks having "wrong" status (has to be checked/tuned)
293     ULong_t status = esdTrack->GetStatus();
294     if ((status & AliESDtrack::kTPCout)   == 0) continue;
295     //     if ((status & AliESDtrack::kTRDout)   == 0) continue;
296     //     if ((status & AliESDtrack::kTRDrefit) == 1) continue;
297     
298     //Continue extrapolation from TPC outer surface
299     const AliExternalTrackParam *outerParam=esdTrack->GetOuterParam();
300     if (!outerParam) continue;
301     AliExternalTrackParam t(*outerParam);
302     
303     t.GetBxByBz(b) ;
304     //Direction to the current PHOS module
305     Double_t phiMod=kAlpha0-kAlpha*mod ;
306     if(!t.Rotate(phiMod))
307       continue ;
308     
309     Double_t y;                       // Some tracks do not reach the PHOS
310     if (!t.GetYAt(rPHOS,bz,y)) continue; //    because of the bending
311     
312     Double_t z; 
313     if(!t.GetZAt(rPHOS,bz,z))
314       continue ;
315     if (TMath::Abs(z) > kZmax) 
316       continue; // Some tracks miss the PHOS in Z
317     if(TMath::Abs(y) < kYmax){
318       t.PropagateToBxByBz(rPHOS,b);        // Propagate to the matching module
319       //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
320       t.GetXYZ(gposTrack) ;
321       TVector3 globalPositionTr(gposTrack) ;
322       TVector3 localPositionTr ;
323       fPHOSGeo->Global2Local(localPositionTr,globalPositionTr,mod) ;
324       Double_t ddx = locpos->X()-localPositionTr.X();
325       Double_t ddz = locpos->Z()-localPositionTr.Z();
326       Double_t d2 = ddx*ddx + ddz*ddz;
327       if(d2 < minDistance) {
328         dx = ddx ;
329         dz = ddz ;
330         minDistance=d2 ;
331         pt=esdTrack->Pt() ;
332         charge=esdTrack->Charge() ;
333       }
334     }
335   } //Scanned all tracks
336   
337 }
338 //____________________________________________________________
339 Float_t AliPHOSTenderSupply::CorrectNonlinearity(Float_t en){
340
341   //For backward compatibility, if no RecoParameters found
342   if(fNonlinearityVersion=="Default"){
343     return 0.0241+1.0504*en+0.000249*en*en ;
344   }
345
346   if(fNonlinearityVersion=="NoCorrection"){
347     return en ;
348   }
349   if(fNonlinearityVersion=="Gustavo2005"){
350     return fNonlinearityParams[0]+fNonlinearityParams[1]*en + fNonlinearityParams[2]*en*en ;
351   }
352   if(fNonlinearityVersion=="Henrik2010"){
353     return en*(fNonlinearityParams[0]+fNonlinearityParams[1]*TMath::Exp(-en*fNonlinearityParams[2]))*(1.+fNonlinearityParams[3]*TMath::Exp(-en*fNonlinearityParams[4]))*(1.+fNonlinearityParams[6]/(en*en+fNonlinearityParams[5])) ;
354   }
355
356   return en ;
357 }
358 //_____________________________________________________________________________
359 Double_t AliPHOSTenderSupply::TestLambda(Double_t pt,Double_t l1,Double_t l2){
360   
361   Double_t l2Mean  = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
362   Double_t l1Mean  = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
363   Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
364   Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
365   Double_t c=-0.35-0.550*TMath::Exp(-0.390730*pt) ;
366   Double_t R2=0.5*(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma + 
367               0.5*(l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
368               0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
369   return (R2<2.5*2.5) ;
370   
371 }
372 //____________________________________________________________________________
373 Double_t AliPHOSTenderSupply::TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge){
374   //Parameterization of LHC10h period
375   //_true if neutral_
376   
377   Double_t meanX=0;
378   Double_t meanZ=0.;
379   Double_t sx=TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
380               6.58365e-01*5.91917e-01*5.91917e-01/((pt-9.61306e-01)*(pt-9.61306e-01)+5.91917e-01*5.91917e-01)+1.59219);
381   Double_t sz=TMath::Min(2.75,4.90341e+02*1.91456e-02*1.91456e-02/(pt*pt+1.91456e-02*1.91456e-02)+1.60) ;
382   AliESDEvent *event=fTender->GetEvent();
383   Double_t mf = event->GetMagneticField(); //Positive for ++ and negative for --
384
385   if(mf<0.){ //field --
386     meanZ = -0.468318 ;
387     if(charge>0)
388       meanX=TMath::Min(7.3, 3.89994*1.20679*1.20679/(pt*pt+1.20679*1.20679)+0.249029+2.49088e+07*TMath::Exp(-pt*3.33650e+01)) ;
389     else
390       meanX=-TMath::Min(7.7,3.86040*0.912499*0.912499/(pt*pt+0.912499*0.912499)+1.23114+4.48277e+05*TMath::Exp(-pt*2.57070e+01)) ;
391   }
392   else{ //Field ++
393     meanZ= -0.468318;
394     if(charge>0)
395       meanX=-TMath::Min(8.0,3.86040*1.31357*1.31357/(pt*pt+1.31357*1.31357)+0.880579+7.56199e+06*TMath::Exp(-pt*3.08451e+01)) ;
396     else
397       meanX= TMath::Min(6.85, 3.89994*1.16240*1.16240/(pt*pt+1.16240*1.16240)-0.120787+2.20275e+05*TMath::Exp(-pt*2.40913e+01)) ;     
398   }
399
400   Double_t rz=(dz-meanZ)/sz ;
401   Double_t rx=(dx-meanX)/sx ;
402   return TMath::Sqrt(rx*rx+rz*rz) ;
403 }
404
405 //________________________________________________________________________
406 Bool_t AliPHOSTenderSupply::IsGoodChannel(Int_t mod, Int_t ix, Int_t iz)
407 {
408   //Check if this channel belogs to the good ones
409   
410   if(mod>4 || mod<1){
411 //    AliError(Form("No bad map for PHOS module %d ",mod)) ;
412     return kTRUE ;
413   }
414   if(!fPHOSBadMap[mod]){
415 //     AliError(Form("No Bad map for PHOS module %d",mod)) ;
416      return kTRUE ;
417   }
418   if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0)
419     return kFALSE ;
420   else
421     return kTRUE ;
422 }
423 //________________________________________________________________________
424 void AliPHOSTenderSupply::ForceUsingBadMap(const char * filename){
425   //Read TH2I histograms with bad maps from local or alien file 
426   TFile * fbm = TFile::Open(filename) ;
427   if(!fbm || !fbm->IsOpen()){
428     AliError(Form("Can not open BadMaps file %s \n",filename)) ;
429     return ;
430   }
431   gROOT->cd() ;
432   char key[55] ;
433   for(Int_t mod=1;mod<4; mod++){
434     snprintf(key,55,"PHOS_BadMap_mod%d",mod) ;
435     TH2I * h = (TH2I*)fbm->Get(key) ;
436     if(h)
437        fPHOSBadMap[mod] = new TH2I(*h) ;
438   }
439   fbm->Close() ;
440   fUsePrivateBadMap=kTRUE ;
441 }
442 //________________________________________________________________________
443 void AliPHOSTenderSupply::ForceUsingCalibration(const char * filename){
444   //Read PHOS recalibration parameters from the file.
445   //We assume that file contains single entry: AliPHOSCalibData
446   TFile * fc = TFile::Open(filename) ;
447   if(!fc || !fc->IsOpen()){
448     AliFatal(Form("Can not open Calibration file %s \n",filename)) ;
449     return ;
450   }
451   fPHOSCalibData = (AliPHOSCalibData*)fc->Get("PHOSCalibration") ;
452   fc->Close() ;
453   fUsePrivateCalib=kTRUE; 
454 }
455
456
457
458
459