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