]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEER/AliGlobalQADataMaker.cxx
First version of HLT component using AliAnalysisManager
[u/mrichter/AliRoot.git] / STEER / STEER / AliGlobalQADataMaker.cxx
1 /*
2  The class for calculating the global (not detector specific) quality assurance.
3  It reuses the following TLists from its base class 
4     AliQADataMaker::fRecPointsQAList (for keeping the track residuals)
5     AliQADataMaker::fESDsQAList      (for keeping global ESD QA data)
6 */
7
8 #include <TPDGCode.h>
9 #include <TH1F.h>
10 #include <TH2F.h>
11
12 #include "AliCDBPath.h"
13 #include "AliCDBEntry.h"
14 #include "AliCDBManager.h"
15 #include "AliDetectorRecoParam.h"
16 #include "AliQAChecker.h"
17 #include "AliGlobalQADataMaker.h"
18 #include "AliGeomManager.h"
19 #include "AliESDEvent.h"
20 #include "AliESDv0.h"
21 #include "AliRawReader.h"
22 #include "AliESDVZERO.h"
23 #include "AliMultiplicity.h" 
24 #include "AliPID.h"
25
26 ClassImp(AliGlobalQADataMaker)
27  
28 //____________________________________________________________________________ 
29 void AliGlobalQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
30 {
31   //Detector specific actions at end of cycle
32   // do the QA checking
33   AliQAChecker::Instance()->Run(AliQAv1::kGLOBAL, task, list) ;  
34 }
35
36 //____________________________________________________________________________ 
37 void AliGlobalQADataMaker::InitRaws()
38 {
39   // create Raws histograms in Raws subdir
40   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
41 }
42
43 //____________________________________________________________________________
44 void AliGlobalQADataMaker::InitRecoParams() 
45 {
46   // Get the recoparam form the OCDB 
47   if (!fRecoParam) {
48     TString name("GRP") ; 
49     AliDebug(AliQAv1::GetQADebugLevel(), Form("Loading reconstruction parameter objects for detector %s", name.Data()));
50     AliCDBPath path(name.Data(),"Calib","RecoParam");
51     AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
52     if(!entry) {
53       fRecoParam = NULL ; 
54       AliDebug(AliQAv1::GetQADebugLevel(), Form("Couldn't find RecoParam entry in OCDB for detector %s",name.Data()));
55     }
56     else {
57       TObject * recoParamObj = entry->GetObject() ; 
58       if ( strcmp(recoParamObj->ClassName(), "TObjArray") == 0 ) {
59           // The detector has only one set of reco parameters
60         AliDebug(AliQAv1::GetQADebugLevel(), Form("Array of reconstruction parameters found for detector %s",name.Data()));
61         TObjArray *recoParamArray = static_cast<TObjArray*>(recoParamObj) ;
62         for (Int_t iRP=0; iRP<recoParamArray->GetEntriesFast(); iRP++) {
63           fRecoParam = static_cast<AliDetectorRecoParam*>(recoParamArray->At(iRP)) ;
64           if (!fRecoParam) 
65             break ; 
66           else if (fRecoParam->IsDefault()) 
67             break ; 
68         }
69       }
70       else if (recoParamObj->InheritsFrom("AliDetectorRecoParam")) {
71           // The detector has only one set of reco parameters
72           // Registering it in AliRecoParam
73         AliDebug(AliQAv1::GetQADebugLevel(), Form("Single set of reconstruction parameters found for detector %s",name.Data()));
74         fRecoParam = static_cast<AliDetectorRecoParam*>(recoParamObj) ;
75         static_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
76       } else { 
77         AliError(Form("No valid RecoParam object found in the OCDB for detector %s",name.Data()));
78       }
79     }
80   }
81 }
82
83 //____________________________________________________________________________ 
84 void AliGlobalQADataMaker::InitRecPoints() {
85   //------------------------------------------------------
86   // This function books the histograms of *track*residuals*
87   // as a part of global QA
88   //------------------------------------------------------
89   static Bool_t first = kTRUE ; 
90   if ( ! first ) return; 
91   const Char_t *name[]={
92     "hGlobalSPD1ResidualsY","SPD1ResidualsZ",
93     "hGlobalSPD2ResidualsY","SPD2ResidualsZ",
94     "hGlobalSDD1ResidualsY","SDD1ResidualsZ",
95     "hGlobalSDD2ResidualsY","SDD2ResidualsZ",
96     "hGlobalSSD1ResidualsY","SSD1ResidualsZ",
97     "hGlobalSSD2ResidualsY","SSD2ResidualsZ",
98     
99     "hGlobalTPC1ResidualsY","TPC1ResidualsZ",
100     "hGlobalTPC2ResidualsY","TPC2ResidualsZ",
101     
102     "hGlobalTRD1ResidualsY","TRD1ResidualsZ",
103     "hGlobalTRD2ResidualsY","TRD2ResidualsZ",
104     "hGlobalTRD3ResidualsY","TRD3ResidualsZ",
105     "hGlobalTRD4ResidualsY","TRD4ResidualsZ",
106     "hGlobalTRD5ResidualsY","TRD5ResidualsZ",
107     "hGlobalTRD6ResidualsY","TRD6ResidualsZ",
108     
109     "hGlobalTOFResidualsY","TOFResidualsZ",
110     
111     "hGlobalPHOS1ResidualsY","PHOS1ResidualsZ",
112     "hGlobalPHOS2ResidualsY","PHOS2ResidualsZ",
113     
114     "hGlobalHMPIDResidualsY","HMPIDResidualsZ",
115     
116     "hGlobalMUONResidualsY","MUONResidualsZ",
117     
118     "hGlobalEMCALResidualsY","EMCALResidualsZ"
119   };
120   const Char_t *title[]={
121     "SPD1 residuals Y","SPD1 residuals Z",
122     "SPD2 residuals Y","SPD2 residuals Z",
123     "SDD1 residuals Y","SDD1 residuals Z",
124     "SDD2 residuals Y","SDD2 residuals Z",
125     "SSD1 residuals Y","SSD1 residuals Z",
126     "SSD2 residuals Y","SSD2 residuals Z",
127     
128     "TPC1 residuals Y","TPC1 residuals Z",
129     "TPC2 residuals Y","TPC2 residuals Z",
130     
131     "TRD1 residuals Y","TRD1 residuals Z",
132     "TRD2 residuals Y","TRD2 residuals Z",
133     "TRD3 residuals Y","TRD3 residuals Z",
134     "TRD4 residuals Y","TRD4 residuals Z",
135     "TRD5 residuals Y","TRD5 residuals Z",
136     "TRD6 residuals Y","TRD6 residuals Z",
137     
138     "TOF residuals Y","TOF residuals Z",
139     
140     "PHOS1 residuals Y","PHOS1 residuals Z",
141     "PHOS2 residuals Y","PHOS2 residuals Z",
142     
143     "HMPID residuals Y","HMPID residuals Z",
144     
145     "MUON residuals Y","MUON residuals Z",
146     
147     "EMCAL residuals Y","EMCAL residuals Z"
148   };
149   
150   for (Int_t m=1; m<AliGeomManager::kLastLayer; m++) {
151     Int_t i=2*m-2;
152     TH1F *h=new TH1F(name[i],title[i],100,-5.,5.);
153     Add2RecPointsList(h,i);    
154     h=new TH1F(name[i+1],title[i+1],100,-5.,5.);
155     Add2RecPointsList(h,i+1);    
156   }
157
158   Add2RecPointsList(
159   new TH1F("hGlobalSSD1AbsoluteResidualsYNegZ",
160            "SSD1 Absolute Residuals Y Neg Z",100,-2.,2.),40);
161   Add2RecPointsList(
162   new TH1F("hGlobalSSD1AbsoluteResidualsZNegZ",
163            "SSD1 Absolute Residuals Z Neg Z",100,-2.,2.),41);
164   Add2RecPointsList(
165   new TH1F("hGlobalSSD1AbsoluteResidualsYPosZ",
166            "SSD1 Absolute Residuals Y Pos Z",100,-2.,2.),42);
167   Add2RecPointsList(
168   new TH1F("hGlobalSSD1AbsoluteResidualsZPosZ",
169            "SSD1 Absolute Residuals Z Pos Z",100,-2.,2.),43);
170
171
172   Add2RecPointsList(
173   new TH1F("hGlobalSSD2AbsoluteResidualsYNegZ",
174            "SSD2 Absolute Residuals Y Neg Z",100,-3.,3.),44);
175   Add2RecPointsList(
176   new TH1F("hGlobalSSD2AbsoluteResidualsZNegZ",
177            "SSD2 Absolute Residuals Z Neg Z",100,-3.,3.),45);
178   Add2RecPointsList(
179   new TH1F("hGlobalSSD2AbsoluteResidualsYPosZ",
180            "SSD2 Absolute Residuals Y Pos Z",100,-3.,3.),46);
181   Add2RecPointsList(
182   new TH1F("hGlobalSSD2AbsoluteResidualsZPosZ",
183            "SSD2Absolute Residuals Z Pos Z",100,-3.,3.),47);
184   
185   first = kFALSE;
186   //
187   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line 
188 }
189
190 //____________________________________________________________________________ 
191 void AliGlobalQADataMaker::InitESDs() {
192   //------------------------------------------------------
193   // This function books the ESD QA histograms
194   // as a part of global QA
195   //------------------------------------------------------
196
197   const Bool_t expert   = kTRUE ; 
198   const Bool_t image    = kTRUE ;
199
200   {// Event related QA
201     const Char_t *name[]={
202       "hGlobalPrimaryVertex"
203     };
204     const Char_t *title[]={
205       "Z-distribution of the primary vertex"
206     };
207     Add2ESDsList(new TH1F(name[0],title[0],100,-20.,20.),kEvt0,!expert,image);
208   }
209  
210   {// Cluster related QA
211     const Char_t *name[]={
212       "hGlobalFractionAssignedClustersITS",
213       "hGlobalFractionAssignedClustersTPC",
214       "hGlobalFractionAssignedClustersTRD",
215       "hGlobalClustersPerITSModule"
216     };
217     const Char_t *title[]={
218       "Fraction of the assigned clusters in ITS",
219       "Fraction of the assigned clusters in TPC",
220       "Fraction of the assigned clusters in TRD",
221       "Number of clusters per an ITS module"
222     };
223     Add2ESDsList(new TH1F(name[0],title[0],100,0.,2.),kClr0, !expert, image);
224     Add2ESDsList(new TH1F(name[1],title[1],100,0.,2.),kClr1, !expert, image);
225     Add2ESDsList(new TH1F(name[2],title[2],100,0.,2.),kClr2, !expert, image);
226     Add2ESDsList(new TH1F(name[3],title[3],2201,-0.5,2200.5),kClr3, !expert, image);
227   }
228
229   {// Track related QA
230     const Char_t *name[]={
231       "hGlobalTrackAzimuthe",                               // kTrk0
232       "hGlobalTrackEta",                                    // kTrk1
233       "hGlobalTPCTrackpT",                                  // kTrk2
234       "hGlobalTPCITSMatchedpT",                             // kTrk3
235       "hGlobalTPCTOFMatchedpT",                             // kTrk4
236       "hGlobalTPCITSMatchingProbability",                   // kTrk5
237       "hGlobalTPCTOFMatchingProbability",                   // kTrk6
238       "hGlobalTPCsideAposDCA",                              // kTrk7
239       "hGlobalTPCsideAnegDCA",                              // kTrk8
240       "hGlobalTPCsideCposDCA",                              // kTrk9
241       "hGlobalTPCsideCnegDCA"                               // kTrk10
242   };
243     const Char_t *title[]={
244       "Track azimuthal distribution (rad)",                   // kTrk0
245       "Track pseudo-rapidity distribution",                   // kTrk1
246       "TPC: track momentum distribution (GeV)",               // kTrk2
247       "TPC-ITS matched: track momentum distribution (GeV)",   // kTrk3
248       "TPC-TOF matched: track momentum distribution (GeV)",   // kTrk4
249       "TPC-ITS track-matching probability",                   // kTrk5
250       "TPC-TOF track-matching probability",                   // kTrk6
251       "TPC side A: DCA for the positive tracks (mm)",         // kTrk7
252       "TPC side A: DCA for the negative tracks (mm)",         // kTrk8
253       "TPC side C: DCA for the positive tracks (mm)",         // kTrk9
254       "TPC side C: DCA for the negative tracks (mm)"          // kTrk10
255     };
256   Add2ESDsList(new TH1F(name[0],title[0],100, 0.,TMath::TwoPi()),kTrk0, !expert, image);
257   Add2ESDsList(new TH1F(name[1],title[1],100,-2.00,2.00),kTrk1, !expert, image);
258   Add2ESDsList(new TH1F(name[2],title[2],50,  0.20,5.00),kTrk2, !expert, image);
259   Add2ESDsList(new TH1F(name[3],title[3],50,  0.20,5.00),kTrk3, !expert, image);
260   Add2ESDsList(new TH1F(name[4],title[4],50,  0.20,5.00),kTrk4, !expert, image);
261   Add2ESDsList(new TH1F(name[5],title[5],50,  0.20,5.00),kTrk5, !expert, image);
262   Add2ESDsList(new TH1F(name[6],title[6],50,  0.20,5.00),kTrk6, !expert, image);
263   Add2ESDsList(new TH1F(name[7],title[7],50, -25.0,25.0),kTrk7, !expert, image);
264   Add2ESDsList(new TH1F(name[8],title[8],50, -25.0,25.0),kTrk8, !expert, image);
265   Add2ESDsList(new TH1F(name[9],title[9],50, -25.0,25.0),kTrk9, !expert, image);
266   Add2ESDsList(new TH1F(name[10],title[10],50, -25.0,25.0),kTrk10, !expert, image);
267   }
268
269   {// V0 related QA
270     const Char_t *name[]={
271       "hGlobalPromptK0sMass",
272       "hGlobalOfflineK0sMass",
273       "hGlobalPromptLambda0Lambda0BarMass",
274       "hGlobalOfflineLambda0Lambda0BarMass"
275     };
276     const Char_t *title[]={
277       "On-the-fly K0s mass (GeV)",
278       "Offline K0s mass (GeV)",
279       "On-the-fly Lambda0 + Lambda0Bar mass (GeV)",
280       "Offline Lambda0 + Lambda0Bar mass (GeV)"
281     };
282     Add2ESDsList(new TH1F(name[0],title[0],50,  0.4477,0.5477),kK0on, !expert, image);
283     Add2ESDsList(new TH1F(name[1],title[1],50,  0.4477,0.5477),kK0off, !expert, image);
284     Add2ESDsList(new TH1F(name[2],title[2],50,  1.0657,1.1657),kL0on, !expert, image);
285     Add2ESDsList(new TH1F(name[3],title[3],50,  1.0657,1.1657),kL0off, !expert, image);
286   }
287
288   {// PID related QA
289   const Char_t *name[]={
290     "hGlobalITSdEdx",
291     "hGlobalTPCdEdx",
292     "hGlobalTOFTrackingvsMeasured",
293     "hGlobalTPCdEdxvsMomentum"
294    };
295     const Char_t *title[]={
296       "ITS: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
297       "TPC: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
298       "TOF: tracking - measured (ps)",
299       "TPC: dEdx (A.U.) vs momentum (GeV)"
300      };
301     Add2ESDsList(new TH1F(name[0],title[0],50,0.00,200.),kPid0, !expert, image);
302     Add2ESDsList(new TH1F(name[1],title[1],50,0.00,100.),kPid1, !expert, image);
303     Add2ESDsList(new TH1F(name[2],title[2],50,-3500.,3500.),kPid2, !expert, image);
304     Add2ESDsList(new TH2F(name[3],title[3],1500,0.05,15.,700,0.,700.),kPid3,!expert,image);
305    }
306   {// Multiplicity related QA
307     const Char_t *name[]={
308       "hGlobalV0AvsITS",
309       "hGlobalV0CvsITS"
310     };
311     const Char_t *title[]={
312       "Multiplicity: V0A vs ITS",
313       "Multiplicity: V0C vs ITS"
314     };
315     TH2F *h0=new TH2F(name[0],title[0],41,-0.5,40.5, 33,-0.5,32.5);
316     Add2ESDsList(h0,kMlt0, !expert, image);
317     TH2F *h1=new TH2F(name[1],title[1],41,-0.5,40.5, 33,-0.5,32.5);
318     Add2ESDsList(h1,kMlt1, !expert, image);
319   }
320   //
321   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
322 }
323
324 //____________________________________________________________________________
325 void AliGlobalQADataMaker::MakeRaws(AliRawReader* rawReader)
326 {
327   //Fill prepared histograms with Raw digit properties
328   rawReader->Reset() ;
329   IncEvCountCycleRaws();
330   IncEvCountTotalRaws();
331   //
332 }
333
334 //____________________________________________________________________________ 
335 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * event) {
336   //-----------------------------------------------------------
337   // This function fills the ESD QA histograms
338   // as a part of global QA
339   //-----------------------------------------------------------
340
341   const AliESDEvent *esd=event;
342
343   // Event related QA
344   const AliESDVertex *vtx=esd->GetPrimaryVertex();
345   if (!vtx->GetStatus()) return;
346
347   Double_t xv=vtx->GetXv();
348   Double_t yv=vtx->GetYv();
349   Double_t zv=vtx->GetZv();
350   FillESDsData(kEvt0,zv);
351
352
353   Int_t ntrk=esd->GetNumberOfTracks() ; 
354   for (Int_t i=0; i<ntrk; i++) {
355     const AliESDtrack *track=esd->GetTrack(i);
356
357     // Cluster related QA
358     if (track->IsOn(AliESDtrack::kITSrefit)) {
359       Int_t n=track->GetITSclusters(0);
360       FillESDsData(kClr0,Float_t(n)/6.); //6 is the number of ITS layers
361     }
362
363     for (Int_t j=0; j<6; ++j) {
364       Int_t idet=0, sts=0;
365       Float_t xloc,zloc;
366       if (!track->GetITSModuleIndexInfo(j,idet,sts,xloc,zloc)) continue;
367       if (j>=2) idet+=240;
368       if (j>=4) idet+=260;
369       if ((sts==1)||(sts==2)||(sts==4)) FillESDsData(kClr3,idet);  
370     }
371
372     if (track->IsOn(AliESDtrack::kTPCrefit)) {
373       Int_t n =track->GetTPCNcls();
374       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
375       if (nf>0) {
376         Double_t val = n*1.0/nf; 
377         FillESDsData(kClr1,val); 
378       }
379     }
380
381     if (track->IsOn(AliESDtrack::kTRDrefit)) {
382       Int_t n=track->GetTRDclusters(0);
383       FillESDsData(kClr2,Float_t(n)/(6*24));//(6*24) is the number of TRD time bins
384     }
385
386     Double_t p=track->GetP();
387
388     // Track related QA
389     if (track->IsOn(AliESDtrack::kTPCrefit)) {
390       Float_t dz[2]; 
391       track->GetDZ(xv,yv,zv,esd->GetMagneticField(),dz); 
392       if ((TMath::Abs(dz[0])<3.) && (TMath::Abs(dz[1])<3.)) { // beam pipe
393         Double_t phi=track->Phi();
394         FillESDsData(kTrk0,phi);
395         Double_t y=track->Eta();
396         FillESDsData(kTrk1,y);
397
398         if (TMath::Abs(y)<0.9) {
399           FillESDsData(kTrk2,p);
400           if (track->IsOn(AliESDtrack::kITSrefit)) FillESDsData(kTrk3,p);
401           //if (track->IsOn(AliESDtrack::kTOFout)) FillESDsData(kTrk4,p);
402           if (track->GetTOFsignal()>0) FillESDsData(kTrk4,p);
403         }
404       }
405     }
406     const AliExternalTrackParam *tpcTrack=track->GetTPCInnerParam();
407     const AliExternalTrackParam *innTrack=track->GetInnerParam();
408     if (tpcTrack)
409     if (innTrack) {
410        Float_t dz[2];
411        tpcTrack->GetDZ(xv,yv,zv,esd->GetMagneticField(),dz);
412        dz[0]*=10.; // in mm
413        if (innTrack->GetZ()  > 0)
414        if (innTrack->GetTgl()> 0) { // TPC side A
415          if (tpcTrack->GetSign() > 0) FillESDsData(kTrk7,dz[0]);
416          else                         FillESDsData(kTrk8,dz[0]);
417        }
418        if (innTrack->GetZ()  < 0)
419        if (innTrack->GetTgl()< 0) { // TPC side C
420          if (tpcTrack->GetSign() > 0) FillESDsData(kTrk9,dz[0]);
421          else                         FillESDsData(kTrk10,dz[0]);
422        }
423     }
424
425     // PID related QA
426     if ((p>0.4) && (p<0.5)) {
427       if (track->IsOn(AliESDtrack::kITSpid)) {
428         Double_t dedx=track->GetITSsignal();
429         FillESDsData(kPid0,dedx);
430       }
431       if (track->IsOn(AliESDtrack::kTPCpid)) {
432         Double_t dedx=track->GetTPCsignal();
433         FillESDsData(kPid1,dedx);
434       }
435     }
436     if (p>1.0) {
437       if (track->IsOn(AliESDtrack::kITSrefit))
438       if (track->IsOn(AliESDtrack::kTPCrefit))
439       if (track->IsOn(AliESDtrack::kTOFout)) {
440          Float_t dz[2];
441          track->GetDZ(xv,yv,zv,esd->GetMagneticField(),dz);
442          if (dz[1]<3.) {
443            Double_t times[AliPID::kSPECIESC];
444            track->GetIntegratedTimes(times);
445            Double_t tof=track->GetTOFsignal()/*-847055 -1771207*/;
446            FillESDsData(kPid2,times[AliPID::kPion]-tof);
447          }
448       }
449     }
450     const AliExternalTrackParam *par=track->GetInnerParam();
451     if (par) {
452       Double_t pp=par->GetP();
453       Double_t dedx=track->GetTPCsignal();
454       FillESDsData(kPid3,pp,dedx);
455     }
456  
457   }
458
459   // Multiplicity related QA
460   AliESDVZERO     *mltV0 =esd->GetVZEROData();
461   const AliMultiplicity *mltITS=esd->GetMultiplicity();
462   if (mltV0)
463     if (mltITS) {
464        Short_t nv0a=mltV0->GetNbPMV0A();
465        Short_t nv0c=mltV0->GetNbPMV0C();
466        Int_t   nits=mltITS->GetNumberOfTracklets();
467        FillESDsData(kMlt0,nits,nv0a);
468        FillESDsData(kMlt1,nits,nv0c);
469     }
470
471   // RS
472   for (int itr = -1; itr<GetNEventTrigClasses(); itr++) {
473     TH1 *tpc = GetMatchingESDsHisto(kTrk2,itr);
474     TH1 *its = GetMatchingESDsHisto(kTrk3,itr);
475     TH1 *tof = GetMatchingESDsHisto(kTrk4,itr);
476     TH1* h5 =  GetMatchingESDsHisto(kTrk5,itr);
477     TH1* h6 =  GetMatchingESDsHisto(kTrk6,itr);
478     if (h5 && h6 && tpc && its && tof) {
479       tpc->Sumw2();
480       its->Sumw2();
481       tof->Sumw2();
482       h5->Divide(its,tpc,1,1.,"b");
483       h6->Divide(tof,tpc,1,1.,"b");
484     }
485   }
486
487   // V0 related QA
488   Int_t nV0=esd->GetNumberOfV0s();
489   for (Int_t i=0; i<nV0; i++) {
490     Double_t mass;
491     AliESDv0 v0(*esd->GetV0(i));
492
493     Int_t nidx=TMath::Abs(v0.GetNindex());
494     AliESDtrack *ntrack1=esd->GetTrack(nidx);
495     if (!ntrack1->IsOn(AliESDtrack::kTPCrefit)) continue;
496
497     Int_t pidx=TMath::Abs(v0.GetPindex());
498     AliESDtrack *ptrack1=esd->GetTrack(pidx);
499     if (!ptrack1->IsOn(AliESDtrack::kTPCrefit)) continue;
500
501     v0.ChangeMassHypothesis(kK0Short);
502     mass=v0.GetEffMass();
503     if (v0.GetOnFlyStatus())
504       FillESDsData(kK0on,mass);
505     else
506       FillESDsData(kK0off,mass);
507
508     v0.ChangeMassHypothesis(kLambda0);
509     mass=v0.GetEffMass();
510     if (v0.GetOnFlyStatus())
511       FillESDsData(kL0on,mass);
512     else
513       FillESDsData(kL0off,mass);
514
515     v0.ChangeMassHypothesis(kLambda0Bar);
516     mass=v0.GetEffMass();
517     if (v0.GetOnFlyStatus())
518       FillESDsData(kL0on,mass);
519     else
520       FillESDsData(kL0off,mass);
521   }
522   //
523   IncEvCountCycleESDs();
524   IncEvCountTotalESDs();
525   //
526 }