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