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