]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliPHOSTenderSupply.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / ANALYSIS / TenderSupplies / AliPHOSTenderSupply.h
1 #ifndef ALIPHOSTENDERSUPPLY_H
2 #define ALIPHOSTENDERSUPPLY_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 ////////////////////////////////////////////////////////////////////////
8 //                                                                    //
9 //  PHOS tender, apply corrections to PHOS clusters                   //
10 //  and do track matching                                             //
11 //  Author : Dmitri Peressounko (RRC KI)                              //
12 //                                                                    //
13 ////////////////////////////////////////////////////////////////////////
14
15 #include <AliTenderSupply.h>
16
17 class TVector3;
18 class AliPHOSGeometry;
19 class AliPHOSCalibData ;
20 class TH2I ;
21 class AliVCluster ;
22 class AliVCaloCells ;
23 class AliAnalysisTaskSE ;
24
25 class AliPHOSTenderSupply: public AliTenderSupply {
26   
27 public:
28   AliPHOSTenderSupply();
29   AliPHOSTenderSupply(const char *name, const AliTender *tender=NULL);
30   virtual ~AliPHOSTenderSupply();
31
32   virtual void   Init(){}
33   virtual void   ProcessEvent();
34   
35   void SetTask(AliAnalysisTaskSE * task){fTask=task;} //if work with AOD and special task
36   
37   void  SetNonlinearityVersion(const char * ver="Gustavo2005"){fNonlinearityVersion=ver;}
38   void  SetNonlinearityParams(Int_t n, const Double_t * par){
39             if(n>10){printf("Only 10 parameters allowed \n"); return ;}
40             for(Int_t i=0;i<n;i++)fNonlinearityParams[i]=par[i]; }
41   void  SetReconstructionPass(Int_t ipass=2){fRecoPass=ipass;}
42   
43   //Use this function to let tender know that it works with MC production
44   //should read calibration from PHOSMCCalibration file and use object for specified production
45   //By defaul real data is assumed.
46   void SetMCProduction(const char * name ="LHC13_b2"){fIsMC=kTRUE ; fMCProduction=name ;}
47   
48   //If you want to override automatic choise of bad maps and calibration
49   void ForceUsingBadMap(const char * filename="alien:///alice/cern.ch/user/p/prsnko/BadMaps/BadMap_LHC10b.root") ;
50   void ForceUsingCalibration(const char * filename="alien:///alice/cern.ch/user/p/prsnko/Recalibrations/LHC10b_pass1.root") ;
51   void   InitTender();
52
53 protected:
54   AliPHOSTenderSupply(const AliPHOSTenderSupply&c);
55   AliPHOSTenderSupply& operator= (const AliPHOSTenderSupply&c);
56   void   FindTrackMatching(Int_t mod,TVector3 *locpos,Double_t &dx, Double_t &dz, Double_t &pttrack, Int_t &charge); 
57   Float_t CorrectNonlinearity(Float_t en) ;
58   Double_t TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge) ;
59   Double_t TestLambda(Double_t pt,Double_t l1,Double_t l2) ;
60   Bool_t IsGoodChannel(Int_t mod, Int_t ix, Int_t iz) ;
61   void   CorrectPHOSMisalignment(TVector3 & globalPos, Int_t module);
62   void   EvalLambdas(AliVCluster * clu, Double_t &m02, Double_t &m20) ;
63   Double_t EvalTOF(AliVCluster * clu,AliVCaloCells * cells); 
64   Double_t CalibrateTOF(Double_t tof, Int_t absId, Bool_t isHG); 
65   void DistanceToBadChannel(Int_t mod, TVector3 * locPos, Double_t &minDist) ;
66
67  
68 private:
69
70   TString fOCDBpass ;                        // Pass to OCDB recalibration object, local or alien
71   TString fNonlinearityVersion;              // Version of non-linearity correction to aaply
72   AliPHOSGeometry   *fPHOSGeo;               // PHOS geometry
73   Double_t fNonlinearityParams[10] ;         // Parameters for non-linearity calculation
74   TH2I * fPHOSBadMap[5] ;                    // Bad channels map
75   Int_t fRecoPass ;                          // Reconstruction pass
76   Bool_t fUsePrivateBadMap ;
77   Bool_t fUsePrivateCalib ;
78   
79   AliPHOSCalibData *fPHOSCalibData;          // PHOS calibration object
80   AliAnalysisTaskSE     *fTask;              // analysis task
81
82   Bool_t fIsMC;                              //True if work with MC data
83   TString fMCProduction ;                    //Name of MC production
84  
85   ClassDef(AliPHOSTenderSupply, 2); // PHOS tender task
86 };
87
88
89 #endif
90
91