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