]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AD/AliAD.cxx
add the possibility to process the event stat file
[u/mrichter/AliRoot.git] / AD / AliAD.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: AliAD.cxx  $ */
17
18 ///////////////////////////////////////////////////////////////////////////
19 //                                                                       //
20 //                  AD (ALICE Diffractive)  Detector                     //
21 //                                                                       //
22 //  This class contains the base procedures for the AD  detector         //
23 //  All comments should be sent to :                                     //
24 //                                                                       //
25 //                                                                       //
26 ///////////////////////////////////////////////////////////////////////////
27
28
29 // --- Standard libraries ---
30 #include <Riostream.h>
31 #include <stdlib.h>
32
33 // --- ROOT libraries ---
34 #include <TNamed.h>
35 #include "TROOT.h"
36 #include "TFile.h"
37 #include "TNetFile.h"
38 #include "TRandom.h"
39 #include "TTree.h"
40 #include "TBranch.h"
41 #include "TClonesArray.h"
42 #include "TGeoGlobalMagField.h"
43 #include "AliMagF.h"
44 #include "TStopwatch.h"
45 #include "TParameter.h"
46 #include "TF1.h"
47
48 // --- AliRoot header files ---
49 #include "AliRun.h"
50 #include "AliMC.h"
51 #include "AliAD.h"
52 #include "AliADhit.h"
53 #include "AliADLoader.h"
54 #include "AliADDigitizer.h"
55 #include "AliADBuffer.h"
56 #include "AliADConst.h"
57 #include "AliDigitizationInput.h"
58 #include "AliADdigit.h"
59 #include "AliADSDigit.h"
60 #include "AliADCalibData.h"
61 #include "AliDAQ.h"
62 #include "AliRawReader.h"
63 #include "AliCDBManager.h"
64 #include "AliCDBEntry.h"
65 #include "AliADReconstructor.h"
66
67 ClassImp(AliAD)
68  //__________________________________________________________________
69 AliAD::AliAD()
70    : AliDetector(),
71      fSetADAToInstalled(0),
72      fSetADCToInstalled(0)
73 {
74         /// Default Constructor
75     
76
77 }
78
79 //_____________________________________________________________________________
80 AliAD::AliAD(const char *name, const char *title)
81    : AliDetector(name,title),
82      fSetADAToInstalled(kTRUE),
83      fSetADCToInstalled(kTRUE)
84
85 {
86   
87    // Standard constructor for AD Detector
88   
89
90 }
91
92 //_____________________________________________________________________________
93 AliAD::~AliAD()
94 {
95    //
96    // Default destructor for AD Detector
97    //
98   
99 }
100
101 //_____________________________________________________________________________
102 void AliAD::CreateMaterials()
103 {
104    //
105    // MATERIALS FOR ADC AND ADA
106    //
107    
108    // Parameters for simulation scope for ADA and ADC (stolen from AliVZEROv7::CreateMaterials )
109    Int_t    fieldType       = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();     // Field type 
110    Double_t maxField        = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();       // Field max.
111    Double_t maxBending      = 10;    // Max Angle
112    Double_t maxStepSize     = 0.01;  // Max step size 
113    Double_t maxEnergyLoss   = 1;     // Max Delta E
114    Double_t precision       = 0.003; // Precision
115    Double_t minStepSize     = 0.003; // Minimum step size 
116    Float_t  density,  as[11], zs[11], ws[11];
117    Double_t radLength, absLength, a_ad, z_ad;
118    Int_t    id;
119    
120    //
121    // Parameters  for AD scintillator: NE-102 (BC400)
122    //
123    // NE-102, has the following properties : (from internet, need better reference)
124    //    Density : ca. 1.03 g/cm3
125    //    Electrons/cm3: 3.39 x 10^23
126    //    H atoms/cm3: 5.28 x 10^22
127    //    C atoms/cm3: 4.78 x 10^22
128    //    Ratio of H to C : 1.104 .
129    //    wavelength of emission : ~4.23 nm.
130    //    Decay time : ~2.4 ns.
131    //    Luminescent efficiency : typically 18% of NaI(Tl)
132    //    Photons/MeV: 2.5 x 10^4 
133    //
134    // H                // C 
135    as[0] = 1.00794;    as[1] = 12.011;
136    zs[0] = 1.;         zs[1] = 6.;
137    ws[0] = 5.23;       ws[1] = 4.74;
138    density = 1.032;
139    id      = 1;
140    AliMixture( id, "NE102", as, zs, density, -2, ws );
141    AliMedium( id, "NE102", id, 1, fieldType, maxField, maxBending, maxStepSize,
142               maxEnergyLoss, precision, minStepSize );
143
144    //
145    // Parameters for lightGuide:  
146    //     TODO check material 
147    // Should be Poly(methyl methacrylate) (PMMA) acrylic 
148    // (C5O2H8)n 
149    // Density  1.18 g/cm3
150    // Mixture PMMA    Aeff=12.3994 Zeff=6.23653 rho=1.18 radlen=34.0677 intlen=63.3073
151    // Element #0 : C  Z=  6.00 A= 12.01 w= 0.600 natoms=5
152    // Element #1 : H  Z=  1.00 A=  1.01 w= 0.081 natoms=8
153    // Element #2 : O  Z=  8.00 A= 16.00 w= 0.320 natoms=2
154
155    // Carbon          Hydrogen          Oxygen
156    as[0] = 12.0107;   as[1] = 1.00794;  as[2] = 15.9994;
157    zs[0] = 6.;        zs[1] = 1.;       zs[2] = 8.;
158    ws[0] = 0.60;      ws[1] = 0.081;    ws[2] = 0.32;
159    density = 1.18;
160    id      = 2;
161    AliMixture( id, "PMMA", as, zs, density, 3, ws );
162    AliMedium( id,"PMMA", id, 1, fieldType, maxField, maxBending, maxStepSize,
163               maxEnergyLoss, precision, minStepSize );
164
165    
166    // mu-metal
167    // Niquel          Iron              Molybdenum        Manganese
168    as[0] = 58.6934;   as[1] = 55.845;   as[2] = 95.94;    as[3] = 54.9380;  
169    zs[0] = 28.;       zs[1] = 26.;      zs[2] = 42.;      zs[3] = 25.;   
170    ws[0] = 0.802;     ws[1] = 0.14079;  ws[2] = 0.0485;   ws[3] = 0.005;
171    // Silicon         Chromium          Cobalt            Aluminium
172    as[4] = 28.0855;   as[5] = 51.9961;  as[6] = 58.9332;  as[7] = 26.981539;   
173    zs[4] = 14.;       zs[5] = 24.;      zs[6] = 27.;      zs[7] = 13.;   
174    ws[4] = 0.003;     ws[5] = 0.0002;   ws[6] = 0.0002;   ws[7] = 0.0001;
175    // Carbon          Phosphorus        Sulfur
176    as[8] = 12.0107;   as[9] = 30.97376; as[10] = 32.066;
177    zs[8] = 6.;        zs[9] = 15.;      zs[10] = 16.;
178    ws[8] = 0.00015;   ws[9] = 0.00005;  ws[10] = 0.00001;
179    density = 8.25;
180    id      = 3;
181    AliMixture( id, "MuMetal", as, zs, density, 11, ws );
182    AliMedium( id,"MuMetal", id, 1, fieldType, maxField, maxBending, maxStepSize,
183               maxEnergyLoss, precision, minStepSize );
184
185    // Parameters for ADCPMA: Aluminium
186    a_ad = 26.98; 
187    z_ad = 13.00;
188    density   = 2.7;
189    radLength = 8.9;
190    absLength = 37.2;
191    id = 4;
192    AliMaterial (id, "Alum",  a_ad, z_ad, density, radLength, absLength, 0, 0 );
193    AliMedium( id, "Alum", id, 1, fieldType, maxField, maxBending, maxStepSize,
194               maxEnergyLoss, precision, minStepSize );
195
196    // Parameters for ADCPMG: Glass for the simulation Aluminium 
197    // TODO fix material
198    a_ad = 26.98; 
199    z_ad = 13.00;
200    density   = 2.7;
201    radLength = 8.9;
202    absLength = 37.2;
203    id = 5;
204    AliMaterial( id, "Glass",  a_ad, z_ad, density, radLength, absLength, 0, 0 );
205    AliMedium( id, "Glass", id, 1, fieldType, maxField, maxBending, maxStepSize,
206               maxEnergyLoss, precision, minStepSize );
207
208
209 }
210 //_____________________________________________________________________________
211 void AliAD::SetTreeAddress()
212 {
213    //
214    // Sets tree address for hits.
215    //
216
217         TBranch *branch;
218         char branchname[20];
219         snprintf(branchname,19,"%s",GetName());
220         // Branch address for hit tree
221         TTree *treeH = fLoader->TreeH();
222         if (treeH ) 
223         {
224                 branch = treeH->GetBranch(branchname);
225                 if (branch) branch->SetAddress(&fHits);
226         }
227 }
228
229
230 //_____________________________________________________________________________
231 void AliAD::MakeBranch(Option_t* opt)
232 {
233         const char* oH = strstr(opt,"H");
234         if (fLoader->TreeH() && oH && (fHits==0x0))
235         {
236                 fHits = new TClonesArray("AliADhit",1000);
237                 fNhits = 0;
238         }
239         AliDetector::MakeBranch(opt);
240 }
241 //_____________________________________________________________________________
242 AliLoader* AliAD::MakeLoader(const char* topfoldername)
243
244  
245    AliDebug(1,Form("Creating AliADLoader, Top folder is %s ",topfoldername));
246    fLoader = new AliADLoader(GetName(),topfoldername);
247    return fLoader;
248 }
249
250 //_____________________________________________________________________________
251 AliDigitizer* AliAD::CreateDigitizer(AliDigitizationInput* digInput) const
252 {
253    //
254    // Creates a digitizer for AD
255    //
256    return new AliADDigitizer(digInput);
257 }
258
259 //_____________________________________________________________________________
260 void AliAD::Hits2Digits(){
261   //
262   // Converts hits to digits
263   //
264   // Creates the AD digitizer 
265   AliADDigitizer* dig = new AliADDigitizer(this,AliADDigitizer::kHits2Digits);
266
267   // Creates the digits
268   dig->Digitize("");
269
270   // deletes the digitizer
271   delete dig;
272 }
273
274 //_____________________________________________________________________________
275 void AliAD::Hits2SDigits(){
276   //
277   // Converts hits to summable digits
278   //
279   // Creates the AD digitizer 
280   AliADDigitizer* dig = new AliADDigitizer(this,AliADDigitizer::kHits2SDigits);
281
282   // Creates the sdigits
283   dig->Digitize("");
284
285   // deletes the digitizer
286   delete dig;
287 }
288
289
290 //_____________________________________________________________________________
291 void AliAD::Digits2Raw()
292 {
293    //
294    //  Converts digits of the current event to raw data
295    //
296    AliAD *fAD = (AliAD*)gAlice->GetDetector("AD");
297    fLoader->LoadDigits();
298    TTree* digits = fLoader->TreeD();
299    if (!digits) {
300       Error("Digits2Raw", "no digits tree");
301       return;
302    }
303    TClonesArray * ADdigits = new TClonesArray("AliADdigit",1000);
304    fAD->SetTreeAddress();               
305    digits->GetBranch("ADDigit")->SetAddress(&ADdigits); 
306   
307    const char *fileName    = AliDAQ::DdlFileName("AD",0);
308    AliADBuffer* buffer  = new AliADBuffer(fileName);
309    //AliADBuffer* buffer  = new AliADBuffer("AD_0.ddl");
310    
311    // Get Trigger information first
312    // Read trigger inputs from trigger-detector object
313    AliDataLoader * dataLoader = fLoader->GetDigitsDataLoader();
314    if( !dataLoader->IsFileOpen() ) 
315         dataLoader->OpenFile( "READ" );
316    AliTriggerDetector* trgdet = (AliTriggerDetector*)dataLoader->GetDirectory()->Get( "Trigger" );
317    UInt_t triggerInfo = 0;
318    if(trgdet) {
319       triggerInfo = trgdet->GetMask() & 0xffff;
320    }
321    else {
322       AliError(Form("There is no trigger object for %s",fLoader->GetName()));
323    }
324
325    buffer->WriteTriggerInfo((UInt_t)triggerInfo); 
326    buffer->WriteTriggerScalers(); 
327    buffer->WriteBunchNumbers(); 
328   
329    // Now retrieve the channel information: charge smaples + time and 
330    // dump it into ADC and Time arrays
331    Int_t nEntries = Int_t(digits->GetEntries());
332    Short_t aADC[16][kNClocks];
333    Float_t aTime[16];
334    Float_t aWidth[16];
335    Bool_t  aIntegrator[16];
336   
337    for (Int_t i = 0; i < nEntries; i++) {
338      fAD->ResetDigits();
339      digits->GetEvent(i);
340      Int_t ndig = ADdigits->GetEntriesFast(); 
341    
342      if(ndig == 0) continue;
343      for(Int_t k=0; k<ndig; k++){
344          AliADdigit* fADDigit = (AliADdigit*) ADdigits->At(k);
345
346          Int_t iChannel       = fADDigit->PMNumber();
347          
348          for(Int_t iClock = 0; iClock < kNClocks; ++iClock) aADC[iChannel][iClock] = fADDigit->ChargeADC(iClock);
349          aTime[iChannel]      = fADDigit->Time(); //divide by resolution
350          aWidth[iChannel]     = fADDigit->Width(); //divide by resolution
351          aIntegrator[iChannel]= fADDigit->Integrator();
352          
353          //AliDebug(1,Form("DDL: %s\tdigit number: %d\tPM number: %d\tADC: %d\tTime: %f",
354                          //fileName,k,fADDigit->PMNumber(),aADC[iChannel][AliADdigit::kNClocks/2],aTime[iChannel])); 
355      }        
356    }
357
358    // Now fill raw data   
359    for (Int_t  iCIU = 0; iCIU < kNCIUBoards; iCIU++) { 
360    // decoding of one Channel Interface Unit numbered iCIU - there are 8 channels per CIU (and 2 CIUs) : 
361       for(Int_t iChannel_Offset = iCIU*8; iChannel_Offset < (iCIU*8)+8; iChannel_Offset=iChannel_Offset+4) { 
362          for(Int_t iChannel = iChannel_Offset; iChannel < iChannel_Offset+4; iChannel++) {
363              buffer->WriteChannel(iChannel, aADC[iChannel], aIntegrator[iChannel]);       
364          }
365          buffer->WriteBeamFlags(); 
366          buffer->WriteMBInfo(); 
367          buffer->WriteMBFlags();   
368          buffer->WriteBeamScalers(); 
369       } 
370
371       for(Int_t iChannel = iCIU*8 + 7; iChannel >= iCIU*8; iChannel--) {
372           buffer->WriteTiming(aTime[iChannel], aWidth[iChannel]); 
373       } // End of decoding of one CIU card    
374   } // end of decoding the eight CIUs
375
376      
377   delete buffer;
378   fLoader->UnloadDigits();  
379
380 }
381
382 //_____________________________________________________________________________
383 Bool_t AliAD::Raw2SDigits(AliRawReader* /*rawReader*/)
384 {
385         // reads raw data to produce digits
386         // for AD not implemented yet (needs detailed hardware info)
387         return kTRUE;
388 }