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