]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - START/AliSTARTReconstructor.cxx
Avoid compiler warnings.
[u/mrichter/AliRoot.git] / START / AliSTARTReconstructor.cxx
... / ...
CommitLineData
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#include <Riostream.h>
19
20#include <TDirectory.h>
21
22#include "AliRunLoader.h"
23#include "AliRun.h"
24#include "AliSTARTLoader.h"
25#include "AliSTARTdigit.h"
26#include "AliSTARTReconstructor.h"
27#include <AliESD.h>
28#include <TClonesArray.h>
29#include "AliSTARTRecPoint.h"
30#include "AliRawReader.h"
31#include "AliSTARTRawReader.h"
32#include "AliLog.h"
33ClassImp(AliSTARTReconstructor)
34
35 void AliSTARTReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
36{
37 //START raw data-> digits conversion
38 // reconstruct time information from raw data
39 AliSTARTRawReader myrawreader(rawReader,digitsTree);
40 myrawreader.NextThing();
41}
42 void AliSTARTReconstructor::Reconstruct(TTree*digitsTree, TTree*clustersTree) const
43{
44// START digits reconstruction
45// STARTRecPoint writing
46
47 AliDebug(1,Form("Start DIGITS reconstruction "));
48 Int_t channelWigth=25; //ps
49 TArrayI* fSumMult = new TArrayI(6);
50 Float_t ph2mV = 150./500.;
51 Float_t mV2channel=200000/(25*25); //5V -> 200ns
52
53 TBranch *brDigits=digitsTree->GetBranch("START");
54 AliSTARTdigit *fDigits = new AliSTARTdigit();
55 if (brDigits) {
56 brDigits->SetAddress(&fDigits);
57 }else{
58 cerr<<"EXEC Branch START digits not found"<<endl;
59 return;
60 }
61 brDigits->GetEntry(0);
62 Int_t besttimeright=channelWigth * (fDigits->BestTimeRight());
63 Int_t besttimeleft=channelWigth * (fDigits->BestTimeLeft());
64
65 //folding with experimental time distribution
66 // Float_t c = 29.9792; // cm/ns
67 Float_t c = 0.0299792; // cm/ps
68 Float_t lenr=TMath::Sqrt(350*350 + 6.5*6.5);
69 Float_t lenl=TMath::Sqrt(69.7*69.7 + 6.5*6.5);
70 Float_t timeDiff=channelWigth * (fDigits->TimeDiff());
71 Int_t meanTime=channelWigth * (fDigits->MeanTime());
72 Float_t ds=(c*(timeDiff)-(lenr-lenl))/2;
73 AliDebug(2,Form(" timediff in ns %f real point%f",timeDiff,ds));
74
75
76 fDigits->GetSumMult(*fSumMult);
77 Int_t multipl[6];
78 for (Int_t i=0; i<6; i++)
79 {
80 Float_t mult=Float_t (fSumMult->At(i));
81 Float_t realMultmV=TMath::Exp(mult/mV2channel);
82 multipl[i]=Int_t ((realMultmV/ph2mV)/500+0.5);
83 }
84 AliDebug(2,Form(" multiplicity Abs side %i multiplicity non-Abs side %i",multipl[1],multipl[2]));
85
86 AliSTARTRecPoint* frecpoints= new AliSTARTRecPoint ();
87 clustersTree->Branch( "START", "AliSTARTRecPoint" ,&frecpoints, 405,1);
88 frecpoints->SetTimeBestRight(besttimeright);
89 frecpoints->SetTimeBestLeft(besttimeleft);
90 frecpoints->SetVertex(ds);
91 frecpoints->SetMeanTime(meanTime);
92 frecpoints->SetMult(multipl[0]);
93 frecpoints->SetMultA(multipl[2]);
94 frecpoints->SetMultC(multipl[1]);
95
96 clustersTree->Fill();
97}
98
99
100void AliSTARTReconstructor::FillESD(AliRunLoader* runLoader, AliESD *pESD) const
101{
102
103 /***************************************************
104 Resonstruct digits to vertex position
105 ****************************************************/
106
107 // Float_t c = 0.3; //speed of light mm/ps
108 // Float_t Zposition=0;
109
110 if (!runLoader) {
111 AliError("Reconstruct >> No run loader");
112 return;
113 }
114
115 AliDebug(1,Form("Start FillESD START"));
116
117 AliSTARTLoader* pStartLoader = (AliSTARTLoader*) runLoader->GetLoader("STARTLoader");
118
119 pStartLoader->LoadRecPoints("READ");
120
121 TTree *treeR = pStartLoader->TreeR();
122
123 AliSTARTRecPoint* frecpoints= new AliSTARTRecPoint ();
124 if (!frecpoints) {
125 AliError("Reconstruct Fill ESD >> no recpoints found");
126 return;
127 }
128
129 AliDebug(1,Form("Start FillESD START"));
130 TBranch *brRec = treeR->GetBranch("START");
131 if (brRec) {
132 brRec->SetAddress(&frecpoints);
133 }else{
134 cerr<<"EXEC Branch START rec not found"<<endl;
135 exit(111);
136 }
137
138 brRec->GetEntry(0);
139 Float_t Zposition=frecpoints->GetVertex();
140
141 pESD->SetT0zVertex(Zposition);
142 pStartLoader->UnloadRecPoints();
143
144} // vertex in 3 sigma
145
146
147
148
149
150