added prediction processor framework and implementations; moved MUON preprocessors
[u/mrichter/AliRoot.git] / HLT / pendolino / TPC / AliHLTPredictionProcessorTPC.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Haavard Helstrup                                      *
8 //*                  for The ALICE HLT 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   AliHLTPredictionProcessorTPC.cxx
20     @author Haavard Helstrup
21     @date   
22     @brief  
23 */
24
25 #include "AliHLTPredictionProcessorTPC.h"
26
27 #include "AliHLTDCSArray.h"
28 #include "TROOT.h"
29 #include "TArrayF.h"
30 #include "TObjArray.h"
31 #include "TTree.h"
32 #include "TString.h"
33 #include "AliCDBEntry.h"
34 #include "AliCDBMetaData.h"
35 #include "AliDCSValue.h"
36
37 #include <TTimeStamp.h>
38 #include <TObjString.h>
39
40
41 ClassImp(AliHLTPredictionProcessorTPC)
42
43 const TString kAmandaString = "TPC_PT_%d_TEMPERATURE";
44
45 //______________________________________________________________________________________________
46
47 AliHLTPredictionProcessorTPC::AliHLTPredictionProcessorTPC(
48            const char* detector, AliHLTPendolino* pendolino) :
49                AliHLTPredictionProcessorInterface(detector, pendolino),
50                fConfigOK(kTRUE),
51                fConfTreeTemp(0),
52                fTemp(0),
53                fPredict(kFALSE),
54                fRun(0),
55                fStartTime(0),
56                fEndTime(0)
57 {
58  // constructor
59 }
60
61 //______________________________________________________________________________________________
62
63 AliHLTPredictionProcessorTPC::~AliHLTPredictionProcessorTPC()
64 {
65
66   // destructor
67         if (fConfTreeTemp) {
68                 delete fConfTreeTemp;
69         }
70         if (fTemp) {
71                 fTemp->Delete();
72                 delete fTemp;
73         }
74
75 }
76
77 //______________________________________________________________________________________________
78
79 UInt_t AliHLTPredictionProcessorTPC::makePrediction(Bool_t doPrediction) {
80 //
81 // Signal whether prediction making is required
82 //
83    fPredict = doPrediction; // own member indicating that prediction making is required
84    return 0;
85 }
86
87 //______________________________________________________________________________________________
88
89 void AliHLTPredictionProcessorTPC::Initialize(Int_t run, UInt_t startTime,
90            UInt_t endTime) {
91 //
92 //  Initialise TPC prediction processor. Read config entry from HCDB
93 //
94    fRun = run;               // storing locally the run number
95    fStartTime = startTime;   // storing locally the start time
96    fEndTime = endTime;       // storing locally the end time
97
98
99    AliCDBEntry* entry = GetFromOCDB("Config", "HLTTemperature");
100    if (entry) fConfTreeTemp = (TTree*) entry->GetObject();
101    if ( fConfTreeTemp==0 ) {
102      Log("AliHLTPredictionProcessorTPC: Temperature Config HCDB entry missing.\n");
103      fConfigOK = kFALSE;
104      return;
105    }
106 }
107
108 //______________________________________________________________________________________________
109
110 UInt_t AliHLTPredictionProcessorTPC::Process(TMap* dcsAliasMap) {
111 //
112 // Event-by-event processing of the PredictionProcessorTPC
113 //
114
115 // test if configuration available and DCS map valid
116
117   if (!dcsAliasMap) return 91;
118   if (dcsAliasMap->GetEntries() == 0 ) return 92;
119   if (!fConfigOK) return 93;
120
121 //
122 // Extract values from DCS maps. Store updated entry to HCDB
123 //
124    UInt_t retVal = 0;
125    Int_t start = 0;
126    TString path2("Calib");        // for the storage path in HCDB
127    TString path3("Temperature");    // for the storage path in HCDB
128
129    //get data from the HCDB
130    AliCDBEntry* entry = GetFromOCDB(path2.Data(), path3.Data());
131    if (entry != 0) {
132        entry->PrintMetaData();  // use the AliCDBEntry
133        fTemp = (TObjArray*)entry->GetObject();
134    } else {
135        Log("Cannot retrieve HCDB entry. New TObjArray generated");
136        fTemp = new TObjArray();
137    }
138
139    // process temperatures
140    
141    UInt_t tempResult = ExtractTemperature(dcsAliasMap); 
142
143    //process dcsAliasMap to ROOT object
144    // and create AliCDBEntry 
145
146    if (tempResult == 0) {
147                 // create meta data entry for HCDB
148         TString comment("HLT temperatures");
149         AliCDBMetaData meta(this->GetName(), 0, "unknownAliRoot", comment.Data());
150
151         // store AliCDBEntry in HCDB
152         if (Store(path2.Data(), path3.Data(), fTemp, &meta, start, kTRUE)) {
153                 Log(" +++ Successfully stored object ;-)");
154         } else {
155                 Log(" *** Storing of OBJECT failed!!");
156                 retVal = 7;
157         }
158     } else {
159         retVal = 94;
160     }
161    
162    return retVal;
163 }
164
165 //______________________________________________________________________________________________
166
167 UInt_t AliHLTPredictionProcessorTPC::ExtractTemperature(TMap* dcsAliasMap)
168 {
169 // Extract temperature values from DCS maps, according to infoamtion from configuration tree
170
171   const Int_t error = 9; 
172   Int_t nentries = fConfTreeTemp->GetEntries();
173   if (nentries<1) return error;
174     
175   TString stringId;
176   Int_t echa=0;
177   fConfTreeTemp->SetBranchAddress("ECha",&echa);
178   fConfTreeTemp->GetEntry(0);
179   stringId = Form(kAmandaString.Data(),echa);
180   UInt_t time = GetCurrentTime(dcsAliasMap,stringId.Data());   
181   AliHLTDCSArray* temperatures = new AliHLTDCSArray(nentries);
182   temperatures->SetTime(time);
183   
184   for (Int_t isensor=0; isensor<nentries; isensor++ ){
185      fConfTreeTemp->GetEntry(isensor);
186      stringId = Form(kAmandaString.Data(),echa);
187      Float_t temp = GetSensorValue(dcsAliasMap,stringId.Data());
188      temperatures->SetValue(isensor,temp);
189   }
190   
191   fTemp->Add(temperatures);
192   return 0;  
193 }
194
195   
196 //______________________________________________________________________________________________
197
198 Float_t AliHLTPredictionProcessorTPC::GetSensorValue(TMap* dcsAliasMap,
199                                      const char* stringId)
200 {
201   // return last value read from sensor specified by stringId
202   Float_t sensorValue=0;
203   TObjArray* valueSet;
204   TPair* pair = (TPair*)dcsAliasMap->FindObject(stringId);
205   if (pair) {
206      valueSet = (TObjArray*)pair->Value();
207          if (valueSet) {
208         Int_t nentriesDCS = valueSet->GetEntriesFast();
209         AliDCSValue *val = (AliDCSValue *)valueSet->At(nentriesDCS - 1);
210         if (val) {
211                         sensorValue=val->GetFloat();
212                 }
213          }
214   }
215   return sensorValue;
216 }
217         
218 //______________________________________________________________________________________________
219
220 UInt_t AliHLTPredictionProcessorTPC::GetCurrentTime(TMap* dcsAliasMap,
221                                      const char* stringId)
222 {
223   // return last time value read from sensor specified by stringId
224   
225   UInt_t time=0;
226   TObjArray* valueSet;
227   TPair* pair = (TPair*) (dcsAliasMap->FindObject(stringId));
228   if (pair) {
229      valueSet = (TObjArray*) (pair->Value());
230          if (valueSet) {
231         Int_t nentriesDCS = valueSet->GetEntriesFast();
232         AliDCSValue *val = (AliDCSValue *) (valueSet->At(nentriesDCS - 1));
233                 if (val) {
234                 time = val->GetTimeStamp();
235                 }
236          }
237   }
238   return time;
239 }
240
241 TMap* AliHLTPredictionProcessorTPC::produceTestData(TString aliasName) {
242     TMap* resultMap = 0;
243
244     // here has to come real dummy data :-)
245     resultMap = new TMap();
246     TTimeStamp tt;
247     TObjString* name = new TObjString("DummyData");
248         Float_t fval = 33.3;
249     AliDCSValue* val = new AliDCSValue(fval, tt.GetTime());
250     TObjArray* arr = new TObjArray();
251     arr->Add(val);
252     resultMap->Add(name, arr);
253
254     return resultMap;
255 }
256
257                              
258