]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFLvHvDataPoints.cxx
Update in clusterFinderV1 algorithm
[u/mrichter/AliRoot.git] / TOF / AliTOFLvHvDataPoints.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 /*
17 $Log: AliTOFLvHvDataPoints.cxx,v $
18 */  
19
20 // AliTOFLvHvDataPoints class
21 // main aim to introduce the aliases for the TOF LV and HV DCS
22 // data points to be then
23 // stored in the OCDB, and to process them. 
24 // Process() method called by TOF preprocessor
25
26 #include "TString.h"
27 #include "TTimeStamp.h"
28 #include "TMap.h"
29 #include "TMath.h"
30 #include "TH1C.h"
31
32 #include "AliDCSValue.h"
33 #include "AliLog.h"
34 #include "AliBitPacking.h"
35
36 #include "AliTOFGeometry.h"
37 #include "AliTOFDCSmaps.h"
38 #include "AliTOFLvHvDataPoints.h"
39
40 class AliCDBMetaData;
41 class TDatime;
42
43 ClassImp(AliTOFLvHvDataPoints)
44
45 //---------------------------------------------------------------
46 AliTOFLvHvDataPoints::AliTOFLvHvDataPoints():
47   TObject(),
48   fRun(0),
49   fStartTime(0),
50   fEndTime(0),
51   fStartTimeDCSQuery(0),
52   fEndTimeDCSQuery(0),
53   fIsProcessed(kFALSE),
54   fFDR(kFALSE),
55   fNumberOfLVdataPoints(0),
56   fNumberOfHVdataPoints(0),
57   fNumberOfHVandLVmaps(0),
58   fStartingLVmap(0x0),
59   fStartingHVmap(0x0),
60   fHisto(0x0)
61 {
62   // main constructor 
63
64 }
65
66 //---------------------------------------------------------------
67 AliTOFLvHvDataPoints::AliTOFLvHvDataPoints(Int_t nRun, UInt_t startTime, UInt_t endTime, UInt_t startTimeDCSQuery, UInt_t endTimeDCSQuery):
68   TObject(),
69   fRun(nRun),
70   fStartTime(startTime),
71   fEndTime(endTime),
72   fStartTimeDCSQuery(startTimeDCSQuery),
73   fEndTimeDCSQuery(endTimeDCSQuery),
74   fIsProcessed(kFALSE),
75   fFDR(kFALSE),
76   fNumberOfLVdataPoints(0),
77   fNumberOfHVdataPoints(0),
78   fNumberOfHVandLVmaps(0),
79   fStartingLVmap(new AliTOFDCSmaps()),
80   fStartingHVmap(new AliTOFDCSmaps()),
81   fHisto(new TH1C("histo","",kNpads,-0.5,kNpads-0.5))
82 {
83
84   // constructor with arguments
85
86   AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s \n\tStartTime DCS Query %s \n\tEndTime DCS Query %s", nRun,
87                TTimeStamp(startTime).AsString(),
88                TTimeStamp(endTime).AsString(), 
89                TTimeStamp(startTimeDCSQuery).AsString(), 
90                TTimeStamp(endTimeDCSQuery).AsString()));
91
92   Init();
93
94 }
95
96 //---------------------------------------------------------------
97
98 AliTOFLvHvDataPoints::AliTOFLvHvDataPoints(const AliTOFLvHvDataPoints & data):
99   TObject(data), 
100   fRun(data.fRun),
101   fStartTime(data.fStartTime),
102   fEndTime(data.fEndTime),
103   fStartTimeDCSQuery(data.fStartTimeDCSQuery),
104   fEndTimeDCSQuery(data.fEndTimeDCSQuery),
105   fIsProcessed(data.fIsProcessed),
106   fFDR(data.fFDR),
107   fNumberOfLVdataPoints(data.fNumberOfLVdataPoints),
108   fNumberOfHVdataPoints(data.fNumberOfHVdataPoints),
109   fNumberOfHVandLVmaps(data.fNumberOfHVandLVmaps),
110   fStartingLVmap(data.fStartingLVmap),
111   fStartingHVmap(data.fStartingHVmap),
112   fHisto(data.fHisto)
113 {
114
115   // copy constructor
116
117   for(int i=0;i<kNddl;i++)
118     fAliasNamesXLVmap[i]=data.fAliasNamesXLVmap[i];
119     
120   for(int i=0;i<kNsectors;i++)
121     for(int j=0;j<kNplates;j++)
122       fAliasNamesXHVmap[i][j]=data.fAliasNamesXHVmap[i][j];
123  
124 }
125 //---------------------------------------------------------------
126
127 AliTOFLvHvDataPoints& AliTOFLvHvDataPoints:: operator=(const AliTOFLvHvDataPoints & data) { 
128
129   // assignment operator
130
131   if (this == &data)
132     return *this;
133
134   TObject::operator=(data);
135   fRun=data.GetRun();
136   fStartTime=data.GetStartTime();
137   fEndTime=data.GetEndTime();
138   fStartTimeDCSQuery=data.GetStartTimeDCSQuery();
139   fEndTimeDCSQuery=data.GetEndTimeDCSQuery();
140
141   fNumberOfLVdataPoints=data.fNumberOfLVdataPoints;
142   fNumberOfHVdataPoints=data.fNumberOfHVdataPoints;
143   fNumberOfHVandLVmaps=data.fNumberOfHVandLVmaps;
144
145   fStartingLVmap=data.fStartingLVmap;
146   fStartingHVmap=data.fStartingHVmap;
147
148   for(int i=0;i<kNddl;i++)
149     fAliasNamesXLVmap[i]=data.fAliasNamesXLVmap[i];
150
151   for(int i=0;i<kNsectors;i++)
152     for(int j=0;j<kNplates;j++)
153       fAliasNamesXHVmap[i][j]=data.fAliasNamesXHVmap[i][j];
154
155   fHisto=data.fHisto;
156
157   return *this;
158 }
159 //---------------------------------------------------------------
160 AliTOFLvHvDataPoints::~AliTOFLvHvDataPoints() {
161
162   // destructor
163
164   delete fStartingLVmap;
165   delete fStartingHVmap;
166
167 }
168
169 //---------------------------------------------------------------
170 Bool_t AliTOFLvHvDataPoints::ProcessData(TMap& aliasMap) {
171   //
172   // method to process the data
173   //
174
175   if(!(fAliasNamesXHVmap[0][0]) || !(fAliasNamesXLVmap[0])) Init();
176
177   AliInfo(Form(" Start Time = %i",fStartTime));
178   AliInfo(Form(" End Time = %i",fEndTime));
179   AliInfo(Form(" Start Time DCS Query= %i",fStartTimeDCSQuery));
180   AliInfo(Form(" End Time DCS Query= %i",fEndTimeDCSQuery));
181
182   if (fEndTime==fStartTime){
183     AliError(Form(" Run with null time length: start time = %d = end time = %d",fStartTime,fEndTime));
184     return kFALSE;
185   }
186
187
188   if (!ReadLVDataPoints(aliasMap)) return kFALSE;
189   AliDebug(1,Form(" Number of LV dp value changes = %d",fNumberOfLVdataPoints));
190
191   if (!ReadHVDataPoints(aliasMap)) return kFALSE;
192   AliDebug(1,Form(" Number of HV dp value changes = %d",fNumberOfHVdataPoints));
193
194   if (!MergeLVmap()) return kFALSE;
195
196   if (!MergeHVmap()) return kFALSE;
197
198   if (!MergeMaps()) return kFALSE;
199
200   fIsProcessed=kTRUE;
201
202   return kTRUE;
203
204 }
205
206 //---------------------------------------------------------------
207 Bool_t AliTOFLvHvDataPoints::MergeMaps() {
208   //
209   // Merge together LV and HV maps
210   //
211
212   Int_t timeMaps[kNmaxDataPoints];
213   for (Int_t ii=0; ii<kNmaxDataPoints; ii++) timeMaps[ii]=0;
214
215   for (Int_t ii=0; ii<fNumberOfHVdataPoints; ii++) {
216     timeMaps[fNumberOfHVandLVmaps++]=fHVDataPoints[ii]->GetTime();
217     //timeMaps[fNumberOfHVandLVmaps]=fHVDataPoints[ii]->GetTime();
218     //fNumberOfHVandLVmaps++;
219   }
220
221   Bool_t check = kTRUE;
222   for (Int_t jj=0; jj<fNumberOfLVdataPoints; jj++) {
223     check = kTRUE;
224     for (Int_t ii=0; ii<fNumberOfHVdataPoints; ii++)
225       check=check&&(fLVDataPoints[jj]->GetTime()!=timeMaps[ii]);
226
227     if (check) {
228       timeMaps[fNumberOfHVandLVmaps++]=fLVDataPoints[jj]->GetTime();
229       //timeMaps[fNumberOfHVandLVmaps]=fLVDataPoints[jj]->GetTime();
230       //fNumberOfHVandLVmaps++;
231     }
232   }
233
234   AliInfo(Form(" TOF HV dps = %d; TOF LV dps = %d; TOF HVandLV dps %d",
235                fNumberOfHVdataPoints, fNumberOfLVdataPoints, fNumberOfHVandLVmaps));
236
237
238   Int_t controller[kNmaxDataPoints];
239   for (Int_t ii=0; ii<kNmaxDataPoints; ii++) controller[ii]=-1;
240   TMath::Sort(fNumberOfHVandLVmaps,timeMaps,controller,kFALSE); // increasing order
241
242   //for (Int_t ii=0; ii<fNumberOfHVandLVmaps; ii++)
243   //AliDebug(3,Form(" time[%d] = %d", controller[ii],timeMaps[controller[ii]]));
244
245   Short_t array[kNpads];
246   for (Int_t iPad=0; iPad<kNpads; iPad++) array[iPad]=-1;
247   Int_t time = 0;
248
249   // HVandLV status map during run
250   for (Int_t index=0; index<fNumberOfHVandLVmaps; index++) {
251     time = timeMaps[controller[index]];
252     AliDebug(2,Form(" time[%d] = %d ", controller[index],time));
253
254     for (Int_t ii=0,jj=0; ii<fNumberOfHVdataPoints && jj<fNumberOfLVdataPoints; ii++,jj++) {
255
256       if ( (fHVDataPoints[ii]->GetTime()==time && fLVDataPoints[jj]->GetTime()<=time) ||
257            (fLVDataPoints[jj]->GetTime()==time && fHVDataPoints[ii]->GetTime()<=time) ) {
258
259         AliDebug(2,Form(" HVdp_time[%d]=%d, LVdp_time[%d]=%d",
260                         ii,fHVDataPoints[ii]->GetTime(),
261                         jj,fLVDataPoints[jj]->GetTime()));
262         for (Int_t iPad=0; iPad<kNpads; iPad++)
263           array[iPad] = fHVDataPoints[ii]->GetCellValue(iPad)*fLVDataPoints[jj]->GetCellValue(iPad);
264         AliTOFDCSmaps *object = new AliTOFDCSmaps(time,array);
265         fMap[index]= object;
266         break;
267
268       }
269       else if ( fHVDataPoints[ii]->GetTime()==time && fLVDataPoints[jj]->GetTime()>time ) {
270
271         AliDebug(2,Form(" HVdp_time[%d]=%d, (no LVdp)",ii,fHVDataPoints[ii]->GetTime()));
272         for (Int_t iPad=0; iPad<kNpads; iPad++)
273           array[iPad] = fHVDataPoints[ii]->GetCellValue(iPad);
274         AliTOFDCSmaps *object = new AliTOFDCSmaps(time,array);
275         fMap[index]= object;
276         break;
277
278       } else if ( fLVDataPoints[jj]->GetTime()==time && fHVDataPoints[ii]->GetTime()>time ) {
279
280         AliDebug(2,Form(" LVdp_time[%d]=%d, (no HVdp)",jj,fLVDataPoints[jj]->GetTime()));
281         for (Int_t iPad=0; iPad<kNpads; iPad++)
282           array[iPad] = fLVDataPoints[jj]->GetCellValue(iPad);
283         AliTOFDCSmaps *object = new AliTOFDCSmaps(time,array);
284         fMap[index]= object;
285         break;
286
287       }
288
289     }
290
291   }
292
293   return kTRUE;
294
295 }
296
297 //---------------------------------------------------------------
298 Bool_t AliTOFLvHvDataPoints::MergeHVmap() {
299   //
300   // Create HV maps from HV dps
301   //
302
303   Bool_t check= kFALSE;
304
305   if (fNumberOfHVdataPoints==0)
306     return check;
307   else {
308
309     // first map construction
310     for (Int_t iPad=0; iPad<kNpads; iPad++) {
311       if (fHVDataPoints[0]->GetCellValue(iPad)==-1)
312         fHVDataPoints[0]->SetCellValue(iPad,fStartingHVmap->GetCellValue(iPad));
313       //else
314       //fHVDataPoints[0]->SetCellValue(iPad,fHVDataPoints[0]->GetCellValue(iPad)*fStartingHVmap->GetCellValue(iPad));
315
316       if (iPad%(96*91)==0)
317         AliDebug(2,Form("HVdp0: channel=%6d -> %1d",iPad,fHVDataPoints[0]->GetCellValue(iPad)));
318     }
319
320     // other maps construction
321     for (Int_t ii=1; ii<fNumberOfHVdataPoints; ii++) {
322       for (Int_t iPad=0; iPad<kNpads; iPad++) {
323         if (fHVDataPoints[ii]->GetCellValue(iPad)==-1)
324           fHVDataPoints[ii]->SetCellValue(iPad,fHVDataPoints[ii-1]->GetCellValue(iPad));
325
326         if (iPad%(96*91)==0)
327           AliDebug(2,Form("HVdp%d: channel=%6d -> %1d",ii,iPad,fHVDataPoints[ii]->GetCellValue(iPad)));
328       }
329     }
330
331     check=kTRUE;
332   }
333
334   return kTRUE;
335
336 }
337
338 //---------------------------------------------------------------
339 Bool_t AliTOFLvHvDataPoints::MergeLVmap() {
340   //
341   // Create LV maps from LV dps
342   //
343
344   Bool_t check= kFALSE;
345
346   if (fNumberOfLVdataPoints==0)
347     return check;
348   else {
349
350     // first map construction
351     for (Int_t iPad=0; iPad<kNpads; iPad++) {
352       if (fLVDataPoints[0]->GetCellValue(iPad)==-1)
353         fLVDataPoints[0]->SetCellValue(iPad,fStartingLVmap->GetCellValue(iPad));
354       //else
355       //fLVDataPoints[0]->SetCellValue(iPad,fLVDataPoints[0]->GetCellValue(iPad)*fStartingLVmap->GetCellValue(iPad));
356
357       if (iPad%(96*91)==0)
358         AliDebug(2,Form("LVdp0: channel=%6d -> %1d",iPad,fLVDataPoints[0]->GetCellValue(iPad)));
359     }
360
361     // other maps construction
362     for (Int_t ii=1; ii<fNumberOfLVdataPoints; ii++) {
363       for (Int_t iPad=0; iPad<kNpads; iPad++) {
364         if (fLVDataPoints[ii]->GetCellValue(iPad)==-1)
365           fLVDataPoints[ii]->SetCellValue(iPad,fLVDataPoints[ii-1]->GetCellValue(iPad));
366
367         if (iPad%(96*91)==0)
368           AliDebug(2,Form("LVdp%d: channel=%6d -> %1d",ii,iPad,fLVDataPoints[ii]->GetCellValue(iPad)));
369       }
370     }
371
372     check=kTRUE;
373   }
374
375   return check;
376
377 }
378
379 //---------------------------------------------------------------
380 Bool_t AliTOFLvHvDataPoints::ReadHVDataPoints(TMap& aliasMap) {
381   //
382   // Read HV dps
383   //
384
385   TObjArray *aliasArr;
386   AliDCSValue* aValue;
387   AliDCSValue* aValuePrev;
388   Int_t val = 0;
389   Int_t time = 0;
390   Int_t nEntries = 0;
391
392   Short_t dummy[kNpads];
393
394   for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
395   // starting loop on aliases
396   for (int i=0; i<kNsectors; i++)
397     for (int j=0; j<kNplates; j++) {
398       aliasArr = (TObjArray*) aliasMap.GetValue(fAliasNamesXHVmap[i][j].Data());
399       if (!aliasArr) {
400         AliError(Form("Alias %s not found!", fAliasNamesXHVmap[i][j].Data()));
401         if (!fFDR)
402           return kFALSE;    // returning only in case we are not in a FDR run
403         else
404           continue;
405       }
406
407       nEntries = aliasArr->GetEntries();
408     
409       if (nEntries<2) AliDebug(1, Form(" NB: number of values for dp %s %d (less than 2)",fAliasNamesXHVmap[i][j].Data(),nEntries));
410
411       if (nEntries==0) {
412         AliError(Form("Alias %s has no entries! Nothing will be stored",
413                       fAliasNamesXHVmap[i][j].Data()));
414         continue;
415       }
416       else {
417
418         // read the first value
419         for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
420         aValue = (AliDCSValue*) aliasArr->At(0);
421         val = aValue->GetInt();
422         time = aValue->GetTimeStamp();
423         AliDebug(1, Form(" at t=%d, 1st value for dp %s = %d", time, fAliasNamesXHVmap[i][j].Data(), val));
424         FillHVarrayPerDataPoint(i,j,val,dummy);
425         AliTOFDCSmaps *object0 = new AliTOFDCSmaps(time,dummy);
426         if (InsertHVDataPoint(object0)!=0) return kTRUE; // to be verified
427
428         // read the values from the second one
429         for (Int_t iEntry=1; iEntry<nEntries; iEntry++) {
430           for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
431           aValue = (AliDCSValue*) aliasArr->At(iEntry);
432           val = aValue->GetInt();
433           time = aValue->GetTimeStamp();
434           AliDebug(1, Form(" at t=%d, %dth value for dp %s = %d", time, iEntry+1, fAliasNamesXHVmap[i][j].Data(), val));
435           aValuePrev = (AliDCSValue*) aliasArr->At(iEntry-1);
436           if (aValuePrev->GetInt()!=val) {
437             FillHVarrayPerDataPoint(i,j,val,dummy);
438             AliTOFDCSmaps *object = new AliTOFDCSmaps(time,dummy);
439             if (InsertHVDataPoint(object)!=0) return kTRUE; // to be verified
440           }
441
442         }
443       }
444     }
445
446   if (fNumberOfHVdataPoints==0) {
447     AliInfo("Valid HV dps not found. By default all HV TOF channels (except the ones in the PHOS holes) switched ON, at SOR.");
448     if (InsertHVDataPoint(fStartingHVmap)!=0) return kTRUE; // to be verified
449     /*
450     AliTOFDCSmaps *object = new AliTOFDCSmaps((Int_t)(TMath::Power(2,31)-1),fStartingHVmap->GetArray());
451     if (InsertHVDataPoint(object)!=0) return kTRUE; // to be verified
452     */
453   }
454   /*
455   else if (fNumberOfHVdataPoints==1) {
456     AliInfo(Form("Only one valid HV dp found (time=%d). This dp will be replicated at EOR (%d).",fHVDataPoints[0]->GetTime(),fEndTimeDCSQuery));
457     AliTOFDCSmaps *object = new AliTOFDCSmaps((Int_t)(TMath::Power(2,31)-1),fHVDataPoints[0]->GetArray());
458     if (InsertHVDataPoint(object)!=0) return kTRUE; // to be verified
459   }
460   */
461   return kTRUE;
462
463 }
464
465 //---------------------------------------------------------------
466 Bool_t AliTOFLvHvDataPoints::ReadLVDataPoints(TMap& aliasMap) {
467   //
468   // Read LV dps
469   //
470
471   TObjArray *aliasArr;
472   AliDCSValue* aValue;
473   AliDCSValue* aValuePrev;
474   Int_t val = 0;
475   Int_t time = 0;
476   Int_t nEntries = 0;
477
478   Short_t dummy[kNpads];
479
480   for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
481   // starting loop on aliases
482   for (int i=0; i<kNddl; i++) {
483     aliasArr = (TObjArray*) aliasMap.GetValue(fAliasNamesXLVmap[i].Data());
484     if (!aliasArr) {
485       AliError(Form("Alias %s not found!", fAliasNamesXLVmap[i].Data()));
486       if (!fFDR)
487         return kFALSE;    // returning only in case we are not in a FDR run
488       else
489         continue;
490     }
491
492     nEntries = aliasArr->GetEntries();
493     
494     if (nEntries<2) AliDebug(1, Form(" NB: number of values for dp %s %d (less than 2)",fAliasNamesXLVmap[i].Data(),nEntries));
495     
496     if (nEntries==0) {
497       AliError(Form("Alias %s has no entries! Nothing will be stored",
498                     fAliasNamesXLVmap[i].Data()));
499       continue;
500     }
501     else {
502
503       // read the first value
504       for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
505       aValue = (AliDCSValue*) aliasArr->At(0);
506       val = aValue->GetInt();
507       time = aValue->GetTimeStamp();
508       AliDebug(1, Form(" at t=%d,1st value for dp %s = %d", time, fAliasNamesXLVmap[i].Data(), val));
509       FillLVarrayPerDataPoint(i,val,dummy);
510       AliTOFDCSmaps *object0 = new AliTOFDCSmaps(time,dummy);
511       if (InsertLVDataPoint(object0)!=0) return kTRUE; // to be verified
512
513       // read the values from the second one
514       for (Int_t iEntry=1; iEntry<nEntries; iEntry++) {
515         for (Int_t iBin=0; iBin<kNpads; iBin++) dummy[iBin]=-1;
516         aValue = (AliDCSValue*) aliasArr->At(iEntry);
517         val = aValue->GetInt();
518         time = aValue->GetTimeStamp();
519         AliDebug(1, Form(" at t=%d,%dth value for dp %s = %d", time, iEntry+1, fAliasNamesXLVmap[i].Data(), val));
520         aValuePrev = (AliDCSValue*) aliasArr->At(iEntry-1);
521         if (aValuePrev->GetInt()!=val) {
522           FillLVarrayPerDataPoint(i,val,dummy);
523           AliTOFDCSmaps *object = new AliTOFDCSmaps(time,dummy);
524           if (InsertLVDataPoint(object)!=0) return kTRUE; // to be verified
525         }
526
527       }
528     }
529   }
530
531   if (fNumberOfLVdataPoints==0) {
532     AliInfo("Valid LV dps not found. By default all LV TOF channels switched ON, at SOR.");
533     if (InsertLVDataPoint(fStartingLVmap)!=0) return kTRUE; // to be verified
534     /*
535     AliTOFDCSmaps *object = new AliTOFDCSmaps((Int_t)(TMath::Power(2,31)-1),fStartingHVmap->GetArray());
536     if (InsertLVDataPoint(object)!=0) return kTRUE; // to be verified
537     */
538   }
539   /*
540   else if (fNumberOfLVdataPoints==1) {
541     AliInfo(Form("Only one valid LV dp found (time=%d). This dp will be replicated at EOR (%d).",fLVDataPoints[0]->GetTime(),fEndTimeDCSQuery));
542     AliTOFDCSmaps *object = new AliTOFDCSmaps((Int_t)(TMath::Power(2,31)-1),fLVDataPoints[0]->GetArray());
543     if (InsertLVDataPoint(object)!=0) return kTRUE; // to be verified
544   }
545   */
546   return kTRUE;
547
548 }
549
550 //---------------------------------------------------------------
551 Int_t AliTOFLvHvDataPoints::InsertHVDataPoint(AliTOFDCSmaps *object)
552 {
553   //
554   // Insert HV dp in the HV dps array.
555   // The HV dps array is sorted according to increasing dp timeStamp value
556   //
557
558   if (fNumberOfHVdataPoints==kNmaxDataPoints) {
559     AliError("Too many HV data points!");
560     return 1;
561   }
562
563   if (fNumberOfHVdataPoints==0) {
564     fHVDataPoints[fNumberOfHVdataPoints++] = object;
565     return 0;
566   }
567
568   for (Int_t index=0; index<fNumberOfHVdataPoints; index++) {
569     if (object->GetTime()==fHVDataPoints[index]->GetTime()) {
570       fHVDataPoints[index]->Update(object);
571       return 0;
572     }
573   }
574
575   Int_t ii = FindHVdpIndex(object->GetTime());
576   memmove(fHVDataPoints+ii+1 ,fHVDataPoints+ii,(fNumberOfHVdataPoints-ii)*sizeof(AliTOFDCSmaps*));
577   fHVDataPoints[ii] = object;
578   fNumberOfHVdataPoints++;
579   
580   return 0;
581
582 }
583
584 //_________________________________________________________________________
585 Int_t AliTOFLvHvDataPoints::FindHVdpIndex(Int_t z) const {
586   //
587   // This function returns the index of the nearest HV DP in time
588   //
589
590   if (fNumberOfHVdataPoints==0) return 0;
591   if (z <= fHVDataPoints[0]->GetTime()) return 0;
592   if (z > fHVDataPoints[fNumberOfHVdataPoints-1]->GetTime()) return fNumberOfHVdataPoints;
593   Int_t b = 0, e = fNumberOfHVdataPoints-1, m = (b+e)/2;
594   for (; b<e; m=(b+e)/2) {
595     if (z > fHVDataPoints[m]->GetTime()) b=m+1;
596     else e=m;
597   }
598
599   return m;
600
601 }
602
603 //---------------------------------------------------------------
604 Int_t AliTOFLvHvDataPoints::InsertLVDataPoint(AliTOFDCSmaps *object)
605 {
606   //
607   // Insert LV dp in the LV dps array.
608   // The LV dps array is sorted according to increasing dp timeStamp value
609   //
610
611   if (fNumberOfLVdataPoints==kNmaxDataPoints) {
612     AliError("Too many LV data points!");
613     return 1;
614   }
615
616   if (fNumberOfLVdataPoints==0) {
617     fLVDataPoints[fNumberOfLVdataPoints++] = object;
618     return 0;
619   }
620
621   for (Int_t index=0; index<fNumberOfLVdataPoints; index++) {
622     if (object->GetTime()==fLVDataPoints[index]->GetTime()) {
623       fLVDataPoints[index]->Update(object);
624       return 0;
625     }
626   }
627
628   Int_t ii = FindLVdpIndex(object->GetTime());
629   memmove(fLVDataPoints+ii+1 ,fLVDataPoints+ii,(fNumberOfLVdataPoints-ii)*sizeof(AliTOFDCSmaps*));
630   fLVDataPoints[ii] = object;
631   fNumberOfLVdataPoints++;
632   
633   return 0;
634
635 }
636
637 //_________________________________________________________________________
638 Int_t AliTOFLvHvDataPoints::FindLVdpIndex(Int_t z) const {
639   //
640   // This function returns the index of the nearest LV DP in time
641   //
642
643   if (fNumberOfLVdataPoints==0) return 0;
644   if (z <= fLVDataPoints[0]->GetTime()) return 0;
645   if (z > fLVDataPoints[fNumberOfLVdataPoints-1]->GetTime()) return fNumberOfLVdataPoints;
646   Int_t b = 0, e = fNumberOfLVdataPoints-1, m = (b+e)/2;
647   for (; b<e; m=(b+e)/2) {
648     if (z > fLVDataPoints[m]->GetTime()) b=m+1;
649     else e=m;
650   }
651
652   return m;
653
654 }
655
656 //---------------------------------------------------------------
657 void AliTOFLvHvDataPoints::Init(){
658   //
659   // Initialize aliases and DCS data
660   //
661
662   TString sindex;
663   for(int i=0;i<kNsectors;i++)
664     for(int j=0;j<kNplates;j++) {
665       fAliasNamesXHVmap[i][j] = "TOF_HVSTATUS_";
666       sindex.Form("SM%02dMOD%1d",i,j);
667       fAliasNamesXHVmap[i][j] += sindex;
668     }
669
670
671   for(int i=0;i<kNddl;i++) {
672     fAliasNamesXLVmap[i] = "TOF_FEACSTATUS_";
673     sindex.Form("%02d",i);
674     fAliasNamesXLVmap[i] += sindex;
675   }
676
677   fStartingLVmap->SetTime(0);
678   for (Int_t iPad=0; iPad<kNpads; iPad++)
679     fStartingLVmap->SetCellValue(iPad,1);
680
681   fStartingHVmap->SetTime(0);
682   for (Int_t iPad=0; iPad<kNpads; iPad++) {
683     Int_t vol[5] = {-1,-1,-1,-1,-1};
684     AliTOFGeometry::GetVolumeIndices(iPad,vol);
685     if ( (vol[0]==13 || vol[0]==14 || vol[0]==15) &&
686          vol[1]==2)
687       fStartingHVmap->SetCellValue(iPad,0);
688     else
689       fStartingHVmap->SetCellValue(iPad,1);
690   }
691
692 }
693
694 //---------------------------------------------------------------
695 void AliTOFLvHvDataPoints::FillHVarrayPerDataPoint(Int_t sector, Int_t plate, UInt_t baseWord, Short_t *array) const
696 {
697   //
698   // Set the status of the TOF pads connected to the HV dp
699   // labelled by sector and plate numbers
700   //
701
702   Int_t det[5] = {sector, plate, -1, -1, -1};
703   UInt_t checkBit = 0;
704
705   Int_t channel = -1;
706
707   for (Int_t iStrip=0; iStrip<AliTOFGeometry::NStrip(plate); iStrip++) {
708     checkBit = AliBitPacking::UnpackWord(baseWord,iStrip,iStrip);
709     for (Int_t iPadZ=0; iPadZ<AliTOFGeometry::NpadZ(); iPadZ++)
710       for (Int_t iPadX=0; iPadX<AliTOFGeometry::NpadX(); iPadX++) {
711         det[2] = iStrip;
712         det[3] = iPadZ;
713         det[4] = iPadX;
714         channel = AliTOFGeometry::GetIndex(det);
715         array[channel]=checkBit;
716     }
717   }
718
719
720 }
721
722 //---------------------------------------------------------------
723 void AliTOFLvHvDataPoints::FillLVarrayPerDataPoint(Int_t nDDL, UInt_t baseWord, Short_t *array) const 
724 {
725   //
726   // Set the status of the TOF pads connected to the LV dp
727   // labelled by TOF crate number
728   //
729
730   Int_t det[5] = {nDDL/4, -1, -1, -1, -1};
731   UInt_t checkBit = 0;
732
733   Int_t iStripXsm[6] = {-1,-1,-1,-1,-1,-1};
734   Int_t firstPadX = -1;
735   Int_t lastPadX = -1;
736   Int_t plate = -1;
737   Int_t strip = -1;
738
739   Int_t channel = -1;
740
741   for (Int_t nFEAC=0; nFEAC<8; nFEAC++) {
742     checkBit = AliBitPacking::UnpackWord(baseWord,nFEAC,nFEAC);
743     firstPadX = -1;
744     lastPadX = -1;
745     GetStripsConnectedToFEAC(nDDL, nFEAC, iStripXsm, firstPadX,lastPadX);
746     for (Int_t index=0; index<6; index++) {
747       if (iStripXsm[index]==-1) continue;
748
749       for (Int_t iPadZ=0; iPadZ<AliTOFGeometry::NpadZ(); iPadZ++)
750         for (Int_t iPadX=firstPadX; iPadX<=lastPadX; iPadX++) {
751           AliTOFGeometry::GetStripAndModule(iStripXsm[index],plate,strip);
752           det[1] = plate;
753           det[2] = strip;
754           det[3] = iPadZ;
755           det[4] = iPadX;
756           channel = AliTOFGeometry::GetIndex(det);
757           array[channel]=checkBit;
758         }
759     }
760   }
761
762
763 }
764
765 //---------------------------------------------------------------
766 void AliTOFLvHvDataPoints::GetStripsConnectedToFEAC(Int_t nDDL, Int_t nFEAC, Int_t *iStrip, Int_t &firstPadX, Int_t &lastPadX) const
767 {
768   //
769   // FEAC-strip mapping:
770   // return the strips and first PadX numbers
771   // connected to the FEAC number nFEAC in the crate number nDDL
772   //
773
774   for (Int_t ii=0; ii<6; ii++) iStrip[ii]=-1;
775
776   if (nDDL<0 || nDDL>=kNddl || nFEAC<0 || nFEAC>=8) return;
777
778   switch (nDDL%4) {
779   case 0:
780     firstPadX =  0;
781     lastPadX  = AliTOFGeometry::NpadX()/2-1;
782
783     if (nFEAC<=2)
784       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC;
785     else if (nFEAC==3)
786       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=ii+6*nFEAC;
787     else if (nFEAC==4)
788       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC-1;
789     else if (nFEAC==5)
790       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=ii+6*nFEAC-1;
791     else if (nFEAC==6)
792       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC-2;
793     else if (nFEAC==7)
794       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=ii+6*nFEAC-2;
795
796     break; 
797   case 1:
798     firstPadX = AliTOFGeometry::NpadX()/2;
799     lastPadX  = AliTOFGeometry::NpadX()-1;
800
801     if (nFEAC<=2)
802       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC;
803     else if (nFEAC==3)
804       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=ii+6*nFEAC;
805     else if (nFEAC==4)
806       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC-1;
807     else if (nFEAC==5)
808       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC-1;
809     else if (nFEAC==6)
810       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=ii+6*nFEAC-1;
811     else if (nFEAC==7)
812       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=ii+6*nFEAC-2;
813
814     break;
815   case 2:
816     firstPadX = AliTOFGeometry::NpadX()/2;
817     lastPadX  = AliTOFGeometry::NpadX()-1;
818
819     if (nFEAC<=2)
820       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC);
821     else if (nFEAC==3)
822       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=90-(ii+6*nFEAC);
823     else if (nFEAC==4)
824       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC-1);
825     else if (nFEAC==5)
826       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=90-(ii+6*nFEAC-1);
827     else if (nFEAC==6)
828       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC-2);
829     else if (nFEAC==7)
830       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=90-(ii+6*nFEAC-2);
831
832     break;
833   case 3:
834     firstPadX =  0;
835     lastPadX  = AliTOFGeometry::NpadX()/2-1;
836
837     if (nFEAC<=2)
838       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC);
839     else if (nFEAC==3)
840       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=90-(ii+6*nFEAC);
841     else if (nFEAC==4)
842       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC-1);
843     else if (nFEAC==5)
844       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC-1);
845     else if (nFEAC==6)
846       for (Int_t ii=0; ii<5; ii++) iStrip[ii]=90-(ii+6*nFEAC-1);
847     else if (nFEAC==7)
848       for (Int_t ii=0; ii<6; ii++) iStrip[ii]=90-(ii+6*nFEAC-2);
849
850     break;
851  }
852
853 }
854
855 //---------------------------------------------------------------
856 void AliTOFLvHvDataPoints::Draw(const Option_t* /*option*/)
857 {
858   //
859   // Draw all histos and graphs
860   //
861
862   if(!fIsProcessed) return;
863   /*
864   TCanvas *ch;
865   TString canvasHistoName="Histos";
866   ch = new TCanvas(canvasHistoName,canvasHistoName,20,20,600,600);
867   ch->cd();
868   */
869   // to be implemented
870
871 }
872
873
874 //---------------------------------------------------------------
875 void AliTOFLvHvDataPoints::DrawHVandLVMap(Int_t index)
876 {
877   //
878   // Draw HV+LV map labelled as index
879   //
880
881   if(!fIsProcessed) return;
882
883   if (index>=fNumberOfHVandLVmaps) return;
884
885   AliTOFDCSmaps *mappa=(AliTOFDCSmaps*)GetHVandLVmap(index);
886
887   char title[100];
888   if (index==0) sprintf(title,"HVandLV map at time %d (%dst map)",mappa->GetTime(),index+1);
889   else if (index==1) sprintf(title,"HVandLV map at time %d (%dnd map)",mappa->GetTime(),index+1);
890   else if (index==2) sprintf(title,"HVandLV map at time %d (%drd map)",mappa->GetTime(),index+1);
891   else if (index>=3) sprintf(title,"HVandLV map at time %d (%dth map)",mappa->GetTime(),index+1);
892   fHisto->Delete();
893   fHisto = new TH1C("histo","",kNpads,-0.5,kNpads-0.5);
894   //fHisto->Clear();
895   fHisto->SetTitle(title);
896
897   for (Int_t ii=0; ii<kNpads; ii++)
898     fHisto->SetBinContent(ii+1,mappa->GetCellValue(ii));
899
900   fHisto->Draw();
901
902 }
903
904 //---------------------------------------------------------------
905 void AliTOFLvHvDataPoints::DrawLVMap(Int_t index)
906 {
907   //
908   // Draw LV map labelled as index
909   //
910
911   if(!fIsProcessed) return;
912
913   if (index>=fNumberOfLVdataPoints) return;
914
915   AliTOFDCSmaps *mappa=(AliTOFDCSmaps*)GetLVmap(index);
916
917   char title[100];
918   if (index==0) sprintf(title,"LV map at time %d (%dst map)",mappa->GetTime(),index+1);
919   else if (index==1) sprintf(title,"LV map at time %d (%dnd map)",mappa->GetTime(),index+1);
920   else if (index==2) sprintf(title,"LV map at time %d (%drd map)",mappa->GetTime(),index+1);
921   else if (index>=3) sprintf(title,"LV map at time %d (%dth map)",mappa->GetTime(),index+1);
922   fHisto->Delete();
923   fHisto = new TH1C("histo","",kNpads,-0.5,kNpads-0.5);
924   //fHisto->Clear();
925   fHisto->SetTitle(title);
926
927   for (Int_t ii=0; ii<kNpads; ii++)
928     fHisto->SetBinContent(ii+1,mappa->GetCellValue(ii));
929
930   fHisto->Draw();
931
932 }
933
934 //---------------------------------------------------------------
935 void AliTOFLvHvDataPoints::DrawHVMap(Int_t index)
936 {
937   //
938   // Draw HV map labelled as index
939   //
940
941   if(!fIsProcessed) return;
942
943   if (index>=fNumberOfHVdataPoints) return;
944
945   AliTOFDCSmaps *mappa=(AliTOFDCSmaps*)GetHVmap(index);
946
947   char title[100];
948   if (index==0) sprintf(title,"HV map at time %d (%dst map)",mappa->GetTime(),index+1);
949   else if (index==1) sprintf(title,"HV map at time %d (%dnd map)",mappa->GetTime(),index+1);
950   else if (index==2) sprintf(title,"HV map at time %d (%drd map)",mappa->GetTime(),index+1);
951   else if (index>=3) sprintf(title,"HV map at time %d (%dth map)",mappa->GetTime(),index+1);
952   fHisto->Delete();
953   fHisto = new TH1C("histo","",kNpads,-0.5,kNpads-0.5);
954   //fHisto->Clear();
955   fHisto->SetTitle(title);
956
957   for (Int_t ii=0; ii<kNpads; ii++)
958     fHisto->SetBinContent(ii+1,mappa->GetCellValue(ii));
959
960   fHisto->Draw();
961
962 }
963
964 //---------------------------------------------------------------
965 AliTOFDCSmaps *AliTOFLvHvDataPoints::GetHVandLVmapAtEOR()
966 {
967   //
968   // Returns HVandLV status map at EOR.
969   // If the end-of-run has been caused by TOF self,
970   // the last but two value of HVandLV status map
971   // will be taken into account.
972   // This last condition is true
973   // if the time interval between the second-last DP and the last one
974   // is less than 60s.
975   //
976
977   AliTOFDCSmaps * lvANDhvMap = 0;
978
979   if (fNumberOfHVandLVmaps<=2) { // nothing changed during the run
980
981     lvANDhvMap = fMap[fNumberOfHVandLVmaps-1];
982
983   }
984   else {
985
986     if (fMap[fNumberOfHVandLVmaps-1]->GetTime()-fMap[fNumberOfHVandLVmaps-2]->GetTime()<=60)
987       lvANDhvMap = (AliTOFDCSmaps*)fMap[fNumberOfHVandLVmaps-3];
988     else
989       lvANDhvMap = (AliTOFDCSmaps*)fMap[fNumberOfHVandLVmaps-1];
990
991   }
992
993   return lvANDhvMap;
994
995 }