]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZERODataDCS.cxx
Fixes to make the preprocessor more robust in case of missing adc aliases.
[u/mrichter/AliRoot.git] / VZERO / AliVZERODataDCS.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
17 #include "AliVZERODataDCS.h"
18
19 #include "AliDCSValue.h"
20 #include "AliLog.h"
21
22 #include <TGraph.h>
23 #include <TAxis.h>
24 #include <TCanvas.h>
25 #include <TTimeStamp.h>
26 #include <TMap.h>
27 #include <TString.h>
28 #include <TObjString.h>
29 #include <TH1F.h>
30
31 class TH2;
32 class AliCDBMetaData;
33 class TDatime;
34
35 // AliVZERODataDCS class
36 // main aim to introduce the aliases for the VZERO DCS
37 // data points to be then
38 // stored in the OCDB, and to process them. 
39 // ProcessData() method called by VZEROPreprocessor
40
41 ClassImp(AliVZERODataDCS)
42
43 //_____________________________________________________________________________
44 AliVZERODataDCS::AliVZERODataDCS():
45         TObject(),
46         fRun(0),
47         fStartTime(0),
48         fEndTime(0),
49     fGraphs("TGraph",kNGraphs),
50         fFEEParameters(NULL),
51         fIsProcessed(kFALSE)
52
53 {
54   // Default constructor
55         for(int i=0;i<kNHvChannel;i++) {
56                 fDeadChannel[i] = kFALSE;
57                 fMeanHV[i]      = 100.0;
58                 fWidthHV[i]     = 0.0; 
59         }
60 }
61
62 //_____________________________________________________________________________
63 AliVZERODataDCS::AliVZERODataDCS(Int_t nRun, UInt_t startTime, UInt_t endTime):
64         TObject(),
65         fRun(nRun),
66         fStartTime(startTime),
67         fEndTime(endTime),
68         fGraphs("TGraph",kNGraphs),
69         fFEEParameters(new TMap()),
70         fIsProcessed(kFALSE)
71
72 {
73
74   // constructor with arguments
75         for(int i=0;i<kNHvChannel;i++) {
76          fDeadChannel[i] = kFALSE;        
77          fMeanHV[i]      = 100.0;
78      fWidthHV[i]     = 0.0; 
79         }
80         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", nRun,
81                 TTimeStamp(startTime).AsString(),
82                 TTimeStamp(endTime).AsString()));
83         
84         fFEEParameters->SetOwnerValue();
85         Init();
86
87 }
88
89 //_____________________________________________________________________________
90 AliVZERODataDCS::~AliVZERODataDCS() {
91
92   // destructor
93   fGraphs.Clear("C");
94   delete fFEEParameters;
95
96 }
97
98 //_____________________________________________________________________________
99 void AliVZERODataDCS::ProcessData(TMap& aliasMap){
100
101   // method to process the data
102
103   if(!(fAliasNames[0])) Init();
104
105   TObjArray *aliasArr;
106   AliDCSValue* aValue;
107
108   // starting loop on aliases
109   for(int iAlias=0; iAlias<kNAliases; iAlias++){
110
111     aliasArr = (TObjArray*) aliasMap.GetValue(fAliasNames[iAlias].Data());
112     if(!aliasArr){
113       AliError(Form("Alias %s not found!", fAliasNames[iAlias].Data()));
114       continue;
115     }
116
117     //Introduce(iAlias, aliasArr);
118     
119     if(aliasArr->GetEntries()<2){
120       AliWarning(Form("Alias %s has just %d entries!",
121                     fAliasNames[iAlias].Data(),aliasArr->GetEntries()));
122     }
123     
124     TIter iterarray(aliasArr);
125         
126     if(iAlias<kNHvChannel){ // Treating HV values
127         Int_t nentries = aliasArr->GetEntries();
128
129         Double_t *times = new Double_t[nentries];
130         Double_t *values = new Double_t[nentries];
131
132         UInt_t iValue=0;
133         Float_t variation = 0.0;
134
135         while((aValue = (AliDCSValue*) iterarray.Next())) {
136                         values[iValue] = aValue->GetFloat();
137                         if(iValue>0) {
138                                 if(values[iValue-1]>0.) variation = TMath::Abs(values[iValue]-values[iValue-1])/values[iValue-1];
139                                 if(variation > 0.10) fDeadChannel[GetOfflineChannel(iAlias)] = kTRUE;
140                         }
141                         times[iValue] = (Double_t) (aValue->GetTimeStamp());
142                         fHv[iAlias]->Fill(values[iValue]);
143                         printf("%s %f Dead=%d\n",fAliasNames[iAlias].Data(),values[iValue],fDeadChannel[GetOfflineChannel(iAlias)]);
144                         iValue++;
145         }      
146         CreateGraph(iAlias, aliasArr->GetEntries(), times, values); // fill graphs 
147
148         // calculate mean and rms of the first two histos
149         // and convert index to aliroot channel
150         Int_t iChannel     = GetOfflineChannel(iAlias); 
151         fMeanHV[iChannel]  = fHv[iAlias]->GetMean();
152         fWidthHV[iChannel] = fHv[iAlias]->GetRMS();
153
154         delete[] values;
155         delete[] times; 
156         } else { // Treating FEE Parameters
157                 AliDCSValue * lastVal = NULL;
158                 while((aValue = (AliDCSValue*) iterarray.Next())) lastVal = aValue; // Take only the last value
159                 fFEEParameters->Add(new TObjString(fAliasNames[iAlias].Data()),lastVal);
160         }      
161   }
162   
163   fIsProcessed=kTRUE;
164 }
165
166 //_____________________________________________________________________________
167 void AliVZERODataDCS::Init(){
168
169   // initialization of aliases and DCS data
170
171   TString sindex;
172   int iAlias = 0;
173   
174   for(int iSide = 0; iSide<2 ; iSide++){
175         for(int iRing = 0; iRing<4 ; iRing++){
176                 for(int iSector = 0; iSector<8 ; iSector++){
177                         if(iSide == 0) fAliasNames[iAlias] = "V00/HV/V0A/SECTOR";
178                         else fAliasNames[iAlias] = "V00/HV/V0C/SECTOR";
179                         sindex.Form("%d/RING%d",iSector,iRing);
180                         fAliasNames[iAlias] += sindex;
181                         
182                         fHv[iAlias] = new TH1F(fAliasNames[iAlias].Data(),fAliasNames[iAlias].Data(), 3000, kHvMin, kHvMax);
183                         fHv[iAlias]->GetXaxis()->SetTitle("Hv");
184                         iAlias++;
185                 }
186         }
187   }
188
189  // Time Resolution Parameters
190         
191         for(int iCIU = 0; iCIU<8 ; iCIU++){
192                 fAliasNames[iAlias++] = Form("V00/FEE/CIU%d/TimeResolution",iCIU);
193                 fAliasNames[iAlias++] = Form("V00/FEE/CIU%d/WidthResolution",iCIU);
194         }
195
196   if(iAlias!=kNAliases) 
197               AliError(Form("Number of DCS Aliases defined not correct"));
198
199 }
200
201 //_____________________________________________________________________________
202 void AliVZERODataDCS::Introduce(UInt_t numAlias, const TObjArray* aliasArr)const
203 {
204
205   // method to introduce new aliases
206
207   int entries=aliasArr->GetEntries();
208   AliInfo(Form("************ Alias: %s **********",fAliasNames[numAlias].Data()));
209   AliInfo(Form("        %d DP values collected",entries));
210
211 }
212
213 //_____________________________________________________________________________
214 void AliVZERODataDCS::CreateGraph(int i, int dim, const Double_t *x, const Double_t *y)
215 {
216
217    // Create graphics
218    
219    TGraph *gr = new(fGraphs[fGraphs.GetEntriesFast()]) TGraph(dim, x, y);
220
221    gr->GetXaxis()->SetTimeDisplay(1);
222    gr->SetTitle(fAliasNames[i].Data());
223
224    AliInfo(Form("Array entries: %d",fGraphs.GetEntriesFast()));
225
226 }
227
228
229 //_____________________________________________________________________________
230 void AliVZERODataDCS::Draw(const Option_t* /*option*/)
231 {
232 // Draw all histos and graphs
233
234   if(!fIsProcessed) return;
235
236   if(fGraphs.GetEntries()==0)  return;
237   
238   TString canvasName;
239   TCanvas *cHV[8];
240   
241   for(int iSide = 0 ;iSide<2;iSide++){
242         for(int iRing=0;iRing<4;iRing++){
243                 if(iSide == 0)  canvasName = "V0A_Ring";
244                 else  canvasName = "V0C_Ring";
245                 canvasName += iRing;
246                 int iCanvas = iSide*4 + iRing;
247                 cHV[iCanvas] = new TCanvas(canvasName,canvasName);
248                 cHV[iCanvas]->Divide(3,3);
249                 for(int iSector=0;iSector<8;iSector++){
250                         cHV[iCanvas]->cd(iSector+1);
251                         int iChannel = iSide*32 + iRing*8 + iSector; 
252                         ((TGraph*) fGraphs.UncheckedAt(iChannel))->SetMarkerStyle(20);
253                         ((TGraph*) fGraphs.UncheckedAt(iChannel))->Draw("ALP");
254
255                 }
256                                 
257         }
258   }
259
260 }
261