]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAODConversionPhoton.cxx
- added new Task for omega analysis using conversion in the 3 pion channel + dependencies
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAODConversionPhoton.cxx
1 #include "AliAODMCParticle.h"
2 #include "AliAODMCHeader.h"
3 #include "AliAODConversionPhoton.h"
4 #include "AliKFConversionPhoton.h"
5
6 using namespace std;
7
8 ClassImp(AliAODConversionPhoton)
9
10 AliAODConversionPhoton::AliAODConversionPhoton() :
11 AliAODConversionParticle(),
12 AliConversionPhotonBase(),
13 fDCArPrimVtx(0),
14 fDCAzPrimVtx(0),
15 fCaloPhoton(0),
16 fNCaloPhotonMCLabels(0),
17 fCaloPhotonMCFlags(0)
18 {
19         // initialize calo photon MC labels
20         for (Int_t i =0; i<50; i++){
21                 fCaloPhotonMCLabels[i]=-1;              
22         }
23   //Standard constructor
24 }
25
26 AliAODConversionPhoton::AliAODConversionPhoton(AliKFConversionPhoton *kfphoton) :
27 AliAODConversionParticle(kfphoton),
28 AliConversionPhotonBase(*((AliConversionPhotonBase*)kfphoton)),
29 fDCArPrimVtx(0),
30 fDCAzPrimVtx(0),
31 fCaloPhoton(0),
32 fNCaloPhotonMCLabels(0),
33 fCaloPhotonMCFlags(0)
34 {
35     //Constructor from kfphoton
36
37     // puts the mass to zero and store dilepton mass
38             SetMass(kfphoton->M());
39
40         //SetE(P());
41                 
42         // initialize calo photon MC labels
43         for (Int_t i =0; i<50; i++){
44                 fCaloPhotonMCLabels[i]=-1;              
45         }
46         
47 }
48
49 AliAODConversionPhoton::AliAODConversionPhoton(TLorentzVector *vec) :
50 AliAODConversionParticle(vec),
51 AliConversionPhotonBase(),
52 fDCArPrimVtx(0),
53 fDCAzPrimVtx(0),
54 fCaloPhoton(0),
55 fNCaloPhotonMCLabels(0),
56 fCaloPhotonMCFlags(0)
57 {
58     //Constructor from TLorentzVector
59
60         // initialize calo photon MC labels
61         for (Int_t i =0; i<50; i++){
62                 fCaloPhotonMCLabels[i]=-1;              
63         }
64 }
65
66
67
68 AliAODConversionPhoton::AliAODConversionPhoton(const AliAODConversionPhoton & original) :
69 AliAODConversionParticle(original),
70 AliConversionPhotonBase(original),
71 fDCArPrimVtx(original.fDCArPrimVtx),
72 fDCAzPrimVtx(original.fDCAzPrimVtx),
73 fCaloPhoton(original.fCaloPhoton),
74 fNCaloPhotonMCLabels(original.fNCaloPhotonMCLabels),
75 fCaloPhotonMCFlags(original.fCaloPhotonMCFlags)
76 {
77         //Copy constructor
78         
79         // initialize calo photon MC labels
80         for (Int_t i =0; i<50; i++){
81                 fCaloPhotonMCLabels[i]=original.fCaloPhotonMCLabels[i];         
82         }
83 }
84
85 AliAODConversionPhoton::~AliAODConversionPhoton()
86 {
87   // empty standard destructor
88 }
89
90 AliAODConversionPhoton & AliAODConversionPhoton::operator = (const AliAODConversionPhoton & /*source*/)
91 {
92   // assignment operator
93   return *this;
94 }
95
96 ///________________________________________________________________________
97 void AliAODConversionPhoton::CalculateDistanceOfClossetApproachToPrimVtx(const AliVVertex* primVertex ){
98
99    Double_t primCo[3] = {primVertex->GetX(),primVertex->GetY(),primVertex->GetZ()};
100
101    Double_t absoluteP = TMath::Sqrt(TMath::Power(GetPx(),2) + TMath::Power(GetPy(),2) + TMath::Power(GetPz(),2));   
102    Double_t p[3] = {GetPx()/absoluteP,GetPy()/absoluteP,GetPz()/absoluteP};
103    Double_t CP[3];
104    
105    CP[0] =  fConversionPoint[0] - primCo[0];
106    CP[1] =  fConversionPoint[1] - primCo[1];
107    CP[2] =  fConversionPoint[2] - primCo[2];
108    
109    Double_t Lambda = - (CP[0]*p[0]+CP[1]*p[1]+CP[2]*p[2])/(p[0]*p[0]+p[1]*p[1]+p[2]*p[2]);
110    
111    Double_t S[3];
112    S[0] = fConversionPoint[0] + p[0]*Lambda;
113    S[1] = fConversionPoint[1] + p[1]*Lambda;
114    S[2] = fConversionPoint[2] + p[2]*Lambda;
115   
116    fDCArPrimVtx = TMath::Sqrt( TMath::Power(primCo[0]-S[0],2) + TMath::Power(primCo[1]-S[1],2));
117    fDCAzPrimVtx = primCo[2]-S[2];
118     
119    return;
120 }
121
122
123 void AliAODConversionPhoton::SetCaloPhotonMCFlags(AliStack *MCStack){
124         
125         Bool_t isPhoton = kFALSE;                                               // largest contribution to cluster is photon 
126         Bool_t isElectron = kFALSE;                                             // largest contribution to cluster is electron
127         Bool_t isConversion = kFALSE;                                   // largest contribution to cluster is converted electron
128         Bool_t isConversionFullyContained = kFALSE;             // largest contribution to cluster is converted electron, second electron has been found in same cluster
129         Bool_t isMerged = kFALSE;                                               // largest contribution to cluster is photon, second photon or electron from dalitz has been found in same cluster 
130         Bool_t isMergedPartConv = kFALSE;                               // cluster contains more than one particle from the same decay and at least one of the particles came from a conversion
131         Bool_t isDalitz = kFALSE;                                               // this cluster was created by a particle stemming from a dality decay
132         Bool_t isDalitzMerged = kFALSE;                                 // this cluster was created by a particle stemming from a dality decay and more than one particle of the dalitz decay is contained in the cluster
133         Bool_t isPhotonWithElecMother = kFALSE;                 // this cluster is from a photon with an electron as mother
134         Bool_t isShower = kFALSE;                                               // this cluster contains as a largest contribution a particle from a shower or radiative process
135         
136         
137         TParticle* Photon;
138         if (fNCaloPhotonMCLabels==0) return;
139         Photon = MCStack->Particle(GetCaloPhotonMCLabel(0));
140         
141         if(Photon == NULL){
142                 return;
143         }
144
145
146         Int_t particleMotherLabel = Photon->GetMother(0);
147         Int_t particleGrandMotherLabel = -1; 
148         Int_t particleMotherPDG = -1; 
149         Int_t particleGrandMotherPDG = -1; 
150         Int_t particleMotherNDaugthers = 0;
151         Int_t particleGrandMotherNDaugthers = 0;
152         if (particleMotherLabel > -1){
153                 particleMotherNDaugthers = MCStack->Particle(Photon->GetMother(0))->GetNDaughters();
154                 particleGrandMotherLabel = MCStack->Particle(Photon->GetMother(0))->GetMother(0);
155                 particleMotherPDG = MCStack->Particle(Photon->GetMother(0))->GetPdgCode();
156                 if (particleGrandMotherLabel > -1){
157                         particleGrandMotherPDG = MCStack->Particle(MCStack->Particle(Photon->GetMother(0))->GetMother(0))->GetPdgCode();
158                         particleGrandMotherNDaugthers = MCStack->Particle(MCStack->Particle(Photon->GetMother(0))->GetMother(0))->GetNDaughters();
159                 }       
160         }
161                 
162         // Check wether the first contribution was photon       
163         if(MCStack->Particle(GetCaloPhotonMCLabel(0))->GetPdgCode() == 22){
164                 isPhoton=kTRUE;
165                 // did it decay via the dalitz channel
166                 if (particleMotherLabel > -1 && particleMotherNDaugthers == 3) isDalitz = kTRUE;
167                 // Test wether particle stems from a shower or radiation
168                 if (abs(particleMotherPDG) == 11){                                                                      // check whether photon stems from electron
169                         isPhotonWithElecMother = kTRUE;
170                         if (particleGrandMotherLabel > -1){                                                                     // test whether first particle has a grandmother
171                                 if (abs(particleGrandMotherPDG) == 22 ) isShower = kTRUE;       // check whether grandmother is a photon (meaning this is most likely a shower)
172                         }       
173                 }                       
174         }
175         // Check wether the first contribution was electron
176         if( abs(MCStack->Particle(GetCaloPhotonMCLabel(0))->GetPdgCode()) == 11 ){
177                 isElectron=kTRUE;       
178                 if (particleMotherLabel > -1) {
179                         // was it a conversion
180                         if (abs(particleMotherPDG) == 22) isConversion = kTRUE;
181                         // did it decay via the dalitz channel
182                         if (particleGrandMotherLabel > -1 && particleGrandMotherNDaugthers == 3 ) isDalitz = kTRUE;
183                 }
184                 if (particleGrandMotherLabel > -1){                                                                             // check whether electron has a grandmother
185                         if (abs(particleGrandMotherPDG) == 11 ||  abs(particleGrandMotherPDG) == 22){   // test wether electron has photon or electron as grandmother (meaning will most likely be a shower)
186                                 isShower = kTRUE; 
187                         }       
188                 }
189         }
190         
191         // check wether there were other contributions to the cluster
192         if (fNCaloPhotonMCLabels>1){
193                 // largest contribution was from photon and is not from shower or electron mother
194                 if (isPhoton && (!isShower || !isPhotonWithElecMother )){
195 //                      cout << "largest contribution is photon, with mother: " <<  particleMotherLabel <<"(" << particleMotherPDG<< ")" << "\t with " << particleMotherNDaugthers << " daugthers and grand mother: " << particleGrandMotherLabel  <<"(" << particleGrandMotherPDG << ")" << endl;
196                         
197                         TParticle* dummyPart =NULL;
198                         for (Int_t i = 1; i< fNCaloPhotonMCLabels; i++){
199                                 if (i > 49) continue;                                                                                                   // abort if more than 20 entries to the cluster have been checked (more are not stored in these objects)
200                                 dummyPart = MCStack->Particle(GetCaloPhotonMCLabel(i));
201                                 Int_t dummyPartMotherLabel = dummyPart->GetMother(0);
202                                 Int_t dummyPartGrandMotherLabel = -1;
203                                 Int_t dummyPartMotherPDG = -1; 
204                                 Int_t dummyPartGrandMotherPDG = -1; 
205                                         
206                                 // check wether this particle has a mother & obtain the pdg code
207                                 if (dummyPartMotherLabel > -1){
208                                         dummyPartGrandMotherLabel = MCStack->Particle(dummyPart->GetMother(0))->GetMother(0);
209                                         dummyPartMotherPDG = MCStack->Particle(dummyPart->GetMother(0))->GetPdgCode(); 
210                                         // check wether this particle has a grandmother & obtain its pdg code
211                                         if (dummyPartGrandMotherLabel > -1){
212                                                 dummyPartGrandMotherPDG = MCStack->Particle(MCStack->Particle(dummyPart->GetMother(0))->GetMother(0))->GetPdgCode();
213                                         }
214                                 }
215                                 
216                                 if (particleMotherLabel > -1){                                                                                                  // test whether first particle has a mother
217                                         if (dummyPartMotherLabel == particleMotherLabel) isMerged = kTRUE;                      // test whether current and first particle have the same mother => i.e. other gamma from decay or dalitz electron
218                                         if (dummyPartGrandMotherLabel > -1){                                                                            // test whether first particle has a grandmother
219                                                 // check wether particle is an electron from a conversion of a photon from the original mother
220                                                 if (abs(dummyPart->GetPdgCode()) == 11 && dummyPartGrandMotherLabel == particleMotherLabel ) isMergedPartConv = kTRUE;  
221                                                 // check wether particle is an electron from a dalitz decay from the original mother
222                                                 if (abs(dummyPart->GetPdgCode()) == 11 && dummyPartMotherLabel == particleMotherLabel ) isDalitzMerged = kTRUE;
223                                         }       
224                                 }       
225 //                              cout << i << "\t" <<GetCaloPhotonMCLabel(i) << "\t" << dummyPart->GetPdgCode() << "\t mother is: " <<  dummyPartMotherLabel <<"(" << dummyPartMotherPDG<< ")"<< "\t grandmother is: " <<  dummyPartGrandMotherLabel<<"(" << dummyPartGrandMotherPDG<< ")"<<  endl;
226                                 
227                         }       
228                 }
229                 // largest contribution was from electron & not a from a shower
230                 if (isElectron && !isShower){
231 //                      cout << "largest contribution is electron, with mother: " <<  particleMotherLabel <<"(" << particleMotherPDG<< ")" << "\t and grand mother: " << particleGrandMotherLabel  <<"(" << particleGrandMotherPDG << ")" << endl;
232                         TParticle* dummyPart =NULL;
233                         for (Int_t i = 1; i< fNCaloPhotonMCLabels; i++){
234                                 if (i > 49) continue;
235                                 
236                                 dummyPart = MCStack->Particle(GetCaloPhotonMCLabel(i));
237                                 Int_t dummyPartMotherLabel = dummyPart->GetMother(0);
238                                 Int_t dummyPartGrandMotherLabel = -1;
239                                 Int_t dummyPartMotherPDG = -1; 
240                                 Int_t dummyPartGrandMotherPDG = -1; 
241                                 
242                                 // check wether this particle has a mother & obtain the pdg code
243                                 if (dummyPartMotherLabel > -1){
244                                         dummyPartGrandMotherLabel = MCStack->Particle(dummyPart->GetMother(0))->GetMother(0);
245                                         dummyPartMotherPDG = MCStack->Particle(dummyPart->GetMother(0))->GetPdgCode();
246                                         // check wether this particle has a grandmother & obtain its pdg code
247                                         if (dummyPartGrandMotherLabel > -1){
248                                                 dummyPartGrandMotherPDG = MCStack->Particle(MCStack->Particle(dummyPart->GetMother(0))->GetMother(0))->GetPdgCode();
249                                         }
250                                 }
251                                 
252                                 if (particleMotherLabel > -1){                                                                                                                                                  // test whether first particle has a mother
253                                         if (isConversion && dummyPartMotherLabel == particleMotherLabel) isConversionFullyContained = kTRUE;            // test whether conversion is fully contained in cluster 
254
255                                         if (particleGrandMotherLabel > -1){                                                                                                                                     // test whether first particle has a grandmother                                                                                                                        
256                                                 if (abs(dummyPart->GetPdgCode()) == 22){                                                                                                        // test whether this particle is a photon
257                                                         // check wether orginal electron and this photon stem from the same particle and electron stems from conversion
258                                                         if( dummyPartMotherLabel == particleGrandMotherLabel && (abs(dummyPartMotherPDG) != 11 ||  abs(dummyPartMotherPDG) != 22 )   ) isMergedPartConv = kTRUE;
259                                                         // check wether orginal electron and this photon stem from the same particle and electron originated in dalitz
260                                                         if( dummyPartMotherLabel == particleMotherLabel && (abs(dummyPartMotherPDG) != 11 ||  abs(dummyPartMotherPDG) != 22 )   ) isDalitzMerged = kTRUE;
261                                                 }
262                                                 if (abs(dummyPart->GetPdgCode()) == 11) {
263                                                         // check wether orginal electron and this electron stem from the same particle and electron stems from conversion
264                                                         if( dummyPartGrandMotherLabel == particleGrandMotherLabel &&  (abs(dummyPartGrandMotherPDG) != 11 ||  abs(dummyPartGrandMotherPDG) != 22 )   ) isMergedPartConv = kTRUE;
265                                                         // check wether orginal electron and this electron stem from the same particle and electron originated in dalitz decay
266                                                         if( dummyPartMotherLabel == particleMotherLabel && abs(particleMotherPDG) != 22 &&  (abs(dummyPartGrandMotherPDG) != 11 ||  abs(dummyPartGrandMotherPDG) != 22 )   ) isDalitzMerged = kTRUE;
267                                                 }
268                                                 
269                                         }       
270                                 }       
271 //                              cout << i << "\t" << GetCaloPhotonMCLabel(i) << "\t" << dummyPart->GetPdgCode() << "\t mother is: " <<  dummyPartMotherLabel <<"(" << dummyPartMotherPDG<< ")"<< "\t grandmother is: " <<  dummyPartGrandMotherLabel<<"(" << dummyPartGrandMotherPDG<< ")"<<  endl;
272                                 
273                         }       
274                 }
275         }               
276         fCaloPhotonMCFlags = isPhoton *1 + isElectron *2 + isConversion*4+ isConversionFullyContained *8 + isMerged *16 + isMergedPartConv*32 + isDalitz *64 + isDalitzMerged *128 + isPhotonWithElecMother *256 + isShower * 512;              
277 //      cout << "isPhoton: \t" << isPhoton << "\t isElectron: \t" << isElectron << "\t isConversion: \t" << isConversion <<  "\t isConversionFullyContained: \t" << isConversionFullyContained << "\t isMerged: \t" << isMerged << "\t isMergedPartConv: \t" << isMergedPartConv << "\t isPhotonWithElecMother: \t" << isPhotonWithElecMother<< "\t isDalitz: \t" << isDalitz << "\t isDalitzMerged: \t" << isDalitzMerged << "\t isShower: \t" << isShower << "\t"<< fCaloPhotonMCFlags<< endl;                
278         
279 }
280
281 void AliAODConversionPhoton::SetCaloPhotonMCFlagsAOD(AliVEvent* event){
282         
283         TClonesArray *AODMCTrackArray = dynamic_cast<TClonesArray*>(event->FindListObject(AliAODMCParticle::StdBranchName()));
284         if (!AODMCTrackArray) return;
285         
286         
287         Bool_t isPhoton = kFALSE;                                               // largest contribution to cluster is photon 
288         Bool_t isElectron = kFALSE;                                             // largest contribution to cluster is electron
289         Bool_t isConversion = kFALSE;                                   // largest contribution to cluster is converted electron
290         Bool_t isConversionFullyContained = kFALSE;             // largest contribution to cluster is converted electron, second electron has been found in same cluster
291         Bool_t isMerged = kFALSE;                                               // largest contribution to cluster is photon, second photon or electron from dalitz has been found in same cluster 
292         Bool_t isMergedPartConv = kFALSE;                               // cluster contains more than one particle from the same decay and at least one of the particles came from a conversion
293         Bool_t isDalitz = kFALSE;                                               // this cluster was created by a particle stemming from a dality decay
294         Bool_t isDalitzMerged = kFALSE;                                 // this cluster was created by a particle stemming from a dality decay and more than one particle of the dalitz decay is contained in the cluster
295         Bool_t isPhotonWithElecMother = kFALSE;                 // this cluster is from a photon with an electron as mother
296         Bool_t isShower = kFALSE;                                               // this cluster contains as a largest contribution a particle from a shower or radiative process
297         
298         
299         AliAODMCParticle* Photon;
300         if (fNCaloPhotonMCLabels==0) return;
301         Photon = (AliAODMCParticle*) AODMCTrackArray->At(GetCaloPhotonMCLabel(0));
302         
303         if(Photon == NULL){
304                 return;
305         }
306
307         AliAODMCParticle* PhotonMother;
308         AliAODMCParticle* PhotonGrandMother;
309         Int_t particleMotherLabel = Photon->GetMother();
310         Int_t particleGrandMotherLabel = -1; 
311         Int_t particleMotherPDG = -1; 
312         Int_t particleGrandMotherPDG = -1; 
313         Int_t particleMotherNDaugthers = 0;
314         Int_t particleGrandMotherNDaugthers = 0;
315         if (particleMotherLabel > -1){
316                 PhotonMother = (AliAODMCParticle*) AODMCTrackArray->At(Photon->GetMother());
317                 particleMotherNDaugthers = PhotonMother->GetNDaughters();
318                 particleGrandMotherLabel = PhotonMother->GetMother();
319                 particleMotherPDG = PhotonMother->GetPdgCode();
320                 if (particleGrandMotherLabel > -1){
321                         PhotonGrandMother = (AliAODMCParticle*) AODMCTrackArray->At(PhotonMother->GetMother());
322                         particleGrandMotherPDG = PhotonGrandMother->GetPdgCode();
323                         particleGrandMotherNDaugthers = PhotonGrandMother->GetNDaughters();
324                 }       
325         }
326                 
327         // Check wether the first contribution was photon       
328         if(abs(Photon->GetPdgCode()) == 22){
329                 isPhoton=kTRUE;
330                 // did it decay via the dalitz channel
331                 if (particleMotherLabel > -1 && particleMotherNDaugthers == 3) isDalitz = kTRUE;
332                 // Test wether particle stems from a shower or radiation
333                 if (abs(particleMotherPDG) == 11){                                                                      // check whether photon stems from electron
334                         isPhotonWithElecMother = kTRUE;
335                         if (particleGrandMotherLabel > -1){                                                                     // test whether first particle has a grandmother
336                                 if (abs(particleGrandMotherPDG) == 22 ) isShower = kTRUE;       // check whether grandmother is a photon (meaning this is most likely a shower)
337                         }       
338                 }                       
339         }
340         // Check wether the first contribution was electron
341         if(abs(Photon->GetPdgCode()) == 11 ){
342                 isElectron=kTRUE;       
343                 if (particleMotherLabel > -1) {
344                         // was it a conversion
345                         if (abs(particleMotherPDG) == 22) isConversion = kTRUE;
346                         // did it decay via the dalitz channel
347                         if (particleGrandMotherLabel > -1 && particleGrandMotherNDaugthers == 3 ) isDalitz = kTRUE;
348                 }
349                 if (particleGrandMotherLabel > -1){                                                                             // check whether electron has a grandmother
350                         if (abs(particleGrandMotherPDG) == 11 ||  abs(particleGrandMotherPDG) == 22){   // test wether electron has photon or electron as grandmother (meaning will most likely be a shower)
351                                 isShower = kTRUE; 
352                         }       
353                 }
354         }
355         
356         // check wether there were other contributions to the cluster
357         if (fNCaloPhotonMCLabels>1){
358                 // largest contribution was from photon and is not from shower or electron mother
359                 if (isPhoton && (!isShower || !isPhotonWithElecMother )){
360 //                      cout << "largest contribution is photon, with mother: " <<  particleMotherLabel <<"(" << particleMotherPDG<< ")" << "\t with " << particleMotherNDaugthers << " daugthers and grand mother: " << particleGrandMotherLabel  <<"(" << particleGrandMotherPDG << ")" << endl;
361                         
362                         AliAODMCParticle* dummyPart =NULL;
363                         AliAODMCParticle* dummyPartMother = NULL;
364                         AliAODMCParticle* dummyPartGrandMother = NULL;
365
366                         for (Int_t i = 1; i< fNCaloPhotonMCLabels; i++){
367                                 if (i > 49) continue;                                                                                                   // abort if more than 20 entries to the cluster have been checked (more are not stored in these objects)
368                                 dummyPart = (AliAODMCParticle*) AODMCTrackArray->At(GetCaloPhotonMCLabel(i));
369                                 Int_t dummyPartMotherLabel = dummyPart->GetMother();
370                                 Int_t dummyPartGrandMotherLabel = -1;
371                                 Int_t dummyPartMotherPDG = -1; 
372                                 Int_t dummyPartGrandMotherPDG = -1; 
373                                         
374                                 // check wether this particle has a mother & obtain the pdg code
375                                 if (dummyPartMotherLabel > -1){
376                                         dummyPartMother = (AliAODMCParticle*) AODMCTrackArray->At(dummyPart->GetMother());
377                                         dummyPartGrandMotherLabel = dummyPartMother->GetMother();
378                                         dummyPartMotherPDG = dummyPartMother->GetPdgCode(); 
379                                         // check wether this particle has a grandmother & obtain its pdg code
380                                         if (dummyPartGrandMotherLabel > -1){
381                                                 dummyPartGrandMother = (AliAODMCParticle*) AODMCTrackArray->At(dummyPartMother->GetMother());
382                                                 dummyPartGrandMotherPDG = dummyPartGrandMother->GetPdgCode();
383                                         }
384                                 }
385                                 
386                                 if (particleMotherLabel > -1){                                                                                                  // test whether first particle has a mother
387                                         if (dummyPartMotherLabel == particleMotherLabel) isMerged = kTRUE;                      // test whether current and first particle have the same mother => i.e. other gamma from decay or dalitz electron
388                                         if (dummyPartGrandMotherLabel > -1){                                                                            // test whether first particle has a grandmother
389                                                 // check wether particle is an electron from a conversion of a photon from the original mother
390                                                 if (abs(dummyPart->GetPdgCode()) == 11 && dummyPartGrandMotherLabel == particleMotherLabel ) isMergedPartConv = kTRUE;  
391                                                 // check wether particle is an electron from a dalitz decay from the original mother
392                                                 if (abs(dummyPart->GetPdgCode()) == 11 && dummyPartMotherLabel == particleMotherLabel ) isDalitzMerged = kTRUE;
393                                         }       
394                                 }       
395 //                              cout << i << "\t" <<GetCaloPhotonMCLabel(i) << "\t" << dummyPart->GetPdgCode() << "\t mother is: " <<  dummyPartMotherLabel <<"(" << dummyPartMotherPDG<< ")"<< "\t grandmother is: " <<  dummyPartGrandMotherLabel<<"(" << dummyPartGrandMotherPDG<< ")"<<  endl;
396                                 
397                         }       
398                 }
399                 // largest contribution was from electron & not a from a shower
400                 if (isElectron && !isShower){
401 //                      cout << "largest contribution is electron, with mother: " <<  particleMotherLabel <<"(" << particleMotherPDG<< ")" << "\t and grand mother: " << particleGrandMotherLabel  <<"(" << particleGrandMotherPDG << ")" << endl;
402                         AliAODMCParticle* dummyPart =NULL;
403                         AliAODMCParticle* dummyPartMother = NULL;
404                         AliAODMCParticle* dummyPartGrandMother = NULL;
405
406                         for (Int_t i = 1; i< fNCaloPhotonMCLabels; i++){
407                                 if (i > 49) continue;
408                                 
409                                 dummyPart = (AliAODMCParticle*)AODMCTrackArray->At(GetCaloPhotonMCLabel(i));
410                                 Int_t dummyPartMotherLabel = dummyPart->GetMother();
411                                 Int_t dummyPartGrandMotherLabel = -1;
412                                 Int_t dummyPartMotherPDG = -1; 
413                                 Int_t dummyPartGrandMotherPDG = -1; 
414                                 
415                                 if (dummyPartMotherLabel > -1){
416                                         dummyPartMother = (AliAODMCParticle*) AODMCTrackArray->At(dummyPart->GetMother());
417                                         dummyPartGrandMotherLabel = dummyPartMother->GetMother();
418                                         dummyPartMotherPDG = dummyPartMother->GetPdgCode(); 
419                                         // check wether this particle has a grandmother & obtain its pdg code
420                                         if (dummyPartGrandMotherLabel > -1){
421                                                 dummyPartGrandMother = (AliAODMCParticle*) AODMCTrackArray->At(dummyPartMother->GetMother());
422                                                 dummyPartGrandMotherPDG = dummyPartGrandMother->GetPdgCode();
423                                         }
424                                 }
425                                 
426                                 if (particleMotherLabel > -1){                                                                                                                                                  // test whether first particle has a mother
427                                         if (isConversion && dummyPartMotherLabel == particleMotherLabel) isConversionFullyContained = kTRUE;            // test whether conversion is fully contained in cluster 
428
429                                         if (particleGrandMotherLabel > -1){                                                                                                                                     // test whether first particle has a grandmother                                                                                                                        
430                                                 if (abs(dummyPart->GetPdgCode()) == 22){                                                                                                        // test whether this particle is a photon
431                                                         // check wether orginal electron and this photon stem from the same particle and electron stems from conversion
432                                                         if( dummyPartMotherLabel == particleGrandMotherLabel && (abs(dummyPartMotherPDG) != 11 ||  abs(dummyPartMotherPDG) != 22 )   ) isMergedPartConv = kTRUE;
433                                                         // check wether orginal electron and this photon stem from the same particle and electron originated in dalitz
434                                                         if( dummyPartMotherLabel == particleMotherLabel && (abs(dummyPartMotherPDG) != 11 ||  abs(dummyPartMotherPDG) != 22 )   ) isDalitzMerged = kTRUE;
435                                                 }
436                                                 if (abs(dummyPart->GetPdgCode()) == 11) {
437                                                         // check wether orginal electron and this electron stem from the same particle and electron stems from conversion
438                                                         if( dummyPartGrandMotherLabel == particleGrandMotherLabel &&  (abs(dummyPartGrandMotherPDG) != 11 ||  abs(dummyPartGrandMotherPDG) != 22 )   ) isMergedPartConv = kTRUE;
439                                                         // check wether orginal electron and this electron stem from the same particle and electron originated in dalitz decay
440                                                         if( dummyPartMotherLabel == particleMotherLabel && abs(particleMotherPDG) != 22 &&  (abs(dummyPartGrandMotherPDG) != 11 ||  abs(dummyPartGrandMotherPDG) != 22 )   ) isDalitzMerged = kTRUE;
441                                                 }
442                                                 
443                                         }       
444                                 }       
445 //                              cout << i << "\t" << GetCaloPhotonMCLabel(i) << "\t" << dummyPart->GetPdgCode() << "\t mother is: " <<  dummyPartMotherLabel <<"(" << dummyPartMotherPDG<< ")"<< "\t grandmother is: " <<  dummyPartGrandMotherLabel<<"(" << dummyPartGrandMotherPDG<< ")"<<  endl;
446                                 
447                         }       
448                 }
449         }               
450
451         fCaloPhotonMCFlags = isPhoton *1 + isElectron *2 + isConversion*4+ isConversionFullyContained *8 + isMerged *16 + isMergedPartConv*32 + isDalitz *64 + isDalitzMerged *128 + isPhotonWithElecMother *256 + isShower * 512;              
452 //      cout << "isPhoton: \t" << isPhoton << "\t isElectron: \t" << isElectron << "\t isConversion: \t" << isConversion <<  "\t isConversionFullyContained: \t" << isConversionFullyContained << "\t isMerged: \t" << isMerged << "\t isMergedPartConv: \t" << isMergedPartConv << "\t isPhotonWithElecMother: \t" << isPhotonWithElecMother<< "\t isDalitz: \t" << isDalitz << "\t isDalitzMerged: \t" << isDalitzMerged << "\t isShower: \t" << isShower << "\t"<< fCaloPhotonMCFlags<< endl;                
453         
454 }