]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDPreprocessor.cxx
VZERO volumes now placed in common mother.
[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 //                                                                        //
25 // Author:                                                                //
26 //   R. Bailhache (R.Bailhache@gsi.de)                                    //
27 //                                                                        //
28 ////////////////////////////////////////////////////////////////////////////
29
30 #include "AliTRDPreprocessor.h"
31
32 #include <TTimeStamp.h>
33 #include <TFile.h>
34 #include <TProfile2D.h>
35 #include <TH2I.h>
36 #include <TStopwatch.h>
37 #include <TObjString.h>
38 #include <TString.h>
39 #include <TList.h>
40 #include <TCollection.h>
41
42 #include "AliCDBMetaData.h"
43 #include "AliDCSValue.h"
44 #include "AliLog.h"
45
46 #include "AliTRDCalibraFit.h"
47 #include "AliTRDCalibraMode.h"
48 #include "Cal/AliTRDCalDet.h"
49
50 ClassImp(AliTRDPreprocessor)
51
52 //______________________________________________________________________________________________
53 AliTRDPreprocessor::AliTRDPreprocessor(AliShuttleInterface *shuttle)
54                    :AliPreprocessor("TRD", shuttle)
55 {
56   //
57   // Constructor
58   //
59
60 }
61
62 //______________________________________________________________________________________________
63 AliTRDPreprocessor::~AliTRDPreprocessor()
64 {
65   //
66   // Destructor
67   //
68
69 }
70
71 //______________________________________________________________________________________________
72 void AliTRDPreprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
73 {
74   //
75   // Initialization routine for the TRD preprocessor
76   //
77
78   AliPreprocessor::Initialize(run,startTime,endTime);
79
80 }
81
82 //______________________________________________________________________________________________
83 UInt_t AliTRDPreprocessor::Process(TMap* /*dcsAliasMap*/)
84 {
85   //
86   // Process the calibration data for the HLT part
87   //
88
89   // How long does it take for the HLT part?
90   TStopwatch timer;
91   timer.Start();
92
93   // Metadata for the reference data
94   AliCDBMetaData metaData;
95   metaData.SetBeamPeriod(1);
96   metaData.SetResponsible("Raphaelle Bailhache");
97   metaData.SetComment("This preprocessor fills reference data.");
98
99   // Take the file from the HLT file exchange server
100   TList *filesources = GetFileSources(kHLT,"GAINDRIFTPRF");
101   if (!filesources) {
102     AliError(Form("No sources found for GAINDRIFTPRF for run %d !",fRun));
103     return 0;
104   }
105   if (filesources->GetSize() != 1) {
106     AliError(Form("More than one source found for GAINDRIFTPRF for run %d!",fRun));
107     return 0;
108   }
109
110   // Call a AliTRDCalibra instance for fit
111   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
112
113   //Choose the fit methods
114   calibra->SetFitChargeNDB(4); //for the relative gain
115   calibra->SetFitMeanWOn();   //weighted mean
116   calibra->SetFitPHNDB(3);    //for the average pulse height
117   calibra->SetFitLagrPolOn(); //LagrPol
118   calibra->SetFitPRFNDB(0);   //for the PRF
119   calibra->SetFitPRFOn();     //gaussian fit
120
121   //Debug mode
122   //calibra->SetDebug(1);       //Debug
123
124   // Init some things
125   AliTRDCalDet *objgaindet          = 0x0; // Object for det average gain factor
126   AliTRDCalDet *objdriftvelocitydet = 0x0; // Object for det average drift velocity
127   AliTRDCalDet *objtime0det         = 0x0; // Object for det average time0 
128   TObject      *objgainpad          = 0x0; // Object for pad (relative to the det) gain factor
129   TObject      *objdriftvelocitypad = 0x0; // Object for pad (relative to the det) drift velocity
130   TObject      *objtime0pad         = 0x0; // Object for pad (relative to the det) time0
131   TObject      *objPRFpad           = 0x0; // Object for pad prf width
132   TH2I         *histogain           = 0x0; // Histogram taken from HLT for gain factor
133   TProfile2D   *histodriftvelocity  = 0x0; // Profile taken from HLT for drift velocity and time0
134   TProfile2D   *histoprf            = 0x0; // Profile taken from HLT for prf
135
136   Int_t    numberfit[3]        = { 0,   0,   0   }; // Number of histos fitted for gain, drift velocity and prf
137   Int_t    numberEnt[3]        = { 0,   0,   0   }; // Number of histos with entries
138   Double_t statisticmean[3]    = { 0.0, 0.0, 0.0 }; // Mean values of the number of entries in these histos
139   Int_t    numbertotalgroup[3] = { 0,   0,   0   }; // Total number of groups
140
141   // Loop over the files taken from the HLT
142   TIter iter(filesources);
143   TObjString *source;
144   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
145     
146     TString filename = GetFile(kHLT,"GAINDRIFTPRF",source->GetName());
147     if (filename.Length() == 0) {
148       AliError(Form("Error retrieving file from source %d failed!", source->GetName()));
149       delete filesources;
150       return 0;
151     }
152
153     // Take the histos
154     TFile *file = TFile::Open(filename);
155     histogain = (TH2I *) file->Get("CH2d");
156     histogain->SetDirectory(0);
157     if (!histogain) {
158       AliError("Error retrieving 2D histos for gain failed!");
159     }
160     histodriftvelocity = (TProfile2D *) file->Get("PH2d");
161     histodriftvelocity->SetDirectory(0);
162     if (!histodriftvelocity) {
163       AliError("Error retrieving 2D Profile for average pulse height failed!");
164     }
165     histoprf = (TProfile2D *) file->Get("PRF2d");
166     histoprf->SetDirectory(0);
167     if (!histoprf) {
168       AliError("Error retrieving 2D Profile for Pad Response Function failed!");
169     }
170     file->Close();
171
172     // Set the mode of calibration from the TObject, store the reference data and try to fit them
173     if (histogain) {
174       calibra->SetModeCalibrationFromTObject((TObject *) histogain,0);
175       StoreReferenceData("HLTData","Gain",(TObject *) histogain,&metaData);
176       AliInfo("Take the CH reference data. Now we will try to fit\n");
177       calibra->SetMinEntries(100); // If there is less than 100 entries in the histo: no fit
178       calibra->FitCHOnline(histogain);
179       numbertotalgroup[0] = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
180                           + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
181       numberfit[0]        = calibra->GetNumberFit();
182       statisticmean[0]    = calibra->GetStatisticMean(); 
183       numberEnt[0]        = calibra->GetNumberEnt();
184       objgaindet          = calibra->CreateDetObjectTree(calibra->GetGain(),0);
185       objgainpad          = calibra->CreatePadObjectTree(calibra->GetGain(),0,objgaindet);
186     }
187     
188     if (histodriftvelocity) {
189       calibra->SetModeCalibrationFromTObject((TObject *) histodriftvelocity,1);
190       StoreReferenceData("HLTData","VdriftT0",(TObject *) histodriftvelocity,&metaData);
191       AliInfo("Take the PH reference data. Now we will try to fit\n");
192       calibra->SetMinEntries(100*20); // If there is less than 2000
193       calibra->FitPHOnline(histodriftvelocity);
194       numbertotalgroup[1] = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
195                           + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
196       numberfit[1]        = calibra->GetNumberFit();
197       statisticmean[1]    = calibra->GetStatisticMean(); 
198       numberEnt[1]        = calibra->GetNumberEnt();
199       objdriftvelocitydet = calibra->CreateDetObjectTree(calibra->GetVdrift(),1);
200       objdriftvelocitypad = calibra->CreatePadObjectTree(calibra->GetVdrift(),1,objdriftvelocitydet);
201       objtime0det         = calibra->CreateDetObjectTree(calibra->GetT0(),3);
202       objtime0pad         = calibra->CreatePadObjectTree(calibra->GetT0(),3,objtime0det);
203     }
204     
205     if (histoprf) {
206       calibra->SetModeCalibrationFromTObject((TObject *) histoprf,2);
207       StoreReferenceData("HLTData","PRF",(TObject *) histoprf,&metaData);
208       AliInfo("Take the PRF reference data. Now we will try to fit\n");
209       calibra->SetMinEntries(100*20); // If there is less than 2000
210       calibra->SetRangeFitPRF(0.5);
211       calibra->FitPRFOnline(histoprf);
212       numbertotalgroup[2] = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
213                            + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
214       numberfit[2]        = calibra->GetNumberFit();
215       statisticmean[2]    = calibra->GetStatisticMean(); 
216       numberEnt[2]        = calibra->GetNumberEnt();
217       objPRFpad           = calibra->CreatePadObjectTree(calibra->GetPRF());
218     }
219
220   }
221   
222   // Bilan of the fit statistic
223   AliInfo(Form("The mean number of entries required for a fit is: %d"
224               ,(Int_t) calibra->GetMinEntries()));
225   AliInfo(Form("FOR THE CH: There is a mean statistic of: %f, with %d fits for %d groups and %d histos with entries"
226               ,statisticmean[0],numberfit[0],numbertotalgroup[0],numberEnt[0]));
227   AliInfo(Form("FOR THE PH: There is a mean statistic of: %f, with %d fits for %d groups and %d histos with entries"
228               ,statisticmean[1],numberfit[1],numbertotalgroup[1],numberEnt[1]));
229   AliInfo(Form("FOR THE PRF: There is a mean statistic of: %f, with %d fits for %d groups and %d histos with entries"
230               ,statisticmean[2],numberfit[2],numbertotalgroup[2],numberEnt[2]));
231   
232   
233   //
234   // Store the coefficients in the grid OCDB if enough statistics
235   //
236   
237   // Store the infos for the detector
238   AliCDBMetaData *md1= new AliCDBMetaData(); 
239   md1->SetObjectClassName("AliTRDCalDet");
240   md1->SetResponsible("Raphaelle Bailhache");
241   md1->SetBeamPeriod(1);
242   md1->SetAliRootVersion("01-10-06"); // root version
243   md1->SetComment("The dummy values in this calibration file are for testing only");
244   if ((numbertotalgroup[0] >                  0) && 
245       (numberfit[0]        >= 0.95*numberEnt[0])) {
246     Store("Calib","ChamberGainFactor",(TObject *) objgaindet         ,md1,0,kTRUE);
247   }
248   if ((numbertotalgroup[1] >                  0) && 
249       (numberfit[1]        >= 0.95*numberEnt[1])) {
250     Store("Calib","ChamberVdrift"    ,(TObject *) objdriftvelocitydet,md1,0,kTRUE);
251     Store("Calib","ChamberT0"        ,(TObject *) objtime0det        ,md1,0,kTRUE);
252   }
253   
254   // Store the infos for the pads
255   AliCDBMetaData *md2= new AliCDBMetaData(); 
256   md2->SetObjectClassName("AliTRDCalPad");
257   md2->SetResponsible("Raphaelle Bailhache");
258   md2->SetBeamPeriod(1);
259   md2->SetAliRootVersion("01-10-06"); //root version
260   md2->SetComment("The dummy values in this calibration file are for testing only");
261   if ((numbertotalgroup[0] >                  0) && 
262       (numberfit[0]        >= 0.95*numberEnt[0])) {
263     Store("Calib","LocalGainFactor"  ,(TObject *) objgainpad         ,md2,0,kTRUE);
264   }
265   if ((numbertotalgroup[1] >                  0) && 
266       (numberfit[1]        >= 0.95*numberEnt[1])) {
267     Store("Calib","LocalVdrift"      ,(TObject *) objdriftvelocitypad,md2,0,kTRUE);
268     Store("Calib","LocalT0"          ,(TObject *) objtime0pad        ,md2,0,kTRUE);
269   }
270   if ((numbertotalgroup[2] >                  0) && 
271       (numberfit[2]        >= 0.95*numberEnt[2])) {
272     Store("Calib","PRFWidth"         ,(TObject *) objPRFpad          ,md2,0,kTRUE);
273   }
274   
275   // End
276   delete filesources;
277   timer.Stop();
278   timer.Print();
279   return 1;  
280
281 }