]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGGA/PHOSTasks/PHOS_embedding/AliPHOSEmbedding.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGGA / PHOSTasks / PHOS_embedding / AliPHOSEmbedding.h
CommitLineData
a65a7e70 1#ifndef AliPHOSEmbedding_h
2#define AliPHOSEmbedding_h
3
4/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
6
7// Class to perform embedding on the AOD level
8// Author: D.Peressounko
9
10class TChain ;
11class TClonesArray ;
12class TH2F ;
13
14class AliPHOSClusterizerv1 ;
15class AliPHOSReconstructor ;
16class AliAODEvent ;
17class AliESDEvent ;
18class AliESDtrack ;
19class AliESDCaloCells ;
20
21#include "AliAnalysisTaskSE.h"
22
23class AliPHOSEmbedding : public AliAnalysisTaskSE {
24public:
25 AliPHOSEmbedding(const char *name = "AliPHOSEmbedding");
26 virtual ~AliPHOSEmbedding() {}
27
28 //Standard methods
29 virtual void UserCreateOutputObjects();
30 virtual void UserExec(Option_t *option);
31 virtual void Terminate(Option_t *){}
32
33 //Chain with signal AOD for embedding
34 void SetSignalChain(TChain * signal){fAODChain =signal;}
35
36 //Calibration used in reconstruction of real data (ESDs)
37 //If not set, assume same calibration as set by default
38 void SetOldCalibration(TH2F **calib) ;
39
40private:
41 AliPHOSEmbedding(const AliPHOSEmbedding&); // not implemented
42 AliPHOSEmbedding& operator=(const AliPHOSEmbedding&); // not implemented
43
44 void Init() ;
45 void InitMF() ; //Mag.Field initialization for track matching
46 void InitGeometry() ;
47
48 AliAODEvent * GetNextSignalEvent(void) ;
49
50 void CopyRecalibrateDigits(void) ;
51 void MakeEmbedding(AliESDEvent * data, AliAODEvent * signal) ;
52 void MakeDigits(AliAODEvent* signal) ;
53
54 //Convert ESD with embedded signal to AOD
55 //First standard stuff
56 void ConvertESDtoAOD(AliESDEvent *esd) ;
57 void ConvertHeader(AliESDEvent &esd) ;
58 void ConvertPrimaryVertices(const AliESDEvent &esd) ;
59 void ConvertCaloClusters(const AliESDEvent &esd) ;
60 void ConvertEMCALCells(const AliESDEvent &esd) ;
61 void ConvertPHOSCells(const AliESDEvent &esd) ;
62
63 //Add new branch
64 void ConvertEmbeddedClusters(const AliESDEvent *esd) ;
65 void ConvertEmbeddedCells(const AliESDEvent *esd) ;
66 void ConvertMCParticles(const AliAODEvent *aod) ;
67
68 Float_t TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge) ;
69
70
71 TChain * fAODChain ; //Signal
72
73 TTree * fDigitsTree ; //! Digits
74 TTree * fClustersTree; //! Clusters
75 TTree * fTreeOut; //Output AOD
76 TClonesArray * fDigitsArr ; //!
77
78 TClonesArray * fEmbeddedClusters ; //!
79 AliAODCaloCells * fEmbeddedCells ; //!
80 AliESDCaloCells * fCellsPHOS ; //! Old PHOS cells
81
82 AliPHOSClusterizerv1 * fClusterizer ; //!
83 AliPHOSReconstructor * fPHOSReconstructor ; //!
84
85 TH2F * fOldPHOSCalibration[5] ; //! Calibration coeff. used in ESD production
86
87 Int_t fNSignal ; // Number of signal evetns processed
88 Int_t fNCaloClustersOld ; //Number of CaloClusters already in ESD
89 Bool_t fInitialized ; //!
90 ClassDef(AliPHOSEmbedding, 1); // PHOS analysis task
91};
92
93#endif