]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RALICE/wa98/Wa98Convert.cxx
16-feb-2005 NvE Support for user selectable split level and buffer size of the output...
[u/mrichter/AliRoot.git] / RALICE / wa98 / Wa98Convert.cxx
CommitLineData
5c9fc915 1// $Id$
2
3///////////////////////////////////////////////////////////////////////////
4// Class Wa98Convert
5// Conversion of Wa98 Hbook ntuple data into Wa98Event physics event structures.
6//
7// Usage example :
8// ---------------
9//
10// gSystem->Load("libHbook");
11// gSystem->Load("ralice");
12// gSystem->Load("rwa98");
13//
14// // Input file with ANALYZ produced Hbook ntuple data
15// THbookFile* f=new THbookFile("pb611258.cwn");
16// TTree* h999=(TTree*)f->Get(999);
17//
18// // Output file for the event structures
19// TFile* ofile=new TFile("run11258.root","RECREATE","WA98 data in RALICE event structure");
20// TTree* otree=new TTree("T","Data of the 1996 Pb+Pb run");
21//
22// Int_t nentries=h999->GetEntries();
23// cout << " Number of entries available : " << nentries << endl;
24// cout << endl;
25//
26// // Limit the number of entries for testing
27// nentries=300;
28//
29// // Print frequency to produce a short summary print every printfreq events
30// Int_t printfreq=10;
31//
aa1cca34 32// // Split level for the output structures
33// Int_t split=2;
34//
35// // Buffer size for the output structures
36// Int_t bsize=32000;
37//
38// Convert q(h999,split,bsize);
5c9fc915 39// q.Loop(otree,nentries,printfreq);
40//
41// otree->Print();
42//
43// // Close output file
44// ofile->Write();
45// ofile->Close();
46//
47//--- Author: Nick van Eijndhoven 06-jul-2004 Utrecht University
48//- Modified: NvE $Date$ Utrecht University
49///////////////////////////////////////////////////////////////////////////
50
51#include "Wa98Convert.h"
52#include "Riostream.h"
53
54ClassImp(Wa98Convert) // Class implementation to enable ROOT I/O
55
aa1cca34 56Wa98Convert::Wa98Convert(TTree* tree,Int_t split,Int_t bsize)
5c9fc915 57{
58// Default constructor.
59// Initialise the input tree (or chain) to be converted.
aa1cca34 60// Also the required split level and buffer size of the output tree
61// can be specified in this constructor.
62// By default tree=0, split=0 and bsize=32000.
5c9fc915 63
64 fChain=tree;
aa1cca34 65 fSplit=split;
66 fBsize=bsize;
67 if (!fChain || split<0 || bsize<0) return;
5c9fc915 68
69 // Link the variables to the branches of the input tree/chain
70 fChain->SetBranchAddress("Jrun",&Jrun);
71 fChain->SetBranchAddress("Jevt",&Jevt);
72 fChain->SetBranchAddress("Jdate",&Jdate);
73 fChain->SetBranchAddress("Jtime",&Jtime);
74 fChain->SetBranchAddress("Jevid",&Jevid);
75 fChain->SetBranchAddress("Jwscal",&Jwscal);
76 fChain->SetBranchAddress("Itword",&Itword);
77 fChain->SetBranchAddress("Zdc",&Zdc);
78 fChain->SetBranchAddress("Emir",&Emir);
79 fChain->SetBranchAddress("Emire",&Emire);
80 fChain->SetBranchAddress("Emirh",&Emirh);
81 fChain->SetBranchAddress("Etm",&Etm);
82 fChain->SetBranchAddress("Etme",&Etme);
83 fChain->SetBranchAddress("Etmh",&Etmh);
84 fChain->SetBranchAddress("Nmod",&Nmod);
85 fChain->SetBranchAddress("Irowl",Irowl);
86 fChain->SetBranchAddress("Icoll",Icoll);
87 fChain->SetBranchAddress("Adcl",Adcl);
88 fChain->SetBranchAddress("Nclu",&Nclu);
89 fChain->SetBranchAddress("Irowc",Irowc);
90 fChain->SetBranchAddress("Icolc",Icolc);
91 fChain->SetBranchAddress("Adcc",Adcc);
92 fChain->SetBranchAddress("Ncluv",&Ncluv);
93 fChain->SetBranchAddress("Iadccv",Iadccv);
94 fChain->SetBranchAddress("Thetacv",Thetacv);
95 fChain->SetBranchAddress("Phicv",Phicv);
96}
97///////////////////////////////////////////////////////////////////////////
98Wa98Convert::~Wa98Convert()
99{
100// Default destructor.
101}
102///////////////////////////////////////////////////////////////////////////
103void Wa98Convert::Loop(TTree* otree,Int_t nentries,Int_t printfreq)
104{
105// Loop over the specified number of entries and convert the
106// ntuple data into the Wa98Event structure.
107// The output will be written on the output tree specified as "otree".
108// If otree=0, a default standard output tree will be created.
109// If nentries<0 (default) all the entries of the input chain
110// will be processed.
111// Every "printfreq" events a short event summary will be printed.
112// The default value is printfreq=1.
113
aa1cca34 114 if (fChain==0 || fSplit<0) return;
5c9fc915 115
116 if (nentries<0) nentries=fChain->GetEntriesFast();
117
118 if (!otree) otree=new TTree("T","Data of the 1996 Pb+Pb run");
119
120 Double_t pi=acos(-1.);
121
5c9fc915 122 Wa98Event* evt=new Wa98Event();
123
124 // Branch in the tree for the event structure
aa1cca34 125 otree->Branch("Wa98Event","Wa98Event",&evt,fBsize,fSplit);
5c9fc915 126
127 // The LEDA specific output data
128 AliCalorimeter* ledaup=new AliCalorimeter(44,144);
129 AliCalorimeter* ledalw=new AliCalorimeter(40,144);
130
131 ledaup->SetName("LedaUp");
132 ledalw->SetName("LedaDown");
133
134 evt->InitLeda(ledaup);
135 evt->InitLeda(ledalw);
136
137 TDatime datim;
138 Float_t pos[3],costh;
139 AliSignal s;
140 s.SetName("CPV signal ADC");
141
142 for (Int_t jentry=0; jentry<nentries; jentry++)
143 {
144 fChain->GetEntry(jentry);
145
146 // Reset the complete Event structure
147 evt->Reset();
148
149 evt->SetRunNumber(Jrun);
150 evt->SetEventNumber(Jevt);
151 datim.Set(Jdate,Jtime);
152 evt->SetDayTime(datim);
153 evt->SetProjectile(207,82,158);
154 evt->SetTarget(207,82,0);
155 evt->SetWeight(Jwscal);
156 evt->SetTrig(Itword);
157 evt->SetZdc(Zdc*1000.);
158 evt->SetMiracE(1000.*Emir,Emire,Emirh);
159 evt->SetMiracEt(Etm,Etme,Etmh);
160
161 ledaup->Reset();
162 ledalw->Reset();
163 // Fill calorimeter with module data
164 for (Int_t i=0; i<Nmod; i++)
165 {
166 if (Adcl[i] > 3) // Adc cut of 3 to remove noise
167 {
168 if (Irowl[i] > 0) ledaup->SetSignal(Irowl[i],Icoll[i],Adcl[i]);
169 if (Irowl[i] < 0) ledalw->SetSignal(-Irowl[i],Icoll[i],Adcl[i]);
170 }
171 }
172
173 // Store associated CPV signals
174 for (Int_t j=0; j<Ncluv; j++)
175 {
176 s.Reset();
177 s.SetSignal(Iadccv[j]);
178 pos[1]=Thetacv[j]*pi/180.;
179 pos[2]=Phicv[j]*pi/180.;
180 costh=cos(pos[1]);
181 pos[0]=0;
182 if (costh) pos[0]=2103./costh;
183 s.SetPosition(pos,"sph");
184 pos[0]=0.4;
185 pos[1]=2.2;
186 pos[2]=0;
187 s.SetPositionErrors(pos,"car");
188 if (Phicv[j]>=0. && Phicv[j]<=180.)
189 {
190 ledaup->AddVetoSignal(s);
191 }
192 else
193 {
194 ledalw->AddVetoSignal(s);
195 }
196 }
197
198 evt->AddDevice(ledaup);
199 evt->AddDevice(ledalw);
200
201 if (!(jentry%printfreq))
202 {
203// cout << " Entry : " << jentry << " Run : " << Jrun << " Event : " << Jevt
204// << " Itword : " << Itword << " Etm : " << Etm << endl;
205// cout << " Jdate : " << Jdate << " Jtime : " << Jtime << endl;
206 cout << " Itword : " << Itword << " Nmod : " << Nmod << " Ncluv : " << Ncluv << endl;
207 evt->HeaderData();
208 }
209
210 // Write the complete structure to the output Tree
211 otree->Fill();
212 }
213
214 if (evt) delete evt;
215 if (ledaup) delete ledaup;
216 if (ledalw) delete ledalw;
217}
218///////////////////////////////////////////////////////////////////////////