]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDPreprocessor.cxx
Split DAQ DA into two tasks
[u/mrichter/AliRoot.git] / TRD / AliTRDPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 // This class is a first implementation for the TRD.                      //
21 // It takes data from HLT and computes the parameters                     //
22 // and stores both reference data and online calibration                  //
23 // parameters in the CDB                                                  //
24 // It alsotakes DCS data, does spline fits                                //
25 // and stores both reference data and spline fits results                 //
26 // in the CDB                                                             //
27 //                                                                        //
28 // Author:                                                                //
29 //   R. Bailhache (R.Bailhache@gsi.de)                                    //
30 //   W. Monange   (w.monange@gsi.de)                                      //
31 //                                                                        //
32 ////////////////////////////////////////////////////////////////////////////
33
34 #include "AliTRDPreprocessor.h"
35
36 #include <TFile.h>
37 #include <TProfile2D.h>
38 #include <TStopwatch.h>
39 #include <TObjString.h>
40 #include <TString.h>
41 #include <TList.h>
42 #include <TCollection.h>
43
44 #include "AliCDBMetaData.h"
45 #include "AliLog.h"
46
47 #include "AliTRDSensorArray.h"
48 #include "AliTRDCalibraFit.h"
49 #include "AliTRDCalibraMode.h"
50 #include "AliTRDCalibPadStatus.h"
51 #include "Cal/AliTRDCalDet.h"
52 #include "Cal/AliTRDCalPadStatus.h"
53
54 ClassImp(AliTRDPreprocessor)
55
56 //______________________________________________________________________________________________
57 AliTRDPreprocessor::AliTRDPreprocessor(AliShuttleInterface *shuttle)
58   :AliPreprocessor("TRD", shuttle),
59   fResult(0),
60   fVdriftHLT(0)
61 {
62   //
63   // Constructor
64   //
65
66 }
67
68 //______________________________________________________________________________________________
69 AliTRDPreprocessor::~AliTRDPreprocessor()
70 {
71   //
72   // Destructor
73   //
74
75 }
76
77 //______________________________________________________________________________________________
78 void AliTRDPreprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
79 {
80   //
81   // Initialization routine for the TRD preprocessor
82   //
83
84   AliPreprocessor::Initialize(run,startTime,endTime);
85
86 }
87
88 //______________________________________________________________________________________________
89 UInt_t AliTRDPreprocessor::Process(TMap* dcsAliasMap)
90 {
91   //
92   // Process DCS and calibration part for HLT
93   //
94
95   fResult = 0;
96
97   //
98   // DCS
99   //
100   
101   ProcessDCS(dcsAliasMap);
102
103   //
104   // Process the calibration data for the HLT and DAQ part
105   //
106
107   TString runType = GetRunType();
108   printf("runtype %s\n",(const char*)runType);
109
110
111   if (strcmp(runType, "PEDESTAL") == 0){
112     ExtractPedestals();
113   }
114   
115   if (strcmp(runType, "PHYSICS") == 0){
116     printf("test\n");
117     //if(GetRunParameter("HLTStatus")==1) ExtractHLT();
118     ExtractHLT();
119     ExtractDriftVelocityDAQ();
120   }
121   printf("result is %d\n",fResult);  
122   return fResult;  
123   
124 }
125
126
127 //______________________________________________________________________________
128 void AliTRDPreprocessor::ProcessDCS(TMap * dcsAliasMap)
129 {
130
131   AliCDBMetaData metaData;
132   metaData.SetBeamPeriod(0);
133   metaData.SetResponsible("Wilfried Monange/Raphaelle Bailhache");
134   metaData.SetComment("TRD calib test");
135         
136         
137   Log ("****** DCS ******\n");
138         
139   TObjArray * list = AliTRDSensorArray::GetList ();
140         
141   if (list == 0x0) {
142     Log ("Error during AliTRDSensorArray::GetList");
143     Log ("DCS will not be processing");
144     return;
145   }
146
147   Int_t nEntries = list->GetEntries ();
148   Log (Form ("%d alias loaded", nEntries));
149                 
150   Bool_t * results=new Bool_t [nEntries];
151   Int_t  * nGraph=new Int_t [nEntries];
152                 
153   for (Int_t iAlias = 0; iAlias < nEntries; iAlias++) {
154                         
155     AliTRDSensorArray * oneTRDDCS = (AliTRDSensorArray *)list->At (iAlias);
156                         
157     oneTRDDCS->SetStartTime (TTimeStamp (fStartTime));
158     oneTRDDCS->SetEndTime (TTimeStamp (fEndTime));
159                         
160     Log (Form("Processing DCS : \"%s\"", oneTRDDCS->GetStoreName ().Data ()));
161                         
162     TMap * map;
163
164     map=oneTRDDCS->ExtractDCS (dcsAliasMap);
165                 
166     nGraph [iAlias] = map->GetEntries ();
167                 
168     if (nGraph [iAlias] == 0) {
169       Log("No TGraph for this dcsDatapointAlias : not stored");
170       results [iAlias] = kFALSE;
171       continue;
172     }
173                 
174     oneTRDDCS->SetGraph(map);
175     results[iAlias]=Store("Calib", oneTRDDCS->GetStoreName().Data(), oneTRDDCS, &metaData, 0, kTRUE); 
176     delete map;         
177
178     //results [iAlias] = StoreReferenceData("Calib", oneTRDDCS->GetStoreName ().Data (), oneTRDDCS, &metaData); 
179
180
181     if (!results[iAlias]) {
182       AliError("Problem during StoreRef DCS"); 
183       fResult|=kEStore;
184     }
185
186
187     //BEGIN TEST (should not be removed ...)
188     /*
189     oneTRDDCS->ClearGraph();
190     oneTRDDCS->ClearFit();
191     oneTRDDCS->SetDiffCut2 (0.1);
192     map=oneTRDDCS->ExtractDCS (dcsAliasMap);
193     oneTRDDCS->SetGraph (map);
194     Store("Calib", ("cut_"+oneTRDDCS->GetStoreName()).Data(), oneTRDDCS, &metaData, 0, kTRUE); 
195     delete map;
196
197
198     if(iAlias==1 || iAlias==19) continue;
199     
200     oneTRDDCS->ClearGraph();
201     oneTRDDCS->ClearFit();
202     oneTRDDCS->SetDiffCut2(0);
203     map=oneTRDDCS->ExtractDCS(dcsAliasMap);
204     oneTRDDCS->MakeSplineFit(map);
205     Store("Calib", ("fit_"+oneTRDDCS->GetStoreName()).Data() , oneTRDDCS, &metaData, 0, kTRUE); 
206     delete map;
207
208      
209     oneTRDDCS->ClearGraph(); 
210     oneTRDDCS->ClearFit();
211     oneTRDDCS->SetDiffCut2 (0.1);
212     map=oneTRDDCS->ExtractDCS (dcsAliasMap);
213     oneTRDDCS->MakeSplineFit(map);
214     Store("Calib", ("cutfit_"+oneTRDDCS->GetStoreName()).Data() , oneTRDDCS, &metaData, 0, kTRUE); 
215     delete map;
216     */    
217     //END TEST
218
219       
220       
221   }
222                 
223   Log ("         Summury of DCS :\n");
224   Log (Form("%30s %10s %10s", "dcsDatapointAlias", "Stored ?", "# graph"));
225   for (Int_t iAlias = 0; iAlias < nEntries; iAlias++) {
226     AliTRDSensorArray * oneTRDDCS = (AliTRDSensorArray *)list->At (iAlias);
227     Log (Form ("%30s %10s %4d", 
228                oneTRDDCS->GetStoreName ().Data (),
229                results[iAlias] ? "ok" : "X",
230                nGraph [iAlias]));
231   }
232   Log ("*********** End of DCS **********");
233   
234   delete results;
235   delete nGraph;
236
237
238 }
239
240 //______________________________________________________________________________________________
241 void AliTRDPreprocessor::ExtractPedestals()
242 {
243   //
244   // Pedestal running on LDCs at the DAQ
245   //
246
247   // Init a AliTRDCalibPadStatus
248   AliTRDCalibPadStatus calPedSum = AliTRDCalibPadStatus();
249
250   AliCDBMetaData metaData;
251   metaData.SetBeamPeriod(0);
252   metaData.SetResponsible("Raphaelle Bailhache");
253   metaData.SetComment("TRD calib test");
254   
255   // Sum the contributions of the LDCs
256   TList * listpad = GetFileSources(kDAQ,"PADSTATUS");
257   if (listpad) {
258     
259     // loop through all files from LDCs
260
261     UInt_t index = 0;
262     while (listpad->At(index)!=NULL) {
263       TObjString* fileNameEntry = (TObjString*) listpad->At(index);
264       if (fileNameEntry != NULL) {
265         TString fileName = GetFile(kDAQ, "PADSTATUS",
266                                    fileNameEntry->GetString().Data());
267         if(fileName.Length() !=0){
268           TFile *f = TFile::Open(fileName);
269           AliTRDCalibPadStatus *calPed;
270           f->GetObject("calibpadstatus",calPed);
271           
272           if(calPed){
273
274             // analyse
275             calPed->AnalyseHisto();
276             
277             // store as reference data
278             TString name("PadStatus");
279             name += (Int_t)index;
280             if(!StoreReferenceData("DAQData",(const char *)name,(TObject *) calPed,&metaData)){
281               Log(Form("Error storing AliTRDCalibPadStatus object %d as reference data",(Int_t)index));
282               fResult |= kEStore;
283             }
284             
285                
286             // Add to the calPedSum
287             for (Int_t idet=0; idet<540; idet++) {
288               AliTRDCalROC *rocMean  = calPed->GetCalRocMean(idet, kFALSE);
289               if ( rocMean )  calPedSum.SetCalRocMean(rocMean,idet);
290               AliTRDCalROC *rocRMS = calPed->GetCalRocRMS(idet, kFALSE);
291               if ( rocRMS )  calPedSum.SetCalRocRMS(rocRMS,idet);
292             }// det loop
293             
294           } // calPed
295         }
296         else{
297           Log(Form("Error by retrieving the file %d for the pedestal",(Int_t)index));
298           fResult |= kEGetFileDAQ;  
299         }
300       }// fileNameEntry length
301       ++index;
302     }// while (list)
303     Log(Form("%d elements found in the list for the pedestal",(Int_t)index));
304     if(index==0){
305       fResult |= kEEmptyListDAQ;
306     }
307     //
308     // Store pedestal entry to OCDB
309     //
310
311
312     // Create Pad Status
313     AliTRDCalPadStatus *calPadStatus = calPedSum.CreateCalPadStatus();
314     AliCDBMetaData md3; 
315     md3.SetObjectClassName("AliTRDCalPadStatus");
316     md3.SetResponsible("Raphaelle Bailhache");
317     md3.SetBeamPeriod(1);
318     md3.SetComment("TRD calib test");
319     if(!Store("Calib","PadStatus"    ,(TObject *)calPadStatus, &md3, 0, kTRUE)){
320       Log("Error storing the pedestal");
321       fResult |= kEStore;    
322     }
323     delete listpad;
324   }
325   else {
326     Log("No list found for the PEDESTRAL Run");
327     fResult |= kEGetFileDAQ;    
328   }
329   
330
331 }
332 //______________________________________________________________________________________________
333 void AliTRDPreprocessor::ExtractDriftVelocityDAQ()
334 {
335   //
336   // Drift velocity DA running on monitoring servers at the DAQ
337   //
338
339
340   // Objects for HLT and DAQ zusammen
341   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
342   AliCDBMetaData metaData;
343   metaData.SetBeamPeriod(0);
344   metaData.SetResponsible("Raphaelle Bailhache");
345   metaData.SetComment("TRD calib test");
346   // Store the infos for the detector
347   AliCDBMetaData md1; 
348   md1.SetObjectClassName("AliTRDCalDet");
349   md1.SetResponsible("Raphaelle Bailhache");
350   md1.SetBeamPeriod(0);
351   md1.SetComment("TRD calib test");
352   // Store the infos for the pads
353   AliCDBMetaData md2; 
354   md2.SetObjectClassName("AliTRDCalPad");
355   md2.SetResponsible("Raphaelle Bailhache");
356   md2.SetBeamPeriod(0);
357   md2.SetComment("TRD calib test");
358
359
360
361
362   // Take the file from the DAQ file exchange server
363   TList *listdaq = GetFileSources(kDAQ,"VDRIFT");
364   if (listdaq) {
365     if(listdaq->GetSize() ==1){
366       TObjString* fileNameEntry = (TObjString*) listdaq->At(0);
367       if(fileNameEntry != NULL){
368         TString fileName = GetFile(kDAQ, "VDRIFT",
369                                    fileNameEntry->GetString().Data());
370         if(fileName.Length() !=0){
371           TFile *filedaq = TFile::Open(fileName);
372           TProfile2D *histodriftvelocity = (TProfile2D *) filedaq->Get("PH2d");
373           if (histodriftvelocity) {
374             
375             // store as reference data
376             if(!StoreReferenceData("DAQData","VdriftT0",(TObject *) histodriftvelocity,&metaData)){
377               Log("Error storing 2D Profile for vdrift from the DAQ");
378               fResult |= kEStore;
379             }
380             
381             // analyse
382             if(!fVdriftHLT){
383               Log("Take the PH reference data. Now we will try to fit\n");
384               calibra->SetMinEntries(2000); // If there is less than 2000
385               calibra->AnalysePH(histodriftvelocity);
386               
387               Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
388                 + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
389               Int_t nbfit        = calibra->GetNumberFit();
390               Int_t nbE        = calibra->GetNumberEnt();
391               
392               // if enough statistics store the results
393               if ((nbtg >                  0) && 
394                   (nbfit        >= 0.95*nbE)) {
395                 // create the cal objects
396                 TObjArray object      = calibra->GetVectorFit();
397                 AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
398                 TObject *objdriftvelocitypad = calibra->CreatePadObjectVdrift();
399                 object              = calibra->GetVectorFit2();
400                 AliTRDCalDet *objtime0det         = calibra->CreateDetObjectT0(&object,kTRUE);
401                 TObject *objtime0pad         = calibra->CreatePadObjectT0();
402                 calibra->ResetVectorFit();
403                 // store
404                 if(!Store("Calib","ChamberVdrift"    ,(TObject *) objdriftvelocitydet,&md1,0,kTRUE)){
405                   Log("Error storing the calibration object for the chamber vdrift (DAQ)");
406                   fResult |= kEStore;    
407                 }
408                 if(!Store("Calib","ChamberT0"        ,(TObject *) objtime0det        ,&md1,0,kTRUE)){
409                   Log("Error storing the calibration object for the chamber t0 (DAQ)");
410                   fResult |= kEStore;    
411                 }
412                 if(!Store("Calib","LocalVdrift"      ,(TObject *) objdriftvelocitypad,&md2,0,kTRUE)){
413                   Log("Error storing the calibration object for the local drift velocity (DAQ)");
414                   fResult |= kEStore;
415                 }
416                 if(!Store("Calib","LocalT0"          ,(TObject *) objtime0pad        ,&md2,0,kTRUE)){
417                   Log("Error storing the calibration object for the local time0 (DAQ)");
418                   fResult |= kEStore;
419                 }
420               }
421               else{
422                 Log("Not enough statistics for the average pulse height (DAQ)");
423               }
424             } // analyse
425           } // histo here
426         }
427         else{
428           Log("Error retrieving the file vdrift (DAQ)");
429           if(!fVdriftHLT) fResult |= kEGetFileDAQ;
430         }
431       }// if fileNameEntry
432     } // size DAQ list
433     else{
434       Log(Form("Problem on the size of the list: %d (DAQ)",listdaq->GetSize()));
435       if(!fVdriftHLT) fResult |= kEEmptyListDAQ;
436     }
437     delete listdaq; 
438   } 
439   else{
440     Log("No list found for vdrift (DAQ)");
441     if(!fVdriftHLT){
442       fResult |= kEGetFileDAQ;
443     }
444   }
445    
446 }
447 //______________________________________________________________________________________________
448 void AliTRDPreprocessor::ExtractHLT()
449 {
450   //
451   // Gain, vdrift and PRF calibration running on HLT
452   // return kFALSE if NULL pointer to the list
453   //
454
455  
456   // Objects for HLT and DAQ zusammen
457   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
458   AliCDBMetaData metaData;
459   metaData.SetBeamPeriod(0);
460   metaData.SetResponsible("Raphaelle Bailhache");
461   metaData.SetComment("TRD calib test");
462   // Store the infos for the detector
463   AliCDBMetaData md1; 
464   md1.SetObjectClassName("AliTRDCalDet");
465   md1.SetResponsible("Raphaelle Bailhache");
466   md1.SetBeamPeriod(0);
467   md1.SetComment("TRD calib test");
468   // Store the infos for the pads
469   AliCDBMetaData md2; 
470   md2.SetObjectClassName("AliTRDCalPad");
471   md2.SetResponsible("Raphaelle Bailhache");
472   md2.SetBeamPeriod(0);
473   md2.SetComment("TRD calib test");
474
475   
476   // Take the file from the HLT file exchange server
477   TList *listhlt = GetFileSources(kHLT,"GAINDRIFTPRF");
478   if (listhlt) {
479     if (listhlt->GetSize() == 1) {
480       TObjString* fileNameEntry = (TObjString*) listhlt->At(0);
481       if(fileNameEntry != NULL){
482         TString fileName = GetFile(kHLT, "GAINDRIFTPRF",
483                                    fileNameEntry->GetString().Data());
484         if(fileName.Length() !=0){
485           // Take the file
486           TFile *filehlt = TFile::Open(fileName);
487           
488           
489           // gain
490           TH2I *histogain = (TH2I *) filehlt->Get("CH2d");
491           histogain->SetDirectory(0);
492           if (histogain) {
493             // store the reference data
494             if(!StoreReferenceData("HLTData","Gain",(TObject *) histogain,&metaData)){
495               Log("Error storing 2D histos for gain");
496               fResult |= kEStore;
497             }
498             // analyse
499             Log("Take the CH reference data. Now we will try to fit\n");
500             calibra->SetMinEntries(1000); // If there is less than 1000 entries in the histo: no fit
501             calibra->AnalyseCH(histogain);
502             Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
503               + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
504             Int_t nbfit       = calibra->GetNumberFit();
505             Int_t nbE         = calibra->GetNumberEnt();
506             // enough statistics
507             if ((nbtg >                  0) && 
508                 (nbfit        >= 0.95*nbE)) {
509               // create the cal objects
510               TObjArray object           = calibra->GetVectorFit();
511               AliTRDCalDet *objgaindet   = calibra->CreateDetObjectGain(&object,calibra->GetScaleFitFactor(),kTRUE);
512               TObject *objgainpad        = calibra->CreatePadObjectGain();
513               // store them
514               if(!Store("Calib","ChamberGainFactor",(TObject *) objgaindet         ,&md1,0,kTRUE)){
515                 Log("Error storing the calibration object for the chamber gain");
516                 fResult |= kEStore;
517               }
518               if(!Store("Calib","LocalGainFactor"  ,(TObject *) objgainpad         ,&md2,0,kTRUE)){
519                 Log("Error storing the calibration object for the local gain factor");
520                 fResult |= kEStore;
521               }
522             }
523             calibra->ResetVectorFit();
524           }// if histogain
525           
526           
527           
528           // vdrift
529           fVdriftHLT = kFALSE;
530           TProfile2D *histodriftvelocity = (TProfile2D *) filehlt->Get("PH2d");
531           histodriftvelocity->SetDirectory(0);
532           if (histodriftvelocity) {
533             // store the reference data
534             if(!StoreReferenceData("HLTData","VdriftT0",(TObject *) histodriftvelocity,&metaData)){
535               Log("Error storing 2D Profile for average pulse height (HLT)");
536               fResult |= kEStore;
537             }
538             // analyse
539             Log("Take the PH reference data. Now we will try to fit\n");
540             calibra->SetMinEntries(1000*20); // If there is less than 20000
541             calibra->AnalysePH(histodriftvelocity);
542             Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
543               + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
544             Int_t nbfit        = calibra->GetNumberFit();
545             Int_t nbE          = calibra->GetNumberEnt();
546             // enough statistics
547             if ((nbtg >                  0) && 
548                 (nbfit        >= 0.95*nbE)) {
549               // create the cal objects
550               TObjArray object  = calibra->GetVectorFit();
551               AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
552               TObject *objdriftvelocitypad      = calibra->CreatePadObjectVdrift();
553               object              = calibra->GetVectorFit2();
554               AliTRDCalDet *objtime0det  = calibra->CreateDetObjectT0(&object,kTRUE);
555               TObject *objtime0pad       = calibra->CreatePadObjectT0();
556               // store them
557               if(!Store("Calib","ChamberVdrift"    ,(TObject *) objdriftvelocitydet,&md1,0,kTRUE)){
558                 Log("Error storing the calibration object for the chamber vdrift (HLT)");
559                 fResult |= kEStore;    
560               }
561               if(!Store("Calib","ChamberT0"        ,(TObject *) objtime0det        ,&md1,0,kTRUE)){
562                 Log("Error storing the calibration object for the chamber t0 (HLT)");
563                 fResult |= kEStore;    
564               }
565               if(!Store("Calib","LocalVdrift"      ,(TObject *) objdriftvelocitypad,&md2,0,kTRUE)){
566                 Log("Error storing the calibration object for the local drift velocity (HLT)");
567                 fResult |= kEStore;
568               }
569               if(!Store("Calib","LocalT0"          ,(TObject *) objtime0pad        ,&md2,0,kTRUE)){
570                 Log("Error storing the calibration object for the local time0 (HLT)");
571                 fResult |= kEStore;
572               }
573               fVdriftHLT = kTRUE;
574             }
575             calibra->ResetVectorFit();
576           }// if TProfile2D
577           
578           
579           // prf
580           TProfile2D *histoprf = (TProfile2D *) filehlt->Get("PRF2d");
581           histoprf->SetDirectory(0);
582           if (histoprf) {
583             // store reference data
584             if(!StoreReferenceData("HLTData","PRF",(TObject *) histoprf,&metaData)){
585               Log("Error storing the 2D Profile for Pad Response Function");
586               fResult |= kEStore;
587             }
588             // analyse
589             Log("Take the PRF reference data. Now we will try to fit\n");
590             calibra->SetMinEntries(600); // If there is less than 20000
591             calibra->AnalysePRFMarianFit(histoprf);
592             Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
593               + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
594             Int_t nbfit        = calibra->GetNumberFit();
595             Int_t nbE          = calibra->GetNumberEnt();
596             // enough statistics
597             if ((nbtg >                  0) && 
598                 (nbfit        >= 0.95*nbE)) {
599               // create cal pad objects 
600               TObjArray object            = calibra->GetVectorFit();
601               TObject *objPRFpad          = calibra->CreatePadObjectPRF(&object);
602               // store them
603               if(!Store("Calib","PRFWidth"         ,(TObject *) objPRFpad          ,&md2,0,kTRUE)){
604                 Log("Error storing the calibration object for the Pad Response Function");
605                 fResult |= kEStore;
606               }
607             }
608             calibra->ResetVectorFit();
609           }// if PRF
610         }
611         else{
612           Log("Error retrieving the file (HLT)");
613           fResult |= kEGetFileHLT;
614         } 
615       } // fileNameEntry
616     } // if HLT size of tlist
617     else{
618       Log(Form("Problem on the size of the list: %d (HLT)",listhlt->GetSize()));
619       fResult |= kEEmptyListHLT;
620     }
621     delete listhlt;
622   } //if HLT tlist
623   else{
624     Log("No list found for the HLT");
625     fResult |= kEGetFileHLT;
626   }
627   
628 }