]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MONITOR/AliMonitorV0s.cxx
Removing warnings (icc), adding more detailed description
[u/mrichter/AliRoot.git] / MONITOR / AliMonitorV0s.cxx
CommitLineData
04fa961a 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"
c4bd737c 29#include "AliRunLoader.h"
30#include <TFolder.h>
31#include <TTree.h>
04fa961a 32#include <TPDGCode.h>
33
34
35ClassImp(AliMonitorV0s)
36
37
38//_____________________________________________________________________________
39AliMonitorV0s::AliMonitorV0s()
40{
41// create a monitor object for V0s
42
43}
44
c4bd737c 45//_____________________________________________________________________________
46AliMonitorV0s::AliMonitorV0s(const AliMonitorV0s& monitor) :
47 AliMonitor(monitor)
48{
49 Fatal("AliMonitorV0s", "copy constructor not implemented");
50}
51
52//_____________________________________________________________________________
53AliMonitorV0s& AliMonitorV0s::operator = (const AliMonitorV0s& /*monitor*/)
54{
55 Fatal("operator =", "assignment operator not implemented");
56 return *this;
57}
58
04fa961a 59
60//_____________________________________________________________________________
61void AliMonitorV0s::CreateHistos(TFolder* folder)
62{
63// create the V0s monitor histograms
64
65 fFolder = folder->AddFolder("V0s", "V0s");
66
67 fRadius = CreateHisto1("Radius", "radius of V0 vertices",
68 90, 0., 3., "r_{xy} [cm]",
69 "#Delta N/N", AliMonitorHisto::kNormEvents);
70
71 fMassK0 = CreateHisto1("MassK0", "invariant mass of K^{0} candidates",
72 50, 0.4, 0.6, "m_{#pi^{+}#pi^{-}} [GeV/c^{2}]",
73 "#Delta N/N", AliMonitorHisto::kNormEvents);
74
75 fMassLambda = CreateHisto1("MassLambda",
76 "invariant mass of #Lambda candidates",
77 50, 1.0, 1.2, "m_{p#pi^{-}} [GeV/c^{2}]",
78 "#Delta N/N", AliMonitorHisto::kNormEvents);
79
80 fMassAntiLambda = CreateHisto1("MassAntiLambda",
81 "invariant mass of #bar{#Lambda} candidates",
82 50, 1.0, 1.2,
83 "m_{#bar{p}#pi^{+}} [GeV/c^{2}]",
84 "#Delta N/N", AliMonitorHisto::kNormEvents);
85}
86
87
88//_____________________________________________________________________________
89void AliMonitorV0s::FillHistos(AliRunLoader* runLoader,
90 AliRawReader*)
91{
92// fill the TPC-ITS correlation monitor histogrms
93
94 AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
95 if (!itsLoader) return;
96
97 itsLoader->LoadV0s();
98 TTree* v0s = itsLoader->TreeV0();
99 if (!v0s) return;
100 AliV0vertex* vertex = new AliV0vertex;
101 v0s->SetBranchAddress("vertices", &vertex);
102
103 for (Int_t i = 0; i < v0s->GetEntries(); i++) {
104 v0s->GetEntry(i);
105 Double_t x, y, z;
106 vertex->GetXYZ(x, y, z);
107 fRadius->Fill(TMath::Sqrt(x*x + y*y));
108 vertex->ChangeMassHypothesis(kK0Short);
109 fMassK0->Fill(vertex->GetEffMass());
110 vertex->ChangeMassHypothesis(kLambda0);
111 fMassLambda->Fill(vertex->GetEffMass());
112 vertex->ChangeMassHypothesis(kLambda0Bar);
113 fMassAntiLambda->Fill(vertex->GetEffMass());
114 }
115
116 delete vertex;
117 itsLoader->UnloadV0s();
118}