Line data Source code
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 "AliDigitizationInput.h"
56 : #include "AliADdigit.h"
57 : #include "AliADSDigit.h"
58 : #include "AliADBuffer.h"
59 : #include "AliDAQ.h"
60 : #include "AliRawReader.h"
61 : #include "AliCDBManager.h"
62 : #include "AliCDBEntry.h"
63 : #include "AliADReconstructor.h"
64 : #include "AliADCalibData.h"
65 : #include "AliADConst.h"
66 :
67 12 : ClassImp(AliAD)
68 : //__________________________________________________________________
69 : AliAD::AliAD()
70 0 : : AliDetector(),
71 0 : fCalibData(NULL),
72 0 : fSetADATwoInstalled(0),
73 0 : fSetADCTwoInstalled(0)
74 0 : {
75 : /// Default Constructor
76 :
77 :
78 0 : }
79 :
80 : //_____________________________________________________________________________
81 : AliAD::AliAD(const char *name, const char *title)
82 0 : : AliDetector(name,title),
83 0 : fCalibData(NULL),
84 0 : fSetADATwoInstalled(kTRUE),
85 0 : fSetADCTwoInstalled(kTRUE)
86 :
87 0 : {
88 :
89 : // Standard constructor for AD Detector
90 :
91 :
92 0 : }
93 :
94 : //_____________________________________________________________________________
95 : AliAD::~AliAD()
96 0 : {
97 : //
98 : // Default destructor for AD Detector
99 : //
100 :
101 0 : }
102 :
103 : //_____________________________________________________________________________
104 : void AliAD::CreateMaterials()
105 : {
106 : //
107 : // MATERIALS FOR ADC AND ADA
108 : //
109 :
110 : // Parameters for simulation scope for ADA and ADC (stolen from AliVZEROv7::CreateMaterials )
111 0 : Int_t fieldType = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ(); // Field type
112 0 : Double_t maxField = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max(); // Field max.
113 : // Int_t isxfld1 = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
114 : // Float_t sxmgmx = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
115 0 : Int_t isxfld2 = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->PrecInteg();
116 : Double_t maxBending = 10; // Max Angle
117 : Double_t maxStepSize = 0.01; // Max step size
118 : Double_t maxEnergyLoss = 1; // Max Delta E
119 : Double_t precision = 0.003; // Precision
120 : Double_t minStepSize = 0.003; // Minimum step size
121 0 : Float_t density, as[11], zs[11], ws[11];
122 : Double_t radLength, absLength, a_ad, z_ad;
123 : Int_t id;
124 :
125 : //
126 : // Air
127 : //
128 0 : Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
129 0 : Float_t zAir[4]={6.,7.,8.,18.};
130 0 : Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
131 : Float_t dAir = 1.20479E-3;
132 : Float_t dAir1 = 1.20479E-10;
133 : // Steel
134 0 : Float_t asteel[4] = { 55.847,51.9961, 58.6934, 28.0855 };
135 0 : Float_t zsteel[4] = { 26., 24., 28., 14. };
136 0 : Float_t wsteel[4] = { .715, .18, .1, .005 };
137 : //
138 : // Cast iron
139 : //
140 0 : Float_t acasti[4] = {55.847,12.011,28.085,54.938};
141 0 : Float_t zcasti[4] = {26.,6.,14.,25.};
142 0 : Float_t wcasti[4] = {0.929,0.035,0.031,0.005};
143 : // --- Define the various materials for GEANT ---
144 : // Aluminum
145 0 : AliMaterial(9, "ALU1 ", 26.98, 13., 2.7, 8.9, 37.2);
146 0 : AliMaterial(29, "ALU2 ", 26.98, 13., 2.7, 8.9, 37.2);
147 0 : AliMaterial(49, "ALU3 ", 26.98, 13., 2.7, 8.9, 37.2);
148 :
149 : // Iron
150 0 : AliMaterial(10, "IRON1 ", 55.85, 26., 7.87, 1.76, 17.1);
151 0 : AliMaterial(30, "IRON2 ", 55.85, 26., 7.87, 1.76, 17.1);
152 0 : AliMaterial(50, "IRON3 ", 55.85, 26., 7.87, 1.76, 17.1);
153 :
154 : //
155 : // Copper
156 0 : AliMaterial(11, "COPPER1 ", 63.55, 29., 8.96, 1.43, 15.1);
157 0 : AliMaterial(31, "COPPER2 ", 63.55, 29., 8.96, 1.43, 15.1);
158 0 : AliMaterial(51, "COPPER3 ", 63.55, 29., 8.96, 1.43, 15.1);
159 :
160 : // Vacuum
161 0 : AliMixture(16, "VACUUM1 ", aAir, zAir, dAir1, 4, wAir);
162 0 : AliMixture(36, "VACUUM2 ", aAir, zAir, dAir1, 4, wAir);
163 0 : AliMixture(56, "VACUUM3 ", aAir, zAir, dAir1, 4, wAir);
164 :
165 : // Stainless Steel
166 0 : AliMixture(19, "STAINLESS STEEL1", asteel, zsteel, 7.88, 4, wsteel);
167 0 : AliMixture(39, "STAINLESS STEEL2", asteel, zsteel, 7.88, 4, wsteel);
168 0 : AliMixture(59, "STAINLESS STEEL3", asteel, zsteel, 7.88, 4, wsteel);
169 :
170 : // Cast iron
171 0 : AliMixture(18, "CAST IRON1", acasti, zcasti, 7.2, 4, wcasti);
172 0 : AliMixture(38, "CAST IRON2", acasti, zcasti, 7.2, 4, wcasti);
173 0 : AliMixture(58, "CAST IRON3", acasti, zcasti, 7.2, 4, wcasti);
174 : // ****************
175 : // Defines tracking media parameters.
176 : // Les valeurs sont commentees pour laisser le defaut
177 : // a GEANT (version 3-21, page CONS200), f.m.
178 : Float_t epsil, stmin, tmaxfd, deemax, stemax;
179 : epsil = .001; // Tracking precision,
180 : stemax = -1.; // Maximum displacement for multiple scat
181 : tmaxfd = -20.; // Maximum angle due to field deflection
182 : deemax = -.3; // Maximum fractional energy loss, DLS
183 : stmin = -.8;
184 : // ***************
185 :
186 : // Aluminum
187 0 : AliMedium(9, "ALU_C0 ", 9, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
188 0 : AliMedium(29, "ALU_C1 ", 29, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
189 0 : AliMedium(49, "ALU_C2 ", 49, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
190 :
191 : // Iron
192 0 : AliMedium(10, "FE_C0 ", 10, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
193 0 : AliMedium(30, "FE_C1 ", 30, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
194 0 : AliMedium(50, "FE_C2 ", 50, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
195 :
196 : // Copper
197 0 : AliMedium(11, "Cu_C0 ", 11, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
198 0 : AliMedium(31, "Cu_C1 ", 31, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
199 0 : AliMedium(51, "Cu_C2 ", 51, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
200 :
201 : // Vacuum
202 0 : AliMedium(16, "VA_C0 ", 16, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
203 0 : AliMedium(36, "VA_C1 ", 36, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
204 0 : AliMedium(56, "VA_C2 ", 56, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
205 :
206 : // Steel
207 0 : AliMedium(19, "ST_C0 ", 19, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
208 0 : AliMedium(39, "ST_C1 ", 39, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
209 0 : AliMedium(59, "ST_C3 ", 59, 0, fieldType, maxField, tmaxfd, stemax, deemax, epsil, stmin);
210 :
211 : //
212 : // Parameters for AD scintillator: NE-102 (BC400)
213 : //
214 : // NE-102, has the following properties : (from internet, need better reference)
215 : // Density : ca. 1.03 g/cm3
216 : // Electrons/cm3: 3.39 x 10^23
217 : // H atoms/cm3: 5.28 x 10^22
218 : // C atoms/cm3: 4.78 x 10^22
219 : // Ratio of H to C : 1.104 .
220 : // wavelength of emission : ~4.23 nm.
221 : // Decay time : ~2.4 ns.
222 : // Luminescent efficiency : typically 18% of NaI(Tl)
223 : // Photons/MeV: 2.5 x 10^4
224 : //
225 : // H // C
226 : // as[0] = 1.00794; as[1] = 12.011;
227 : // zs[0] = 1.; zs[1] = 6.;
228 : // ws[0] = 5.23; ws[1] = 4.74;
229 : // density = 1.032;
230 : // id = 1;
231 : // AliMixture( id, "NE102", as, zs, density, -2, ws );
232 : // AliMedium( id, "NE102", id, 1, fieldType, maxField, maxBending, maxStepSize,
233 : // maxEnergyLoss, precision, minStepSize );
234 :
235 : // ecalvovi@cern.ch
236 : // Parameters for AD scintillator: BC404
237 : //
238 : // NE-102, has the following properties : (from internet, need better reference)
239 : // Density : ca. 1.032 g/cm3
240 : // Electrons/cm3: 3.37 x 10^23
241 : // H atoms/cm3: 5.21 x 10^22
242 : // C atoms/cm3: 4.74 x 10^22
243 : // Ratio of H to C : 1.100
244 : // wavelength of emission : 408 nm.
245 : // Decay time : 1.8 ns.
246 : // Luminescent efficiency : typically 18% of NaI(Tl)
247 : // Photons/MeV: ??
248 : //
249 : // H // C
250 0 : as[0] = 1.00794; as[1] = 12.011;
251 0 : zs[0] = 1.; zs[1] = 6.;
252 0 : ws[0] = 5.21; ws[1] = 4.74;
253 : density = 1.032;
254 : id = 1;
255 0 : AliMixture( id, "BC404", as, zs, density, -2, ws );
256 0 : AliMedium ( id, "BC404", id, 1, fieldType, maxField, maxBending, maxStepSize,
257 : maxEnergyLoss, precision, minStepSize );
258 : // parameters AliMedium: numed name nmat isvol ifield fieldm tmaxfd stemax deemax epsil stmin
259 : // ...
260 : // isvol sensitive volume if isvol!=0
261 : // ifield magnetic field flag (see below)
262 : // fieldm maximum magnetic field
263 : // ...
264 : // ifield = 0 no magnetic field
265 : // = -1 user decision in guswim
266 : // = 1 tracking performed with Runge Kutta
267 : // = 2 tracking performed with helix
268 : // = 3 constant magnetic field along z
269 :
270 :
271 : //
272 : // Parameters for lightGuide:
273 : // TODO check material
274 : // Should be Poly(methyl methacrylate) (PMMA) acrylic
275 : // (C5O2H8)n
276 : // Density 1.18 g/cm3
277 : // Mixture PMMA Aeff=12.3994 Zeff=6.23653 rho=1.18 radlen=34.0677 intlen=63.3073
278 : // Element #0 : C Z= 6.00 A= 12.01 w= 0.600 natoms=5
279 : // Element #1 : H Z= 1.00 A= 1.01 w= 0.081 natoms=8
280 : // Element #2 : O Z= 8.00 A= 16.00 w= 0.320 natoms=2
281 :
282 : // Carbon Hydrogen Oxygen
283 0 : as[0] = 12.0107; as[1] = 1.00794; as[2] = 15.9994;
284 0 : zs[0] = 6.; zs[1] = 1.; zs[2] = 8.;
285 0 : ws[0] = 0.60; ws[1] = 0.081; ws[2] = 0.32;
286 : density = 1.18;
287 : id = 2;
288 0 : AliMixture( id, "PMMA", as, zs, density, 3, ws );
289 0 : AliMedium( id,"PMMA", id, 1, fieldType, maxField, maxBending, maxStepSize,
290 : maxEnergyLoss, precision, minStepSize );
291 :
292 :
293 : // mu-metal
294 : // Niquel Iron Molybdenum Manganese
295 0 : as[0] = 58.6934; as[1] = 55.845; as[2] = 95.94; as[3] = 54.9380;
296 0 : zs[0] = 28.; zs[1] = 26.; zs[2] = 42.; zs[3] = 25.;
297 0 : ws[0] = 0.802; ws[1] = 0.14079; ws[2] = 0.0485; ws[3] = 0.005;
298 : // Silicon Chromium Cobalt Aluminium
299 0 : as[4] = 28.0855; as[5] = 51.9961; as[6] = 58.9332; as[7] = 26.981539;
300 0 : zs[4] = 14.; zs[5] = 24.; zs[6] = 27.; zs[7] = 13.;
301 0 : ws[4] = 0.003; ws[5] = 0.0002; ws[6] = 0.0002; ws[7] = 0.0001;
302 : // Carbon Phosphorus Sulfur
303 0 : as[8] = 12.0107; as[9] = 30.97376; as[10] = 32.066;
304 0 : zs[8] = 6.; zs[9] = 15.; zs[10] = 16.;
305 0 : ws[8] = 0.00015; ws[9] = 0.00005; ws[10] = 0.00001;
306 : density = 8.25;
307 : id = 3;
308 0 : AliMixture( id, "MuMetal", as, zs, density, 11, ws );
309 0 : AliMedium( id,"MuMetal", id, 1, fieldType, maxField, maxBending, maxStepSize,
310 : maxEnergyLoss, precision, minStepSize );
311 :
312 : // Parameters for ADCPMA: Aluminium
313 : a_ad = 26.98;
314 : z_ad = 13.00;
315 : density = 2.7;
316 : radLength = 8.9;
317 : absLength = 37.2;
318 : id = 4;
319 0 : AliMaterial (id, "Alum", a_ad, z_ad, density, radLength, absLength, 0, 0 );
320 0 : AliMedium( id, "Alum", id, 1, fieldType, maxField, maxBending, maxStepSize,
321 : maxEnergyLoss, precision, minStepSize );
322 :
323 : // Parameters for ADCPMG: Glass for the simulation Aluminium
324 : // TODO fix material
325 : a_ad = 26.98;
326 : z_ad = 13.00;
327 : density = 2.7;
328 : radLength = 8.9;
329 : absLength = 37.2;
330 : id = 5;
331 0 : AliMaterial( id, "Glass", a_ad, z_ad, density, radLength, absLength, 0, 0 );
332 0 : AliMedium( id, "Glass", id, 1, fieldType, maxField, maxBending, maxStepSize,
333 : maxEnergyLoss, precision, minStepSize );
334 :
335 :
336 0 : }
337 : //_____________________________________________________________________________
338 : void AliAD::SetTreeAddress()
339 : {
340 : //
341 : // Sets tree address for hits.
342 : //
343 :
344 : TBranch *branch;
345 0 : char branchname[20];
346 0 : snprintf(branchname,19,"%s",GetName());
347 : // Branch address for hit tree
348 0 : TTree *treeH = fLoader->TreeH();
349 0 : if (treeH )
350 : {
351 0 : branch = treeH->GetBranch(branchname);
352 0 : if (branch) branch->SetAddress(&fHits);
353 : }
354 0 : }
355 :
356 :
357 : //_____________________________________________________________________________
358 : void AliAD::MakeBranch(Option_t* opt)
359 : {
360 0 : const char* oH = strstr(opt,"H");
361 0 : if (fLoader->TreeH() && oH && (fHits==0x0))
362 : {
363 0 : fHits = new TClonesArray("AliADhit",1000);
364 0 : fNhits = 0;
365 0 : }
366 0 : AliDetector::MakeBranch(opt);
367 0 : }
368 : //_____________________________________________________________________________
369 : AliLoader* AliAD::MakeLoader(const char* topfoldername)
370 : {
371 :
372 0 : AliDebug(1,Form("Creating AliADLoader, Top folder is %s ",topfoldername));
373 0 : fLoader = new AliADLoader(GetName(),topfoldername);
374 0 : return fLoader;
375 0 : }
376 :
377 : //_____________________________________________________________________________
378 : AliDigitizer* AliAD::CreateDigitizer(AliDigitizationInput* digInput) const
379 : {
380 : //
381 : // Creates a digitizer for AD
382 : //
383 0 : return new AliADDigitizer(digInput);
384 0 : }
385 : //_____________________________________________________________________________
386 : void AliAD::Hits2Digits(){
387 : //
388 : // Converts hits to digits
389 : //
390 : // Creates the AD digitizer
391 0 : AliADDigitizer* dig = new AliADDigitizer(this,AliADDigitizer::kHits2Digits);
392 :
393 : // Creates the digits
394 0 : dig->Digitize("");
395 :
396 : // deletes the digitizer
397 0 : delete dig;
398 0 : }
399 :
400 : //_____________________________________________________________________________
401 : void AliAD::Hits2SDigits(){
402 : //
403 : // Converts hits to summable digits
404 : //
405 : // Creates the AD digitizer
406 0 : AliADDigitizer* dig = new AliADDigitizer(this,AliADDigitizer::kHits2SDigits);
407 :
408 : // Creates the sdigits
409 0 : dig->Digitize("");
410 :
411 : // deletes the digitizer
412 0 : delete dig;
413 0 : }
414 :
415 : //_____________________________________________________________________________
416 :
417 : void AliAD::Digits2Raw()
418 : {
419 : //
420 : // Converts digits of the current event to raw data
421 : //
422 0 : GetCalibData();
423 :
424 0 : AliAD *fAD = (AliAD*)gAlice->GetDetector("AD");
425 0 : fLoader->LoadDigits();
426 0 : TTree* digits = fLoader->TreeD();
427 0 : if (!digits) {
428 0 : Error("Digits2Raw", "no digits tree");
429 0 : return;
430 : }
431 0 : TClonesArray * ADdigits = new TClonesArray("AliADdigit",1000);
432 0 : fAD->SetTreeAddress();
433 0 : digits->GetBranch("ADDigit")->SetAddress(&ADdigits);
434 :
435 0 : const char *fileName = AliDAQ::DdlFileName("AD",0);
436 0 : AliADBuffer* buffer = new AliADBuffer(fileName);
437 :
438 : // Get Trigger information first
439 : // Read trigger inputs from trigger-detector object
440 0 : AliDataLoader * dataLoader = fLoader->GetDigitsDataLoader();
441 0 : if( !dataLoader->IsFileOpen() )
442 0 : dataLoader->OpenFile("READ");
443 0 : AliTriggerDetector* trgdet = (AliTriggerDetector*)dataLoader->GetDirectory()->Get("Trigger");
444 : UInt_t triggerInfo = 0;
445 0 : if(trgdet) {
446 0 : triggerInfo = trgdet->GetMask() & 0xffff;
447 0 : }
448 : else {
449 0 : AliError(Form("There is no trigger object for %s",fLoader->GetName()));
450 : }
451 :
452 0 : buffer->WriteTriggerInfo((UInt_t)triggerInfo);
453 0 : buffer->WriteTriggerScalers();
454 0 : buffer->WriteBunchNumbers();
455 :
456 : // Now retrieve the channel information: charge smaples + time and
457 : // dump it into ADC and Time arrays
458 0 : Int_t nEntries = Int_t(digits->GetEntries());
459 0 : Short_t aADC[16][kADNClocks];
460 0 : Short_t aTime[16];
461 0 : Short_t aWidth[16];
462 0 : Bool_t aIntegrator[16];
463 0 : Bool_t aBBflag[16];
464 0 : Bool_t aBGflag[16];
465 :
466 0 : for (Int_t i = 0; i < nEntries; i++) {
467 0 : fAD->ResetDigits();
468 0 : digits->GetEvent(i);
469 0 : Int_t ndig = ADdigits->GetEntriesFast();
470 :
471 0 : if(ndig == 0) continue;
472 0 : for(Int_t k=0; k<ndig; k++){
473 0 : AliADdigit* fADDigit = (AliADdigit*) ADdigits->At(k);
474 :
475 0 : Int_t iChannel = kOnlineChannel[fADDigit->PMNumber()];
476 :
477 0 : for(Int_t iClock = 0; iClock < kADNClocks; ++iClock) aADC[iChannel][iClock] = fADDigit->ChargeADC(20-iClock);
478 0 : Int_t board = AliADCalibData::GetBoardNumber(iChannel);
479 0 : aTime[iChannel] = TMath::Nint(fADDigit->Time() / fCalibData->GetTimeResolution(board));
480 0 : aWidth[iChannel] = TMath::Nint(fADDigit->Width() / fCalibData->GetWidthResolution(board));
481 0 : aIntegrator[iChannel]= fADDigit->Integrator();
482 0 : aBBflag[iChannel] = fADDigit->GetBBflag();
483 0 : aBGflag[iChannel] = fADDigit->GetBGflag();
484 :
485 : //AliDebug(1,Form("DDL: %s\tdigit number: %d\tPM number: %d\tADC: %d\tTime: %f",
486 : //fileName,k,fADDigit->PMNumber(),aADC[iChannel][AliADdigit::kADNClocks/2],aTime[iChannel]));
487 : }
488 0 : }
489 :
490 : // Now fill raw data
491 : Int_t iCIU=0;
492 0 : for (Int_t iV0CIU = 0; iV0CIU < 8; iV0CIU++) {
493 :
494 0 : if(iV0CIU != 2 && iV0CIU != 5) {
495 0 : buffer->WriteEmptyCIU();
496 0 : continue;
497 : }
498 : // decoding of one Channel Interface Unit numbered iCIU - there are 8 channels per CIU (and 2 CIUs) :
499 0 : for(Int_t iChannel_Offset = iCIU*8; iChannel_Offset < (iCIU*8)+8; iChannel_Offset=iChannel_Offset+4) {
500 0 : for(Int_t iChannel = iChannel_Offset; iChannel < iChannel_Offset+4; iChannel++) {
501 0 : buffer->WriteChannel(iChannel, aADC[iChannel], aIntegrator[iChannel]);
502 : }
503 0 : buffer->WriteBeamFlags(&aBBflag[iChannel_Offset],&aBGflag[iChannel_Offset]);
504 0 : buffer->WriteMBInfo();
505 0 : buffer->WriteMBFlags();
506 0 : buffer->WriteBeamScalers();
507 : }
508 :
509 0 : for(Int_t iChannel = iCIU*8 + 7; iChannel >= iCIU*8; iChannel--) {
510 0 : buffer->WriteTiming(aTime[iChannel], aWidth[iChannel]);
511 : } // End of decoding of one CIU card
512 0 : iCIU++;
513 0 : } // end of decoding the eight CIUs
514 :
515 :
516 0 : delete buffer;
517 0 : fLoader->UnloadDigits();
518 0 : }
519 :
520 : //_____________________________________________________________________________
521 :
522 : Bool_t AliAD::Raw2SDigits(AliRawReader* rawReader)
523 : {
524 : // reads raw data to produce sdigits
525 : // for AD not implemented yet
526 0 : return kTRUE;
527 : }
528 :
529 : //_____________________________________________________________________________
530 : void AliAD::GetCalibData()
531 : {
532 : // Gets calibration object for AD set
533 : // Do nothing in case it is already loaded
534 0 : if (fCalibData) return;
535 :
536 0 : AliCDBEntry *entry = AliCDBManager::Instance()->Get("AD/Calib/Data");
537 0 : if (entry) fCalibData = (AliADCalibData*) entry->GetObject();
538 0 : if (!fCalibData) AliFatal("No calibration data from calibration database !");
539 :
540 : return;
541 0 : }
|