]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDCalibrationComponent.cxx
bugfixes by Theodor and Raphaelle
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDCalibrationComponent.cxx
1 // $Id$
2
3 /**************************************************************************
4  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Authors: Matthias Richter <Matthias.Richter@ift.uib.no>                *
7  *          Timm Steinbeck <timm@kip.uni-heidelberg.de>                   *
8  *          for The ALICE Off-line Project.                               *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19 /** @file   AliHLTTRDCalibrationComponent.cxx
20     @author Timm Steinbeck, Matthias Richter
21     @date
22     @brief  A TRDCalibration processing component for the HLT. */
23
24 #if __GNUC__ >= 3
25 using namespace std;
26 #endif
27
28 #include "TTree.h"
29 #include "TFile.h"
30 #include "TBranch.h"
31 #include "TH2I.h"
32 #include "TH2.h"
33 #include "TProfile2D.h"
34
35 #include "AliHLTReadoutList.h"
36
37 #include "AliHLTTRDCalibrationComponent.h"
38 #include "AliHLTTRDDefinitions.h"
39 #include "AliHLTTRDUtils.h"
40
41 #include "AliCDBManager.h"
42 #include "AliCDBStorage.h"
43 #include "AliRawReaderMemory.h"
44
45 #include "AliTRDCalPad.h"
46 #include "AliTRDCalDet.h"
47
48 #include "AliTRDCalibraFillHisto.h"
49 #include "AliTRDtrackV1.h"
50
51 #include "AliTRDCalibraFit.h"
52 #include "AliTRDCalibraMode.h"
53 #include "AliTRDCalibraVector.h"
54 #include "AliTRDCalibraVdriftLinearFit.h"
55
56 #include <cstdlib>
57 #include <cerrno>
58 #include <string>
59
60 ClassImp(AliHLTTRDCalibrationComponent);
61
62 AliHLTTRDCalibrationComponent::AliHLTTRDCalibrationComponent()
63   : AliHLTCalibrationProcessor(),
64     fTRDCalibraFillHisto(NULL),
65     fOutputSize(500000),
66     fTracksArray(NULL),
67     fOutArray(NULL),
68     fAfterRunArray(NULL),
69     fDisplayArray(NULL),
70     fNevent(0),
71     feveryNevent(1000),
72     fRecievedTimeBins(kFALSE),
73     fTrgStrings(NULL),
74     fAccRejTrg(0)
75 {
76   // Default constructor
77 }
78
79 AliHLTTRDCalibrationComponent::~AliHLTTRDCalibrationComponent()
80 {
81   // Destructor
82 }
83
84 const char* AliHLTTRDCalibrationComponent::GetComponentID()
85 {
86   // Return the component ID const char *
87   return "TRDCalibration"; // The ID of this component
88 }
89
90 void AliHLTTRDCalibrationComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
91 {
92   // Get the list of input data
93   list.clear(); // We do not have any requirements for our input data type(s).
94   list.push_back(AliHLTTRDDefinitions::fgkTRDSATracksDataType);
95 }
96
97 AliHLTComponentDataType AliHLTTRDCalibrationComponent::GetOutputDataType()
98 {
99   // Get the output data type
100   return kAliHLTMultipleDataType;
101   //  return AliHLTTRDDefinitions::fgkCalibrationDataType;
102  
103 }
104
105 int AliHLTTRDCalibrationComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
106 {
107   // Get the output data type
108   tgtList.clear();
109   tgtList.push_back(AliHLTTRDDefinitions::fgkCalibrationDataType);
110   tgtList.push_back(AliHLTTRDDefinitions::fgkEORCalibrationDataType);
111   return tgtList.size();
112 }
113
114 void AliHLTTRDCalibrationComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
115 {
116   // Get the output data size
117   constBase = fOutputSize;
118   inputMultiplier = 0;
119 }
120
121 AliHLTComponent* AliHLTTRDCalibrationComponent::Spawn()
122 {
123   // Spawn function, return new instance of this class
124   return new AliHLTTRDCalibrationComponent;
125 };
126
127 Int_t AliHLTTRDCalibrationComponent::ScanArgument( int argc, const char** argv )
128 {
129   // perform initialization. We check whether our relative output size is specified in the arguments.
130   int i = 0;
131   char* cpErr;
132   if(!fTrgStrings)
133     fTrgStrings = new TObjArray();
134
135   while ( i < argc )
136     {
137       HLTDebug("argv[%d] == %s", i, argv[i] );
138       if ( !strcmp( argv[i], "output_size" ) )
139         {
140           if ( i+1>=argc )
141             {
142               HLTError("Missing output_size parameter");
143               return ENOTSUP;
144             }
145           HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
146           fOutputSize = strtoul( argv[i+1], &cpErr, 0 );
147           if ( *cpErr )
148             {
149               HLTError("Cannot convert output_size parameter '%s'", argv[i+1] );
150               return EINVAL;
151             }
152           HLTInfo("Output size set to %lu %%", fOutputSize );
153           i += 2;
154           continue;
155         }
156       if ( !strcmp( argv[i], "-everyNevent" ) )
157         {
158           if ( i+1>=argc )
159             {
160               HLTError("Missing everyNevent parameter");
161               return ENOTSUP;
162             }
163           HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
164           fOutputSize = strtoul( argv[i+1], &cpErr, 0 );
165           if ( *cpErr )
166             {
167               HLTError("Cannot convert everyNevent parameter '%s'", argv[i+1] );
168               return EINVAL;
169             }
170           HLTInfo("Pushing back every %d event", feveryNevent);
171           i += 2;
172           continue;
173         }
174       if ( !strcmp( argv[i], "-TrgStr" ) )
175         {
176           if ( i+1>=argc )
177             {
178               HLTError("Missing parameter for mbTriggerString");
179               return ENOTSUP;
180             }
181           HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
182           fTrgStrings->Add(new TObjString(argv[i+1]));
183           i += 2;
184           continue;
185         }
186
187       if ( !strcmp( argv[i], "-acceptTrgStr" ) )
188         {
189           fAccRejTrg=1;
190           i += 1;
191           continue;
192         }
193       if ( !strcmp( argv[i], "-rejectTrgStr" ) )
194         {
195           fAccRejTrg=-1;
196           i += 1;
197           continue;
198         }
199
200       else {
201         HLTError("Unknown option '%s'", argv[i] );
202         return EINVAL;
203       }
204     }
205   return i;
206 }
207
208 Int_t AliHLTTRDCalibrationComponent::InitCalibration()
209 {
210
211   if(!fTrgStrings)
212     fTrgStrings = new TObjArray();
213
214   if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
215     HLTError("DefaultStorage is not set in CDBManager");
216     return -EINVAL;
217   }
218   if(AliCDBManager::Instance()->GetRun()<0){
219     HLTError("Run Number is not set in CDBManager");
220     return -EINVAL;
221   }
222   HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
223
224   if(fTrgStrings->GetEntriesFast()>0 && !fAccRejTrg){
225     HLTError("Trigger string(s) given, but acceptTrgStr or rejectTrgStr not selected");
226     return -EINVAL;
227   }
228
229   fTRDCalibraFillHisto = AliTRDCalibraFillHisto::Instance();
230   fTRDCalibraFillHisto->SetIsHLT(kTRUE);
231   fTRDCalibraFillHisto->SetHisto2d(); // choose to use histograms
232   fTRDCalibraFillHisto->SetCH2dOn();  // choose to calibrate the gain
233   fTRDCalibraFillHisto->SetPH2dOn();  // choose to calibrate the drift velocity
234   fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
235   fTRDCalibraFillHisto->SetIsHLT(); // per detector
236   //fTRDCalibraFillHisto->SetDebugLevel(1);// debug
237   fTRDCalibraFillHisto->SetFillWithZero(kTRUE);
238   fTRDCalibraFillHisto->SetLinearFitterOn(kTRUE);
239   fTRDCalibraFillHisto->SetNumberBinCharge(100);
240   
241   fTracksArray = new TClonesArray("AliTRDtrackV1");
242   fOutArray = new TObjArray(4);
243   fAfterRunArray=new TObjArray(5);
244   fDisplayArray=new TObjArray(4);
245
246   HLTDebug("run SetupCTPData");
247   SetupCTPData();
248
249   return 0;
250 }
251
252 Int_t AliHLTTRDCalibrationComponent::DeinitCalibration()
253 {
254   
255   // Deinitialization of the component
256   
257   HLTDebug("DeinitCalibration");
258   delete fTracksArray; fTracksArray=0;
259   //fTRDCalibraFillHisto->Destroy();
260   //fOutArray->Delete();
261   delete fOutArray; fOutArray=0;
262   fAfterRunArray->Delete();
263   delete fAfterRunArray; fAfterRunArray=0;
264   fDisplayArray->Delete();
265   delete fDisplayArray; fDisplayArray=0;
266   fTrgStrings->Delete();
267   delete fTrgStrings; fTrgStrings=0;
268   return 0;
269 }
270
271 Int_t AliHLTTRDCalibrationComponent::ProcessCalibration(const AliHLTComponent_EventData& evtData,
272                                                         const AliHLTComponent_BlockData* blocks,
273                                                         AliHLTComponent_TriggerData& trigData,
274                                                         AliHLTUInt8_t* /*outputPtr*/,
275                                                         AliHLTUInt32_t& /*size*/,
276                                                         vector<AliHLTComponent_BlockData>& /*outputBlocks*/)
277 {
278   HLTDebug("NofBlocks %lu", evtData.fBlockCnt );
279   // Process an event
280         
281   
282   // Loop over all input blocks in the event
283   vector<AliHLTComponent_DataType> expectedDataTypes;
284   GetInputDataTypes(expectedDataTypes);
285   for ( unsigned long iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ )
286     {
287       const AliHLTComponentBlockData &block = blocks[iBlock];
288       AliHLTComponentDataType inputDataType = block.fDataType;
289       Bool_t correctDataType = kFALSE;
290
291       for(UInt_t i = 0; i < expectedDataTypes.size(); i++)
292         if( expectedDataTypes.at(i) == inputDataType)
293           correctDataType = kTRUE;
294       if (!correctDataType) {
295         HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - Skipping",
296                   iBlock, evtData.fBlockCnt,
297                   evtData.fEventID, evtData.fEventID,
298                   DataType2Text(inputDataType).c_str());
299         continue;
300       }
301       else {
302         HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
303                  iBlock, evtData.fBlockCnt-1,
304                  evtData.fEventID, evtData.fEventID,
305                  DataType2Text(inputDataType).c_str(),
306                  block.fSize);
307       }
308
309       Int_t nTimeBins;
310       AliHLTTRDUtils::ReadTracks(fTracksArray, block.fPtr, block.fSize, &nTimeBins);
311       
312           
313       if(!fRecievedTimeBins){
314         HLTDebug("Reading number of time bins from input block. Value is: %d", nTimeBins);
315         fTRDCalibraFillHisto->Init2Dhistos(nTimeBins); // initialise the histos
316         fTRDCalibraFillHisto->SetNumberClusters(0); // At least 1 clusters
317         fTRDCalibraFillHisto->SetNumberClustersf(nTimeBins); // Not more than %d  clusters
318         fRecievedTimeBins=kTRUE;
319       }
320       
321         
322       Bool_t TriggerPassed=kFALSE;
323                 
324       if(fAccRejTrg){
325         if(fAccRejTrg>0){
326           TriggerPassed=kFALSE;
327           for(int i = 0; i < fTrgStrings->GetEntriesFast(); i++){
328             const TObjString *const obString=(TObjString*)fTrgStrings->At(i);
329             const TString tString=obString->GetString();
330             //printf("Trigger Output: %i\n",EvaluateCTPTriggerClass(tString.Data(),trigData));
331             if(EvaluateCTPTriggerClass(tString.Data(),trigData)){TriggerPassed=kTRUE; break;}
332           }
333         }
334         else{
335           TriggerPassed=kTRUE;
336           for(int i = 0; i < fTrgStrings->GetEntriesFast(); i++){
337             const TObjString *const obString=(TObjString*)fTrgStrings->At(i);
338             const TString tString=obString->GetString();
339             if(EvaluateCTPTriggerClass(tString.Data(),trigData)){TriggerPassed=kFALSE; break;}
340           }
341         }
342       }
343       
344       
345       Int_t nbEntries = fTracksArray->GetEntries();
346       HLTDebug(" %i TRDtracks in tracksArray", nbEntries);
347       AliTRDtrackV1* trdTrack = 0x0;
348       for (Int_t i = 0; i < nbEntries; i++){
349         HLTDebug("%i/%i: ", i+1, nbEntries);
350         trdTrack = (AliTRDtrackV1*)fTracksArray->At(i);
351         //trdTrack->Print();
352         fTRDCalibraFillHisto->SetCH2dOn(TriggerPassed);
353         fTRDCalibraFillHisto->UpdateHistogramsV1(trdTrack);
354         fTRDCalibraFillHisto->SetCH2dOn(kTRUE);
355       }
356       
357
358       if(!fOutArray->At(0))FormOutput(0);
359       if(!fDisplayArray->At(0))FormOutput(1);
360       if (fNevent%feveryNevent==0 && fOutArray) {
361         PushBack(fDisplayArray, AliHLTTRDDefinitions::fgkCalibrationDataType);
362       }
363
364       fTracksArray->Delete();
365       fNevent++;
366
367     }
368  
369   return 0;
370
371 }
372
373 /**
374  * Form output array of histrograms
375  */
376 //============================================================================
377 void AliHLTTRDCalibrationComponent::FormOutput(Int_t param)
378 {
379   // gain histo
380   TH2I *hCH2d = fTRDCalibraFillHisto->GetCH2d();
381   if(!param)fOutArray->Add(hCH2d);
382   else fDisplayArray->Add(hCH2d);
383
384   // drift velocity histo
385   TProfile2D *hPH2d = fTRDCalibraFillHisto->GetPH2d();
386   if(!param)fOutArray->Add(hPH2d);
387   else fDisplayArray->Add(hPH2d);
388
389   // PRF histo
390   TProfile2D *hPRF2d = fTRDCalibraFillHisto->GetPRF2d();
391   if(!param)fOutArray->Add(hPRF2d);
392   else fDisplayArray->Add(hPRF2d);
393
394   // Vdrift Linear Fit
395   if(!param){
396     AliTRDCalibraVdriftLinearFit *hVdriftLinearFitOne=(AliTRDCalibraVdriftLinearFit *)fTRDCalibraFillHisto->GetVdriftLinearFit();
397     fOutArray->Add(hVdriftLinearFitOne);
398   }
399   else{
400     TH2S *hVdriftLinearFitOne = (TH2S *)(((AliTRDCalibraVdriftLinearFit *)fTRDCalibraFillHisto->GetVdriftLinearFit())->GetLinearFitterHisto(10,kTRUE)); 
401     fDisplayArray->Add(hVdriftLinearFitOne);
402   }
403
404   HLTDebug("GetCH2d = 0x%x; NEntries = %i; size = %i", hCH2d, hCH2d->GetEntries(), sizeof(*hCH2d));
405   hCH2d->Print();
406   HLTDebug("GetPH2d = 0x%x; NEntries = %i; size = %i", hPH2d, hPH2d->GetEntries(), sizeof(*hPH2d));
407   hPH2d->Print();
408   HLTDebug("GetPRF2d = 0x%x; NEntries = %i; size = %i", hPRF2d, hPRF2d->GetEntries(), sizeof(*hPRF2d));
409   hPRF2d->Print();
410   //HLTDebug("GetVdriftLinearFit = 0x%x; size = %i", hVdriftLinearFitOne, sizeof(hVdriftLinearFitOne)); 
411   
412   HLTDebug("output Array: pointer = 0x%x; NEntries = %i; size = %i", fOutArray, fOutArray->GetEntries(), sizeof(fOutArray));
413    
414 }
415
416 Int_t AliHLTTRDCalibrationComponent::ShipDataToFXS(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
417 {
418   //fTRDCalibraFillHisto->DestroyDebugStreamer();
419
420   AliHLTReadoutList rdList(AliHLTReadoutList::kTRD);
421
422   EORCalibration();
423   
424   fOutArray->Remove(fOutArray->FindObject("AliTRDCalibraVdriftLinearFit"));
425   //fOutArray->Remove(fOutArray->FindObject("PRF2d"));
426   //fOutArray->Remove(fOutArray->FindObject("PH2d"));
427   //fOutArray->Remove(fOutArray->FindObject("CH2d"));
428
429   //if(!(fOutArray->FindObject("CH2d"))) {
430   //  TH2I * ch2d = new TH2I("CH2d","Nz0Nrphi0",100,0.0,300.0,540,0,540);
431   //  fOutArray->Add(ch2d);
432   //}
433
434   //if(!(fOutArray->FindObject("PH2d"))) {
435   //  TProfile2D * ph2d = new TProfile2D("PH2d","Nz0Nrphi0",30,-0.05,2.95,540,0,540);
436   //  fOutArray->Add(ph2d);
437   //}
438
439   //if(!(fOutArray->FindObject("PRF2d"))) {
440   //  TProfile2D * prf2d = new TProfile2D("PRF2d","Nz0Nrphi0Ngp3",60,-9.0,9.0,540,0,540);
441   //  fOutArray->Add(prf2d);
442   //}
443
444
445   HLTDebug("Size of the fOutArray is %d\n",fOutArray->GetEntriesFast());
446
447   /*
448   TString fileName="$ALIHLT_TOPDIR/build-debug/output/CalibHistoDump_run";
449   fileName+=".root";
450   HLTInfo("Dumping Histogram file to %s",fileName.Data());
451   TFile* file = TFile::Open(fileName, "RECREATE");
452   //fAfterRunArray->Write();
453   fOutArray->Write();
454   file->Close();
455   HLTInfo("Histogram file dumped");
456   */
457   
458   PushToFXS((TObject*)fOutArray, "TRD", "GAINDRIFTPRF", rdList.Buffer() );
459   //PushToFXS((TObject*)fOutArray->FindObject("CH2d"), "TRD", "GAINDRIFTPRF", rdList.Buffer() );
460
461         
462   return 0;
463 }
464 Int_t AliHLTTRDCalibrationComponent::EORCalibration()
465 {
466   //Also Fill histograms for the online display
467   TH2I *hCH2d=(TH2I*)fOutArray->FindObject("CH2d");
468   TProfile2D *hPH2d=(TProfile2D*)fOutArray->FindObject("PH2d");
469   TProfile2D *hPRF2d= (TProfile2D*)fOutArray->FindObject("PRF2d");
470   AliTRDCalibraVdriftLinearFit* hVdriftLinearFit = (AliTRDCalibraVdriftLinearFit*)fOutArray->FindObject("AliTRDCalibraVdriftLinearFit");
471  
472
473   if(!hCH2d || !hPH2d || !hPRF2d || !hVdriftLinearFit) return 0; 
474
475   //Fit
476   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
477
478   //Gain
479   calibra->SetMinEntries(100);
480   calibra->AnalyseCH(hCH2d);
481   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
482     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
483   Int_t nbfit       = calibra->GetNumberFit();
484   Int_t nbE         = calibra->GetNumberEnt();
485   TH1F *coefgain = 0x0;
486   // enough statistics
487   //if ((nbtg >                  0) && 
488   //   (nbfit        >= 0.2*nbE)) {
489   // create the cal objects
490   //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
491   TObjArray object           = calibra->GetVectorFit();
492   AliTRDCalDet *objgaindet   = calibra->CreateDetObjectGain(&object,kFALSE);
493   coefgain                   = objgaindet->MakeHisto1DAsFunctionOfDet();
494   //}
495   calibra->ResetVectorFit();
496
497   // vdrift second method
498   calibra->SetMinEntries(100); // If there is less than 100
499   hVdriftLinearFit->FillPEArray();
500   calibra->AnalyseLinearFitters(hVdriftLinearFit);
501   nbtg = 540;
502   nbfit = calibra->GetNumberFit();
503   nbE   = calibra->GetNumberEnt();
504   TH1F *coefdriftsecond = 0x0;
505   // enough statistics
506   //if ((nbtg >                  0) && 
507   // (nbfit        >= 0.1*nbE)) {
508   // create the cal objects
509   //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
510   object  = calibra->GetVectorFit();
511   AliTRDCalDet *objdriftvelocitydetsecond = calibra->CreateDetObjectVdrift(&object,kTRUE);
512   objdriftvelocitydetsecond->SetTitle("secondmethodvdrift");
513   coefdriftsecond  = objdriftvelocitydetsecond->MakeHisto1DAsFunctionOfDet();
514   //}
515   calibra->ResetVectorFit();
516   
517   // vdrift first method
518   calibra->SetMinEntries(100*20); // If there is less than 20000
519   calibra->AnalysePH(hPH2d);
520   nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
521     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
522   nbfit        = calibra->GetNumberFit();
523   nbE          = calibra->GetNumberEnt();
524   TH1F *coefdrift = 0x0;
525   TH1F *coeft0 = 0x0;
526   // enough statistics
527   //if ((nbtg >                  0) && 
528   // (nbfit        >= 0.2*nbE)) {
529   // create the cal objects
530   //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
531   //calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
532   object  = calibra->GetVectorFit();
533   AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
534   coefdrift        = objdriftvelocitydet->MakeHisto1DAsFunctionOfDet();
535   object              = calibra->GetVectorFit2();
536   AliTRDCalDet *objtime0det  = calibra->CreateDetObjectT0(&object,kTRUE);
537   coeft0        = objtime0det->MakeHisto1DAsFunctionOfDet();
538   //}
539   calibra->ResetVectorFit();
540            
541
542   //PRF
543   calibra->SetMinEntries(200); 
544   calibra->AnalysePRFMarianFit(hPRF2d);
545   nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
546     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
547   nbfit        = calibra->GetNumberFit();
548   nbE          = calibra->GetNumberEnt();
549   TH1F *coefprf = 0x0;
550   // enough statistics
551   //if ((nbtg >                  0) && 
552   //  (nbfit        >= 0.95*nbE)) {
553   // create cal pad objects 
554   object            = calibra->GetVectorFit();
555   TObject *objPRFpad          = calibra->CreatePadObjectPRF(&object);
556   coefprf                     = ((AliTRDCalPad *) objPRFpad)->MakeHisto1D();
557   //}
558   calibra->ResetVectorFit();
559
560
561   coefgain->SetName("coefgain");
562   coefprf->SetName("coefprf");
563   coefdrift->SetName("coefdrift");
564   coefdriftsecond->SetName("coefdriftsecond");
565   coeft0->SetName("coeft0");
566   if(coefgain) fAfterRunArray->Add(coefgain);
567   if(coefprf) fAfterRunArray->Add(coefprf);
568   if(coefdrift) fAfterRunArray->Add(coefdrift);
569   if(coefdriftsecond) fAfterRunArray->Add(coefdriftsecond);
570   if(coeft0) fAfterRunArray->Add(coeft0);
571   
572
573   if(coefgain||coefprf||coefdrift||coeft0||coefdriftsecond) {
574     PushBack(fAfterRunArray, AliHLTTRDDefinitions::fgkEORCalibrationDataType);
575   }
576
577   /*
578   TString fileName="$ALIHLT_TOPDIR/build-debug/output/CalibHistoDump_run";
579   fileName+=".root";
580   HLTInfo("Dumping Histogram file to %s",fileName.Data());
581   TFile* file = TFile::Open(fileName, "RECREATE");
582   fAfterRunArray->Write();
583   fOutArray->Write();
584   file->Close();
585   HLTInfo("Histogram file dumped");
586   */
587
588   return 0;
589 }       
590