/**************************************************** AliACORDEDataDCS class create to make a pointer to the ACORDE data DCS points Author: Pedro Gonzalez (CIEMAT, Madrid) ACORDE-DCS creator: Mario Ivan Martinez Hdez Last update: June 13th 2014 from Mario Rodríguez Cahuantzi (CINVESTAV, mrodriguez@fis.cinvestav.mx) ==> wrong ACORDE aliases for AMANDA fixed (https://alice.its.cern.ch/jira/browse/ALIROOT-5479) -- old alias: ACO_HV_MODULE[0..59]_VMON -- new alias: ACO_HV_MODULE[0..59]_INSIDE_VMON ()/ ACO_HV_MODULE[0..59]_OUTSIDE_VMON Last update: Fix of coding violations Mario Rodriguez C. (FCFM-BUAP) *****************************************************/ #include "AliACORDEDataDCS.h" #include "AliCDBMetaData.h" #include "AliDCSValue.h" #include "AliLog.h" #include #include #include #include #include #include #include #include ClassImp(AliACORDEDataDCS) //--------------------------------------------------------------- AliACORDEDataDCS::AliACORDEDataDCS(): TObject(), fRun(0), fStartTime(0), fEndTime(0), fGraphs("TGraph",kNGraphs), fFunc(0), fIsProcessed(kFALSE) { for(int i=0;iGetEntries()<2) { AliError(Form("Alias %s has just %d entries!", fAliasNames[j].Data(),aliasArr->GetEntries())); continue; } TIter iterarray(aliasArr); Double_t *time = new Double_t[aliasArr->GetEntries()]; Double_t *val = new Double_t[aliasArr->GetEntries()]; UInt_t ne=0; while ((aValue = (AliDCSValue*) iterarray.Next())) { val[ne] = aValue->GetFloat(); time[ne] = (Double_t) (aValue->GetTimeStamp()); fHv[j]->Fill(val[ne]); ne++; } CreateGraph(j, aliasArr->GetEntries(), time, val); delete[] val; delete[] time; } // calculate mean and rms of the first two histos for(int i=0;iGetMean(); fWidth[i] = fHv[i]->GetRMS(); } fIsProcessed=kTRUE; } //--------------------------------------------------------------- void AliACORDEDataDCS::Init() { // Init of AliACORDEDatDCS procedure // Loop over the aliases TH1::AddDirectory(kFALSE); fGraphs.SetOwner(1); TString aliasName; for(int i=0;iGetXaxis()->SetTitle("Hv"); } } //--------------------------------------------------------------- void AliACORDEDataDCS::Introduce(UInt_t numAlias, const TObjArray* aliasArr) { int entries=aliasArr->GetEntries(); AliInfo(Form("************ Alias: %s **********",fAliasNames[numAlias].Data())); AliInfo(Form(" %d DP values collected",entries)); } //--------------------------------------------------------------- void AliACORDEDataDCS::CreateGraph(int i, int dim, const Double_t *x, const Double_t *y) { // Create the plots for the ACORDE DCS TGraph *gr = new(fGraphs[fGraphs.GetEntriesFast()]) TGraph(dim, x, y); gr->GetXaxis()->SetTimeDisplay(1); gr->SetTitle(fAliasNames[i].Data()); AliInfo(Form("Array entries: %d",fGraphs.GetEntriesFast())); } //--------------------------------------------------------------- void AliACORDEDataDCS::Draw(const Option_t* /*option*/) { // Draw all histos and graphs if(!fIsProcessed) return; TString canvasHistoName; TCanvas *ch[10]; for (int i=0;i<10;i++) { canvasHistoName.Form("ACO_HV_MODULE"); ch[i]=new TCanvas(canvasHistoName,canvasHistoName,20,20,600,600); ch[i]->Divide(2,3); for(int j=0;j<6;j++) { ch[i]->cd(j+1); ((TGraph*) fGraphs.UncheckedAt(i*6+j))->SetMarkerStyle(20); ((TGraph*) fGraphs.UncheckedAt(i*6+j))->Draw("alp"); } } }