/* - Contact: - prino@to.infn.it - Link: - http://www.to.infn.it/~prino/alice/RawData/run11161.date - Run Type: - PEDESTAL_RUN - DA Type: - LDC - Number of events needed: > 10 - Input Files: - - Output Files: - SDDbase_step1_ddl*c*_sid*.data SDDbase_step2_ddl*c*_sid*.data - Trigger types used: */ ////////////////////////////////////////////////////////////////////////////// // Detector Algorithm for analysis of SDD baseline runs. // // // // Produces ASCII and ROOT output files with: // // 1. anode quality bit // // 1. Baseline values // // 2. Raw noise // // 3. Common mode coefficients // // 4. Noise corrected for common mode // // Files are stored locally. // // The next DAQ-DA step on Test Pulse run writes files to FXS // // // // Author: F. Prino (prino@to.infn.it) // // // ////////////////////////////////////////////////////////////////////////////// extern "C" { #include "daqDA.h" } #include "event.h" #include "monitor.h" #include #include // ROOT includes #include #include #include #include #include #include #include #include #include // AliRoot includes #include "AliRawReaderDate.h" #include "AliITSOnlineSDDBase.h" #include "AliITSOnlineSDDCMN.h" #include "AliITSRawStreamSDD.h" #include "AliITSRawStreamSDDCompressed.h" #ifdef ALI_AMORE #include #endif /* Main routine Arguments: list of DATE raw data files */ int main(int argc, char **argv) { int status = 0; // line added to solve IO problems gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", "*", "TStreamerInfo", "RIO", "TStreamerInfo()"); /* log start of process */ printf("ITS SDD BASELINE+NOISE algorithm program started\n"); /* check that we got some arguments = list of files */ if (argc<2) { printf("Wrong number of arguments\n"); return -1; } Int_t maxNEvents=10; // maximum number of events to be analyzed const Int_t kTotDDL=24; const Int_t kModPerDDL=12; const Int_t kSides=2; UInt_t amSamplFreq=40; gSystem->Exec("rm -f SDDbase_*.data"); gSystem->Exec("rm -f SDDbase_step2_LDC.tar"); AliITSOnlineSDDBase **base=new AliITSOnlineSDDBase*[kTotDDL*kModPerDDL*kSides]; AliITSOnlineSDDCMN **corr=new AliITSOnlineSDDCMN*[kTotDDL*kModPerDDL*kSides]; TH2F **histo=new TH2F*[kTotDDL*kModPerDDL*kSides]; Bool_t isFilled[kTotDDL*kModPerDDL*kSides]; Bool_t writtenoutput=kFALSE; Char_t hisnam[20]; for(Int_t iddl=0; iddleventType; switch (event->eventType){ /* START OF RUN */ case START_OF_RUN: break; /* END OF RUN */ case END_OF_RUN: break; // case PHYSICS_EVENT: // comment this line for test raw data // break; // comment this line for test raw data case CALIBRATION_EVENT: break; // uncomment this line for test raw data case PHYSICS_EVENT: // uncomment this line for test raw data printf(" event number = %i \n",iev); ievPed++; AliRawReader *rawReader = new AliRawReaderDate((void*)event); rawReader->Reset(); UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rawReader); amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr); AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader,cdhAttr); if(!writtenoutput){ printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq); writtenoutput=kTRUE; } for(Int_t iddl=0; iddlReset(); } } } while(s->Next()){ Int_t iDDL=rawReader->GetDDLID(); Int_t iCarlos=s->GetCarlosId(); if(s->IsCompletedModule()) continue; if(s->IsCompletedDDL()) continue; if(iDDL>=0 && iDDLGetChannel(); histo[index]->Fill(s->GetCoord2(),s->GetCoord1(),s->GetSignal()); isFilled[index]=1; } } delete s; delete rawReader; for(Int_t iddl=0; iddlSetLastGoodTB(126); else base[index]->SetLastGoodTB(254); base[index]->AddEvent(histo[index]); } if(iStep==1){ if(amSamplFreq==20) corr[index]->SetLastGoodTB(126); else corr[index]->SetLastGoodTB(254); corr[index]->AddEvent(histo[index]); } } } } /* free resources */ free(event); } if(ievPed>=maxNEvents) break; } } if(iStep==0){ for(Int_t iddl=0; iddlValidateAnodes(); base[index]->WriteToASCII(); } } } } } /* write report */ TDatime time; TObjString timeinfo(Form("%02d%02d%02d%02d%02d%02d",time.GetYear()-2000,time.GetMonth(),time.GetDay(),time.GetHour(),time.GetMinute(),time.GetSecond())); printf("Run #%s, received %d calibration events, time %s\n",getenv("DATE_RUN_NUMBER"),ievPed,timeinfo.GetString().Data()); /* report progress */ daqDA_progressReport(90); TObjArray* basHistos=new TObjArray(); TObjArray* noiseHistos=new TObjArray(); TObjArray* cmnHistos=new TObjArray(); TObjArray* corrnHistos=new TObjArray(); TObjArray* statusHistos=new TObjArray(); Char_t filnam[100],command[150]; for(Int_t iddl=0; iddlValidateAnodes(); corr[index]->WriteToASCII(); if(isFilled[index]){ basHistos->AddLast(corr[index]->GetBaselineAnodeHisto()); noiseHistos->AddLast(corr[index]->GetRawNoiseAnodeHisto()); cmnHistos->AddLast(corr[index]->GetCMNCoefAnodeHisto()); corrnHistos->AddLast(corr[index]->GetCorrNoiseAnodeHisto()); statusHistos->AddLast(corr[index]->GetStatusAnodeHisto()); sprintf(filnam,"SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,imod,isid); sprintf(command,"tar -rf SDDbase_step2_LDC.tar %s",filnam); gSystem->Exec(command); } } } } #ifdef ALI_AMORE amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender); Int_t statusamore =0; statusamore += amoreDA.Send("TimeInfoPedestal",&timeinfo); statusamore += amoreDA.Send("Baselines",basHistos); statusamore += amoreDA.Send("RawNoise",noiseHistos); statusamore += amoreDA.Send("CommonMode",cmnHistos); statusamore += amoreDA.Send("CorrectedNoise",corrnHistos); statusamore += amoreDA.Send("NoisyChannels",statusHistos); if ( statusamore ) printf("Warning: Failed to write Arrays in the AMORE database\n"); else printf("amoreDA.Send() OK\n"); #else printf("Warning: SDDBAS DA not compiled with AMORE support\n"); #endif TFile *fh=new TFile("SDDbaseHistos.root","RECREATE"); basHistos->Write(); noiseHistos->Write(); cmnHistos->Write(); corrnHistos->Write(); statusHistos->Write(); fh->Close(); /* report progress */ daqDA_progressReport(100); return status; }