]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
Include pressure sensors. Read Config entry from OCDB (Haavard)
[u/mrichter/AliRoot.git] / TPC / AliTPCPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, 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 #include "AliTPCPreprocessor.h"
18
19 #include "AliCDBMetaData.h"
20 #include "AliDCSValue.h"
21 #include "AliLog.h"
22 #include "AliTPCSensorTempArray.h"
23 #include "AliTPCDBPressure.h"
24
25 #include <TTimeStamp.h>
26
27 const Int_t kValCutTemp = 100;         // discard temperatures > 100 degrees
28 const Int_t kDiffCutTemp = 5;      // discard temperature differences > 5 degrees
29
30 //
31 // This class is the SHUTTLE preprocessor for the TPC detector.
32 // It contains several components, this far the part containing 
33 // temperatures is implemented
34 //
35
36 ClassImp(AliTPCPreprocessor)
37
38 //______________________________________________________________________________________________
39 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
40   AliPreprocessor("TPC",shuttle),
41   fTemp(0), fPressure(0), fConfigOK(kTRUE)
42 {
43   // constructor
44 }
45 //______________________________________________________________________________________________
46 // AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor& org) :
47 //   AliPreprocessor(org),
48 //   fTemp(0), fPressure(0), fConfigOK(kTRUE)
49 // {
50 //   // copy constructor not implemented
51 //   //   -- missing underlying copy constructor in AliPreprocessor 
52 // 
53 //   Fatal("AliTPCPreprocessor", "copy constructor not implemented");
54 //   
55 // //  fTemp = new AliTPCSensorTempArray(*(org.fTemp)); 
56 // }
57
58 //______________________________________________________________________________________________
59 AliTPCPreprocessor::~AliTPCPreprocessor()
60 {
61   // destructor
62   
63   delete fTemp;
64   delete fPressure;
65 }
66 //______________________________________________________________________________________________
67 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
68 {
69   Fatal("operator =", "assignment operator not implemented");
70   return *this;
71 }
72
73
74 //______________________________________________________________________________________________
75 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
76         UInt_t endTime)
77 {
78   // Creates AliTestDataDCS object
79
80   AliPreprocessor::Initialize(run, startTime, endTime);
81
82         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
83                 TTimeStamp(startTime).AsString(),
84                 TTimeStamp(endTime).AsString()));
85
86   // Temperature sensors
87
88         AliCDBEntry* entry = GetFromOCDB("Config", "Temperature"); 
89         TTree *confTree = (TTree*) entry->GetObject();
90         if ( confTree==0 ) {
91            AliError(Form("Temperature Config OCDB entry missing.\n"));
92            Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
93            fConfigOK = kFALSE;
94         }
95         fTemp = new AliTPCSensorTempArray(fStartTime, fEndTime, confTree);
96         fTemp->SetValCut(kValCutTemp);
97         fTemp->SetDiffCut(kDiffCutTemp);
98         confTree->Delete(); delete confTree; confTree=0;
99         entry->Delete(); delete entry; entry=0;
100   
101   // Pressure sensors
102  
103         entry = GetFromOCDB("Config", "Pressure"); 
104         confTree = (TTree*) entry->GetObject();
105         if ( confTree==0 ) {
106            AliError(Form("Pressure Config OCDB entry missing.\n"));
107            Log("AliTPCPreprocsessor: Pressure Config OCDB entry missing.\n");
108            fConfigOK = kFALSE;
109         }
110         fPressure = new AliDCSSensorArray(fStartTime, fEndTime, confTree);
111         confTree->Delete(); delete confTree; confTree=0;
112         entry->Delete(); delete entry; entry=0;
113
114 }
115
116 //______________________________________________________________________________________________
117 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
118 {
119   // Fills data into TPC calibrations objects
120
121    if (!dcsAliasMap) return 9;
122    if (!fConfigOK) return 9;
123
124   // Amanda servers provide information directly through dcsAliasMap
125
126   // Temperature sensors are processed by AliTPCCalTemp
127
128
129   UInt_t tempResult = MapTemperature(dcsAliasMap);
130   UInt_t result=tempResult;
131
132   // Pressure sensors
133
134   UInt_t pressureResult = MapPressure(dcsAliasMap);
135   result += pressureResult;
136   
137   // Other calibration information will be retrieved through FXS files
138   //  examples: 
139   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
140   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
141   //
142   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
143   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
144
145
146   return result;
147 }
148 //______________________________________________________________________________________________
149 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
150 {
151
152    // extract DCS temperature maps. Perform fits to save space
153
154   UInt_t result=0;
155   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
156   if (map) {
157     fTemp->MakeSplineFit(map);
158     AliInfo(Form("Temperature values extracted, fits performed.\n"));
159   } else {
160     AliError(Form("No temperature map extracted.\n"));
161     Log("AliTPCPreprocsessor: no temperature map extracted. \n");
162     result=9;
163   }
164   delete map;
165   // Now store the final CDB file
166   
167   if ( result == 0 ) { 
168         AliCDBMetaData metaData;
169         metaData.SetBeamPeriod(0);
170         metaData.SetResponsible("Haavard Helstrup");
171         metaData.SetComment("Preprocessor AliTPC data base entries.");
172
173         result = Store("Calib", "Temperature", fTemp, &metaData, 0, 0);
174         if ( result == 1 ) {                  
175           result = 0;
176         } else {
177           result = 1;
178         }                      // revert to new return code conventions
179    }
180
181    return result;
182 }
183 //______________________________________________________________________________________________
184 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
185 {
186
187    // extract DCS temperature maps. Perform fits to save space
188
189   UInt_t result=0;
190   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
191   if (map) {
192     fPressure->MakeSplineFit(map);
193     AliInfo(Form("Pressure values extracted, fits performed.\n"));
194   } else {
195     AliError(Form("No atmospheric pressure map extracted.\n"));
196     Log("AliTPCPreprocsessor: no atmospheric pressure map extracted. \n");
197     result=9;
198   }
199   delete map;
200   // Now store the final CDB file
201   
202   if ( result == 0 ) { 
203         AliCDBMetaData metaData;
204         metaData.SetBeamPeriod(0);
205         metaData.SetResponsible("Haavard Helstrup");
206         metaData.SetComment("Preprocessor AliTPC data base entries.");
207
208         result = Store("Calib", "Pressure", fPressure, &metaData, 0, 0);
209         if ( result == 1 ) {                  
210           result = 0;
211         } else {
212           result = 1;
213         }                      // revert to new return code conventions
214    }
215
216    return result;
217 }