]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0.cxx
More changes to compile on Windows/Cygwin
[u/mrichter/AliRoot.git] / T0 / AliT0.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 //  T0 (T-Zero) Detector                                            //
21 //  This class contains the base procedures for the T0     //
22 //  detector                                                                 //
23 //                                                                           //
24 //Begin_Html
25 /*
26 <img src="gif/AliT0Class.gif">
27 </pre>
28 <br clear=left>
29 <font size=+2 color=red>
30 <p>The responsible person for this module is
31 <a href="mailto:Alla.Maevskaia@cern.ch">Alla Maevskaia</a>.
32 </font>
33 <pre>
34 */
35 //End_Html
36 //                                                                           //
37 //                                                                           //
38 ///////////////////////////////////////////////////////////////////////////////
39
40 //#include <Riostream.h>
41
42 //#include <TFile.h>
43 #include <TGeometry.h>
44 //#include <TMath.h>
45 #include <TNode.h>
46 //#include <TParticle.h>
47 //#include <TRandom.h>
48 #include <TTUBE.h>
49 //#include <TVirtualMC.h>
50
51 #include "AliLog.h"
52 #include "AliMC.h"
53 #include "AliLoader.h"
54 #include "AliRun.h"
55 #include "TClonesArray.h"
56 #include "AliT0.h"
57 //#include "AliT0Loader.h"
58 #include "AliT0digit.h"
59 #include "AliT0hit.h"
60 #include "AliT0Digitizer.h"
61 #include "AliT0RawData.h"
62 #include "AliT0RecPoint.h"
63 //#include "AliT0Parameters.h"
64 #include "AliLog.h"
65
66 ClassImp(AliT0)
67
68   //static  AliT0digit *digits; 
69
70 //_____________________________________________________________________________
71 AliT0::AliT0()
72   : AliDetector(), fIdSens(0), fDigits(NULL), fRecPoints(NULL)
73 {
74   //
75   // Default constructor for class AliT0
76   //
77   fIshunt   = 1;
78   fHits     = 0;
79   fDigits   = 0;
80   fRecPoints = 0;
81 }
82  
83 //_____________________________________________________________________________
84 AliT0::AliT0(const char *name, const char *title)
85   : AliDetector(name,title), fIdSens(0), fDigits(new AliT0digit()), fRecPoints(new AliT0RecPoint())
86 {
87   //
88   // Standard constructor for T0 Detector
89   //
90
91   
92   //
93   // Initialise Hit array
94   fHits       = new TClonesArray("AliT0hit",  405);
95   gAlice->GetMCApp()->AddHitList(fHits);
96   //  fDigits    = new AliT0digit();
97   //  fRecPoints = new AliT0RecPoint();
98   fIshunt     =  1;
99   //  fIdSens   =  0;
100   //PH  SetMarkerColor(kRed);
101 }
102
103 //_____________________________________________________________________________
104 AliT0::~AliT0() {
105   
106   //destructor
107   if (fHits) {
108     fHits->Delete();
109     delete fHits;
110   }
111   /*
112   if (fDigits) {
113     fDigits->Delete();
114     delete fDigits;
115     cout<<" delete fDigits; "<<endl;
116   }
117   if (fRecPoints) {
118    fRecPoints ->Delete();
119     delete fRecPoints;
120     cout<<" delete fRecPoints; "<<endl;
121   }
122   */ 
123 }
124
125 //_____________________________________________________________________________
126 void AliT0::AddHit(Int_t track, Int_t *vol, Float_t *hits)
127 {
128   //
129   // Add a T0 hit
130   //
131   TClonesArray &lhits = *fHits;
132   new(lhits[fNhits++]) AliT0hit(fIshunt,track,vol,hits);
133 }
134
135
136 //_____________________________________________________________________________
137
138 void AliT0::AddDigit(Int_t besttimeright, Int_t besttimeleft, Int_t meantime, 
139                      Int_t timediff, Int_t sumMult, Int_t refpoint,
140                         TArrayI *timeCFD, TArrayI *qt0, TArrayI *timeLED, TArrayI *qt1)
141 {
142   
143   //  Add a T0 digit to the list.
144  //
145   
146   if (!fDigits) {
147     fDigits = new AliT0digit();
148   }
149   fDigits-> SetTimeBestA(besttimeright);
150   fDigits->SetTimeBestC(besttimeleft);
151   fDigits-> SetMeanTime(meantime);
152   fDigits-> SetDiffTime(timediff);
153   fDigits-> SetSumMult(sumMult);
154   fDigits->SetTimeCFD(*timeCFD);
155   fDigits->SetTimeLED(*timeLED);
156   fDigits->SetQT0(*qt0);
157   fDigits->SetQT1(*qt1);
158   fDigits->SetRefPoint(refpoint);
159 }
160
161
162 //_____________________________________________________________________________
163 void AliT0::BuildGeometry()
164 {
165   //
166   // Build simple ROOT TNode geometry for event display
167   //
168   TNode *node, *top;
169   const int kColorT0  = 19;
170
171   top=gAlice->GetGeometry()->GetNode("alice");
172
173   // T0 define the different volumes
174   new TRotMatrix("rotx999","rot999",  90,0,90,90,180,0);
175
176   new TTUBE("S_0ST1","T0  volume 1","void",5.,10.7,5.3);
177   top->cd();
178   node = new TNode("0ST1","0ST01","S_0ST1",0,0,-69.7,"");
179   node->SetLineColor(kColorT0);
180   fNodes->Add(node);
181
182   new TTUBE("S_0ST2","T0 volume 2","void",5.,10.7,5.3);
183   top->cd();
184   node = new TNode("0ST2","0ST2","S_0ST2",0,0,350,"rotx999");
185   node->SetLineColor(kColorT0);
186   fNodes->Add(node);
187 }
188  
189 //_____________________________________________________________________________
190 Int_t AliT0::DistanceToPrimitive(Int_t /*px*/, Int_t /*py*/)
191 {
192   //
193   // Calculate the distance from the mouse to the T0 on the screen
194   // Dummy routine
195   //
196   return 9999;
197 }
198  
199 //-------------------------------------------------------------------------
200 void AliT0::Init()
201 {
202   //
203   // Initialis the T0 after it has been built
204   Int_t i;
205   //
206   if(AliLog::GetGlobalDebugLevel()>0) {
207     printf("\n%s: ",ClassName());
208     for(i=0;i<35;i++) printf("*");
209     printf(" T0_INIT ");
210     for(i=0;i<35;i++) printf("*");
211     printf("\n%s: ",ClassName());
212     //
213     // Here the T0 initialisation code (if any!)
214     for(i=0;i<80;i++) printf("*");
215     printf("\n");
216   }
217 }
218
219 //---------------------------------------------------------------------------
220 void AliT0::MakeBranch(Option_t* option)
221 {
222   //
223 // Create Tree branches for the T0.
224
225  // Options:
226   //
227   //    H          Make a branch of TClonesArray of AliT0Hit's
228   //    D          Make a branch of TClonesArray of AliT0Digit's
229   //
230   //    R         Make a branch of  AliT0RecPoints
231   //
232   char branchname[20];
233   sprintf(branchname,"%s",GetName());
234
235   const char *cH = strstr(option,"H");
236   const char *cD = strstr(option,"D");
237   const char *cR = strstr(option,"R");
238
239     if (cH && fLoader->TreeH())
240   {
241      if (fHits == 0x0) fHits  = new TClonesArray("AliT0hit",  405);
242      AliDetector::MakeBranch(option);
243   } 
244     
245     
246   if (cD && fLoader->TreeD())
247     {
248       if (fDigits == 0x0) fDigits  = new AliT0digit();
249       //     MakeBranchInTree(fLoader->TreeD(), branchname,
250       //                       &fDigits, 405, 0);
251       fLoader->TreeD()->Branch(branchname,"AliT0digit",&fDigits,405,1);
252       //   fLoader->TreeD()->Print();
253     } 
254   if (cR && fLoader->TreeR())
255     {
256       if (fRecPoints == 0x0) fRecPoints  = new AliT0RecPoint();
257       MakeBranchInTree(fLoader->TreeR(), branchname,
258                        &fRecPoints, 405, 0);
259     } 
260   
261 }    
262
263 //_____________________________________________________________________________
264 void AliT0::ResetHits()
265 {
266   //
267   //reset hits
268   //
269   AliDetector::ResetHits();
270   
271 }
272 //____________________________________________________________________
273 void AliT0::ResetDigits()
274 {
275   //
276   // Reset number of digits and the digits array for this detector
277   //
278   if (fDigits) fDigits->Clear();
279 }
280
281 //_____________________________________________________________________________
282 void AliT0::SetTreeAddress()
283 {
284
285   TTree    *treeH;
286   treeH = TreeH();
287   
288   if (treeH)
289     {
290       if (fHits == 0x0) fHits  = new TClonesArray("AliT0hit",  405);
291     }
292     
293   AliDetector::SetTreeAddress();
294   TTree *treeD = fLoader->TreeD();
295   if (treeD) {
296     if (fDigits == 0x0)  fDigits  = new AliT0digit();
297     TBranch* branch = treeD->GetBranch ("T0");
298     if (branch) branch->SetAddress(&fDigits);
299   }
300
301   TTree *treeR = fLoader->TreeR();
302   if (treeR) {
303     if (fRecPoints == 0x0) fRecPoints  = new  AliT0RecPoint()  ;
304     TBranch* branch = treeR->GetBranch ("T0");
305     if (branch) branch->SetAddress(&fRecPoints);
306   }
307  
308 }
309
310
311 //_____________________________________________________________________________
312 void AliT0::MakeBranchInTreeD(TTree *treeD, const char *file)
313 {
314     //
315     // Create TreeD branches for the FMD
316     //
317     const Int_t kBufferSize = 4000;
318     char branchname[20];
319     sprintf(branchname,"%s",GetName());
320     if(treeD)
321      {
322        MakeBranchInTree(treeD,  branchname,&fDigits, kBufferSize, file);
323      }
324 }
325
326 //_____________________________________________________________________________
327 AliDigitizer* AliT0::CreateDigitizer(AliRunDigitizer* manager) const
328 {
329   return new AliT0Digitizer(manager);
330 }
331 //____________________________________________________________________________
332 void AliT0::Digits2Raw()
333 {
334 //
335 // Starting from the T0 digits, writes the Raw Data objects
336 //
337 //  AliT0Loader* pStartLoader = (AliT0Loader*)fLoader;
338   fLoader ->LoadDigits("read");
339   TTree* treeD = fLoader->TreeD();
340   if (!treeD) {
341     AliError("no digits tree");
342     return;
343   }
344   if (fDigits == 0x0)  fDigits  = new AliT0digit();
345   
346   TBranch *branch = treeD->GetBranch("T0");
347   if (branch) {
348     branch->SetAddress(&fDigits);
349   }else{
350     AliError("Branch T0 DIGIT not found");
351     exit(111);
352   } 
353   AliT0RawData rawWriter;
354   rawWriter.SetVerbose(0);
355   
356   AliDebug(2,Form(" Formatting raw data for T0 "));
357   branch->GetEntry(0);
358   //  rawWriter.RawDataT0(treeD->GetBranch("T0"));
359   rawWriter.RawDataT0(fDigits);
360   
361   
362   fLoader->UnloadDigits();
363   
364 }
365
366 //____________________________________________________________________________
367 void AliT0::Raw2Digits(AliRawReader *rawReader,TTree* digitsTree)
368 {
369
370  //T0 raw data-> digits conversion
371  // reconstruct time information from raw data
372  // cout<<"  AliT0::Raw2Digits(AliRawReader *rawReader,TTree* digitsTree) "<<
373   // rawReader<<" "<<digitsTree<<endl;
374
375  
376   //  AliT0RawReader myrawreader(rawReader,digitsTree);
377    AliT0RawReader myrawreader(rawReader);
378    if (!myrawreader.Next())
379      AliDebug(1,Form(" no raw data found!! %i", myrawreader.Next()));
380    Int_t allData[110][5];
381    for (Int_t i=0; i<110; i++) {
382      allData[i][0]=myrawreader.GetData(i,0);
383    }
384
385
386    AliT0digit* fDigits = new AliT0digit();
387    digitsTree->Branch("T0","AliT0digit",&fDigits,405,1);
388    
389    
390    TArrayI *timeLED = new TArrayI(24);
391    TArrayI * timeCFD = new TArrayI(24);
392    TArrayI *chargeQT0 = new TArrayI(24);
393    TArrayI *chargeQT1 = new TArrayI(24);
394    
395    for (Int_t in=0; in<24; in++)
396      {
397        timeLED->AddAt(allData[in+1][0],in);
398        timeCFD->AddAt(allData[in+25][0],in);
399        chargeQT0->AddAt(allData[in+55][0],in);
400        chargeQT1->AddAt(allData[in+79][0],in);
401        AliDebug(2, Form(" readed Raw %i %i %i %i %i", in, timeLED->At(in),timeCFD->At(in),chargeQT0->At(in),chargeQT1->At(in)));
402      }
403   
404    fDigits->SetTimeCFD(*timeCFD);
405    fDigits->SetQT0(*chargeQT1);
406
407    fDigits->SetTimeLED(*timeLED);
408    fDigits->SetQT1(*chargeQT1);
409
410    fDigits->SetMeanTime(allData[49][0]);
411    fDigits->SetDiffTime(allData[50][0]);
412    fDigits->SetTimeBestA(allData[51][0]);
413    fDigits->SetTimeBestC(allData[52][0]);
414    digitsTree->Fill();
415    fDigits->Write();
416  
417    delete timeCFD ;
418    delete chargeQT0;
419    delete timeLED ;
420    delete chargeQT1;
421
422
423 }