]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/AliMonitorV0s.cxx
Removing unused include
[u/mrichter/AliRoot.git] / MONITOR / AliMonitorV0s.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  This class creates and fills the monitor histograms for V0s              //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24
25 #include "AliMonitorV0s.h"
26 #include "AliMonitorHisto.h"
27 #include "AliITSLoader.h"
28 #include "AliV0vertex.h"
29 #include <TPDGCode.h>
30
31
32 ClassImp(AliMonitorV0s) 
33
34
35 //_____________________________________________________________________________
36 AliMonitorV0s::AliMonitorV0s()
37 {
38 // create a monitor object for V0s
39
40 }
41
42
43 //_____________________________________________________________________________
44 void AliMonitorV0s::CreateHistos(TFolder* folder)
45 {
46 // create the V0s monitor histograms
47
48   fFolder = folder->AddFolder("V0s", "V0s");
49
50   fRadius = CreateHisto1("Radius", "radius of V0 vertices", 
51                          90, 0., 3., "r_{xy} [cm]", 
52                          "#Delta N/N", AliMonitorHisto::kNormEvents);
53
54   fMassK0 = CreateHisto1("MassK0", "invariant mass of K^{0} candidates", 
55                          50, 0.4, 0.6, "m_{#pi^{+}#pi^{-}} [GeV/c^{2}]", 
56                          "#Delta N/N", AliMonitorHisto::kNormEvents);
57
58   fMassLambda = CreateHisto1("MassLambda", 
59                              "invariant mass of #Lambda candidates", 
60                              50, 1.0, 1.2, "m_{p#pi^{-}} [GeV/c^{2}]", 
61                              "#Delta N/N", AliMonitorHisto::kNormEvents);
62
63   fMassAntiLambda = CreateHisto1("MassAntiLambda", 
64                                  "invariant mass of #bar{#Lambda} candidates", 
65                                  50, 1.0, 1.2, 
66                                  "m_{#bar{p}#pi^{+}} [GeV/c^{2}]", 
67                                  "#Delta N/N", AliMonitorHisto::kNormEvents);
68 }
69
70
71 //_____________________________________________________________________________
72 void AliMonitorV0s::FillHistos(AliRunLoader* runLoader, 
73                                AliRawReader*)
74 {
75 // fill the TPC-ITS correlation monitor histogrms
76
77   AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
78   if (!itsLoader) return;
79
80   itsLoader->LoadV0s();
81   TTree* v0s = itsLoader->TreeV0();
82   if (!v0s) return;
83   AliV0vertex* vertex = new AliV0vertex;
84   v0s->SetBranchAddress("vertices", &vertex);
85
86   for (Int_t i = 0; i < v0s->GetEntries(); i++) {
87     v0s->GetEntry(i);
88     Double_t x, y, z;
89     vertex->GetXYZ(x, y, z);
90     fRadius->Fill(TMath::Sqrt(x*x + y*y));
91     vertex->ChangeMassHypothesis(kK0Short); 
92     fMassK0->Fill(vertex->GetEffMass());
93     vertex->ChangeMassHypothesis(kLambda0); 
94     fMassLambda->Fill(vertex->GetEffMass());
95     vertex->ChangeMassHypothesis(kLambda0Bar); 
96     fMassAntiLambda->Fill(vertex->GetEffMass());
97   }
98
99   delete vertex;
100   itsLoader->UnloadV0s();
101 }