2 Contact: Jean-Luc Charvet <jean-luc.charvet@cern.ch>
3 Link: http://aliceinfo.cern.ch/static/Offline/dimuon/muon_html/README_mchda.html
4 Reference Runs: (station 3)
19 Number of events needed: 400 events for each calibration run (11)
20 Input Files: mutrkcalibvalues and config_ldc-MTRK-S3-0 in path : /afs/cern.ch/user/j/jcharvet/public/DA_validation
21 Output Files: local dir (not persistent) -> MUONTRKGAINda.par FXS -> run<#>_MCH_<ldc>_GAINS
25 /**************************************************************************
26 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
28 * Author: The ALICE Off-line Project. *
29 * Contributors are mentioned in the code where appropriate. *
31 * Permission to use, copy, modify and distribute this software and its *
32 * documentation strictly for non-commercial purposes is hereby granted *
33 * without fee, provided that the above copyright notice appears in all *
34 * copies and that both the copyright notice and this permission notice *
35 * appear in the supporting documentation. The authors make no claims *
36 * about the suitability of this software for any purpose. It is *
37 * provided "as is" without express or implied warranty. *
38 **************************************************************************/
43 -------------------------------------------------------------------------
44 2012-02-29 New version: MUONTRKGAINda.cxx,v 1.7
45 -------------------------------------------------------------------------
47 Version for MUONTRKGAINda MUON tracking
48 (A. Baldisseri, J.-L. Charvet)
51 Rem: AliMUON2DMap stores all channels, even those which are not connected
52 if pedMean == -1, channel not connected to a pad
63 #include <Riostream.h>
70 #include "AliMUONRawStreamTrackerHP.h"
71 #include "AliRawReader.h"
72 #include "AliMUONVStore.h"
73 #include "AliMUON2DMap.h"
74 #include "AliMUONCalibParamND.h"
75 #include "AliMpIntPair.h"
76 #include "AliMpConstants.h"
77 #include "AliRawDataErrorLog.h"
78 #include "AliMUONTrackerIO.h"
80 #include "AliMUONDspHeader.h"
89 #include "TStopwatch.h"
91 #include "TTimeStamp.h"
92 #include "TGraphErrors.h"
95 #include "TPluginManager.h"
97 #include "TObjString.h"
98 #include "THashTable.h"
99 #include <THashList.h>
107 #include "AliMUONPedestal.h"
108 #include "AliMUONGain.h"
109 #include "AliMUONErrorCounter.h"
112 int main(Int_t argc, const char** argv)
114 Int_t status=0 , status1=0 ;
118 const char* prefixDA = "MUONTRKGAINda"; // program prefix
119 const char* prefixLDC = getenv("DATE_ROLE_NAME"); // LDC name
120 if(prefixLDC == NULL) prefixLDC ="MCH" ;
121 printf("%s : -------- Begin execution : %s -------- \n",prefixLDC,prefixDA);
123 // const char* prefixDA = "MUONTRKGAINda"; // program prefix
124 // printf(" ######## Begin execution : %s ######## \n\n",prefixDA);
127 // decode the input line
128 if (argc!=2) { printf("Wrong number of arguments\n"); return -1; }
132 // needed for streamer application
133 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
138 // needed for Minuit plugin
139 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
143 "TMinuitMinimizer(const char*)");
144 Int_t skipEvents = 0;
145 Int_t maxEvents = 1000000;
146 Int_t maxDateEvents = 1000000;
148 Int_t nDateEvents = 0;
149 Int_t nDateRejEvents = 0;
150 Int_t nGlitchErrors= 0;
151 Int_t nParityErrors= 0;
152 Int_t nPaddingErrors= 0;
153 Int_t nTokenlostErrors= 0;
154 Int_t nEventsRecovered = 0;
156 Int_t nEvthreshold = 50; //below this nb_evt the mean of the charge is not calculated and forced to 4085 (sigma)
157 Int_t statusDA = 0 ; // DA return code
159 TString logOutputFile;
161 Char_t flatFile[256]="";
164 UInt_t runNumber = 0;
169 Int_t nEntries = daqDA_ECS_getTotalIteration(); // usually = 11 = Nb of calibration runs
170 Int_t nInit=0; // = 0 all DAC values ; = 1 DAC=0 excluded (default=1)
171 Int_t nbpf1=4; // nb of points for linear fit (default=6)
172 Int_t printLevel = 0; // printout (default=0, =1 =>.ped , => .peak & .param)
173 Int_t plotLevel = 1; // plotout (default=1 => tree , =2 tree+Tgraph+fit)
177 // Reading current iteration
178 nIndex = daqDA_ECS_getCurrentIteration();
179 if(nIndex<0 || nIndex>nEntries) {printf("\n Failed: nIndex = %d\n",nIndex); return -1 ;}
182 AliMUONGain* muonGain = new AliMUONGain();
183 muonGain->SetprefixDA(prefixDA);
184 muonGain->SetprefixLDC(prefixLDC);
185 muonGain->SetAliRootDataFileName(); // MUONTRKGAINda_data.root
186 muonGain->SetStatusDA(statusDA);
188 // Output log file initialisations
189 sprintf(flatFile,"%s.log",prefixDA);
190 logOutputFile=flatFile;
191 AliLog::SetStreamOutput(&filcout); // Print details on logfile
192 filcout.open(logOutputFile.Data());
193 filcout<<"//=================================================" << endl;
194 filcout<<"//" << prefixLDC << " " << prefixDA << endl;
195 filcout<<"//=================================================" << endl;
196 filcout<<"// * Date : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
198 muonGain->SetAlifilcout(&filcout);
199 cout<<prefixLDC << " : Date: " << muonGain->GetDate()->AsString("l") << "\n" << endl;
207 // DAC values stored in array vDAC (reading dbfile in DETDB)
208 // Int_t vDAC[11] = {0,200,400,800,1200,1600,2000,2500,3000,3500,4000}
209 Int_t nConfig = 1; // flag to read or not configuration ascii file in detDB
210 Int_t vDAC[11]; // DAC values
211 Char_t dbfile[256]="";
214 sprintf(dbfile,"mutrkcalibvalues");
215 status=daqDA_DB_getFile(dbfile,dbfile);
216 if(status) {printf(" Failed : input file %s is missing, status = %d\n",dbfile,status); return -1; }
218 ifstream filein(dbfile,ios::in); Int_t k=0, kk;
219 while (k<nEntries ) { filein >> kk >> vDAC[k] ; k++; }
220 injCharge=vDAC[nIndex-1];
222 filein >> nInit >> line ; cout << "mutrkcalibvalues: " << line << "=" << nInit << " " ; // = 0 all DAC values fitted ; = 1 DAC=0 excluded (default=1)
223 filein >> nbpf1 >> line ; cout << line << "=" << nbpf1 << " " ; // nb of points for linear fit (default=6)
224 filein >> printLevel >> line; cout << line << "=" << printLevel << " " ; // printout (default=0, =1 =>.ped /run, =2 => .peak & .param)
225 filein >> plotLevel >> line; cout << line << "=" << plotLevel << " " ; // plotout (default=1 => tree , =2 tree+Tgraph+fit)
226 filein >> nConfig >> line; cout << line << "=" << nConfig << " " ; //nConfig (default=1 => read config in DetDB, otherwise =0)
227 filein >> nEvthres >> line ;
228 if(nEvthres !=0)nEvthreshold=nEvthres; cout << line << "=" << nEvthreshold << " " ; // (default = 0 <=> 50) below nEvthreshold calibration not performed
229 filein >> nbev >> line; // Nb of events to read (default = 0 => reading all events)
230 if(nbev !=0){maxEvents=nbev; cout << line << "=" << maxEvents << " " ;}
233 muonGain->SetAliPrintLevel(printLevel);
234 muonGain->SetAliPlotLevel(plotLevel);
235 muonGain->SetconfigDA(nConfig);
236 muonGain->SetnEvthreshold(nEvthreshold);
237 // muonGain->SetStatusDA(statusDA);
241 // Laod configuration ascii file from DetDB
242 sprintf(dbfile,"config_%s",getenv("DATE_ROLE_NAME"));
243 status=daqDA_DB_getFile(dbfile,dbfile);
244 if(status) {printf(" Failed : Configuration file %s is missing, status = %d\n",dbfile,status); return -1; }
245 // else printf(" *** Copy ascii config file: %s from DetDB to working directory and reading ...*** \n",dbfile);
246 muonGain->LoadConfig(dbfile);
249 // Rawdeader, RawStreamHP
250 AliRawReader* rawReader = AliRawReader::Create(inputFile.Data());
251 AliMUONRawStreamTrackerHP* rawStream = new AliMUONRawStreamTrackerHP(rawReader);
252 // rawStream->DisableWarnings();
253 rawStream->EnabbleErrorLogger();
255 // kLowErrorDetail, /// Logs minimal information in the error messages.
256 // kMediumErrorDetail, /// Logs a medium level of detail in the error messages.
257 // kHighErrorDetail /// Logs maximum information in the error messages.
258 // rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kLowErrorDetail);
259 rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kMediumErrorDetail);
260 // rawStream->SetLoggingDetailLevel(AliMUONRawStreamTrackerHP::kHighErrorDetail);
262 printf("\n%s : Reading data from file %s\n",prefixLDC,inputFile.Data());
264 Int_t tabTokenError[20][14];
265 for ( Int_t i=0 ; i<20 ; i++) { for ( Int_t j=0 ; j<14 ; j++) { tabTokenError[i][j]=0;} }
267 while (rawReader->NextEvent())
269 Int_t eventType = rawReader->GetType();
270 runNumber = rawReader->GetRunNumber();
271 if(nDateEvents==0) { filcout<<"// ----> RUN = " << runNumber << "\n" << endl;}
272 if (nDateEvents >= maxDateEvents) break;
273 if (nEvents >= maxEvents) break;
274 if (nDateEvents>0 && nDateEvents % 100 == 0)
275 cout<< prefixLDC << " : DATE events = " << nDateEvents << " Used events = " << nEvents << endl;
277 // check shutdown condition
278 if (daqDA_checkShutdown())
283 rawReader->NextEvent();
287 if (eventType != PHYSICS_EVENT)
288 continue; // for the moment
290 const char* detail = "";
291 // First lopp over DDL's to find good events
292 // Error counters per event (counters in the decoding lib are for each DDL)
293 Bool_t eventIsErrorMessage = kFALSE;
294 int eventGlitchErrors = 0;
295 int eventParityErrors = 0;
296 int eventPaddingErrors = 0;
297 int eventTokenlostErrors = 0;
301 if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
302 eventGlitchErrors += rawStream->GetGlitchErrors();
303 eventParityErrors += rawStream->GetParityErrors();
304 eventPaddingErrors += rawStream->GetPaddingErrors();
305 eventTokenlostErrors += rawStream->GetTokenLostErrors();
306 if (rawStream->GetTokenLostErrors())
309 const AliMUONRawStreamTrackerHP::AliBlockHeader* blkHeader = 0x0;
310 const AliMUONRawStreamTrackerHP::AliDspHeader* dspHeader = 0x0;
311 Int_t nBlock = rawStream->GetBlockCount();
312 for(Int_t iBlock = 0; iBlock < nBlock ;iBlock++)
314 blkHeader = rawStream->GetBlockHeader(iBlock);
315 // printf("Block %d Total length %d\n",iBlock,blkHeader->GetTotalLength());
316 Int_t nDsp = rawStream->GetDspCount(iBlock);
317 // printf("Block %d DSP %d\n",iBlock,nDsp);
318 for(Int_t iDsp = 0; iDsp < nDsp ;iDsp++)
320 dspHeader = blkHeader->GetDspHeader(iDsp);
321 // printf("Dsp %d Add %X\n",iDsp,dspHeader);
322 if (dspHeader->GetErrorWord())
324 Int_t ddl = rawStream->GetDDL() ;
325 // Int_t ddl = AliDAQ::DdlID("MUONTRK", rawStream->GetDDL()) - 2560 ; // format 2560 + ddl
326 Int_t frt = (dspHeader->GetErrorWord() & 0xFFFF0000) >> 16 ; // 4*4bits right shift
327 tabTokenError[ddl][frt]++;
328 // printf(" DDL %d error word %X %d %d\n",ddl,dspHeader->GetErrorWord(),frt,tabTokenError[8][4]);
334 } while(rawStream->NextDDL());
336 AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
337 if (!eventIsErrorMessage)
339 // Good events (no error) -> compute pedestal for all channels
342 muonGain->SetAliNCurrentEvents(nEvents);
343 while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next()))
345 for(int i = 0; i < busPatch->GetLength(); ++i)
347 busPatch->GetData(i, manuId, channelId, charge);
348 muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
354 // Events with errors
355 if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
357 filcout << " ----------- Date Event recovered = " << nDateEvents << " ----------------" << endl;
358 // Recover parity errors -> compute pedestal for all good buspatches
359 if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
362 filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
363 <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
364 <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;
368 filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
369 <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
370 <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
374 muonGain->SetAliNCurrentEvents(nEvents);
375 while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next()))
377 // Check the buspatch -> if error not use it in the pedestal calculation
379 for(int i = 0; i < busPatch->GetLength(); ++i)
381 if (!busPatch->IsParityOk(i)) errorCount++;
386 for(int i = 0; i < busPatch->GetLength(); ++i)
388 busPatch->GetData(i, manuId, channelId, charge);
389 muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
395 AliMUONErrorCounter* errorCounter;
396 // Bad buspatch -> not used (just print)
397 filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
398 <<" parity errors "<<errorCount<<endl;
399 // Number of events where this buspatch is missing
400 sprintf(bpname,"bp%d",busPatch->GetBusPatchId());
401 if (!(errorCounter = (AliMUONErrorCounter*) (muonGain->GetErrorBuspatchTable()->FindObject(bpname))))
404 errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
405 errorCounter->SetName(bpname);
406 muonGain->GetErrorBuspatchTable()->Add(errorCounter);
411 errorCounter->Increment();
413 // errorCounter->Print();
414 } // end of if (!errorCount)
415 } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
418 } //end of if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
421 // Fatal errors reject the event
422 detail = Form(" ----------- Date Event rejected = %d ----------------",nDateEvents);
424 filcout << detail << endl;
425 if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
428 filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
429 <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
430 <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;
434 filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
435 <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
436 <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
439 } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
440 filcout<<"Number of errors : Glitch "<<eventGlitchErrors
441 <<" Parity "<<eventParityErrors
442 <<" Padding "<<eventPaddingErrors
443 <<" Token lost "<<eventTokenlostErrors<<endl;
445 } // end of if (!rawStream->IsErrorMessage())
447 if (eventGlitchErrors) nGlitchErrors++;
448 if (eventParityErrors) nParityErrors++;
449 if (eventPaddingErrors) nPaddingErrors++;
450 //if (eventTokenlostErrors) nTokenlostErrors++;
452 } // while (rawReader->NextEvent())
456 // process and store mean peak values in .root file
457 sprintf(flatFile,"%s.par",prefixDA);
458 if(shuttleFile.IsNull())shuttleFile=flatFile;
459 injCharge=vDAC[nIndex-1];
460 muonGain->SetAliIndex(nIndex); // fIndex
461 muonGain->SetAliInjCharge(injCharge);
462 muonGain->SetAliNEvents(nEvents);
463 muonGain->SetAliRunNumber(runNumber);
464 muonGain->MakePedStoreForGain(shuttleFile);
466 // writing some counters
467 detail=Form("\n%s : Nb of DATE events = %d",prefixLDC,nDateEvents) ; cout << detail; filcout << detail ;
468 detail=Form("\n%s : Nb of Glitch errors = %d",prefixLDC,nGlitchErrors) ; cout << detail; filcout << detail ;
469 detail=Form("\n%s : Nb of Parity errors = %d",prefixLDC,nParityErrors) ; cout << detail; filcout << detail ;
470 detail=Form("\n%s : Nb of Token lost errors = %d",prefixLDC,nTokenlostErrors) ; cout << detail; filcout << detail ;
471 detail=Form("\n%s : Nb of Rejected DATE events = %d",prefixLDC,nDateRejEvents) ; cout << detail; filcout << detail ;
472 detail=Form("\n%s : Nb of recovered events = %d",prefixLDC,nEventsRecovered) ; cout << detail; filcout << detail ;
473 detail=Form("\n%s : Nb of events without errors = %d",prefixLDC,nEvents-nEventsRecovered) ; cout << detail; filcout << detail ;
474 detail=Form("\n%s : Nb of used events = %d (threshold= %d)\n\n",prefixLDC,nEvents,nEvthreshold); cout << detail; filcout << detail ;
476 // Writing Token Error table
479 char* detail=Form("\nWarning: Token Lost occurence \n");
482 for ( Int_t i=0 ; i<20 ; i++)
484 for ( Int_t j=4 ; j<14 ; j++)
486 if(tabTokenError[i][j]>0)
488 Int_t tab=tabTokenError[i][j];
490 Int_t station = i/4 +1;
491 if( j % 2 == 0)detail=Form(" in DDL= %d (station %d) and FRT%d ( Up ) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
492 else detail=Form(" in DDL= %d (station %d) and FRT%d (Down) => %d Token errors (address = 0x%X0000)",2560+i,station,frt,tab,j);
493 printf("%s\n",detail);
494 filcout << detail << endl;
503 muonGain->SetAliInit(nInit); // fnInit
504 muonGain->SetAliEntries(nEntries); // fnEntries
505 muonGain->SetAliNbpf1(nbpf1); // fnbpf1
506 muonGain->MakeGainStore(shuttleFile);
507 status = muonGain->GetStatusDA() ;
511 detail=Form("%s : Root data file : %s\n",prefixLDC,muonGain->GetRootDataFileName()); filcout << detail ; // cout << detail;
512 detail=Form("%s : Output logfile : %s\n",prefixLDC,logOutputFile.Data()); filcout << detail ; // cout << detail;
513 detail=Form("%s : Gain Histo file : %s\n",prefixLDC,muonGain->GetHistoFileName()); filcout << detail ; // cout << detail;
514 detail=Form("%s : Gain file (to SHUTTLE) : %s\n",prefixLDC,shuttleFile.Data()); filcout << detail ; // cout << detail;
516 // Copying files to local DB folder defined by DAQ_DETDB_LOCAL
518 unsigned int nLastVersions=50;
519 dir= getenv("DAQ_DETDB_LOCAL");
521 unsigned int nLastVersions=50;
522 printf("\n%s : *** Local DataBase: %s (Max= %d) ***\n",prefixLDC,dir,nLastVersions);
523 status1 = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
527 status1 = daqDA_localDB_storeFile(muonGain->GetRootDataFileName(),nLastVersions);
528 status1 = daqDA_localDB_storeFile(muonGain->GetHistoFileName(),nLastVersions);
529 status1 = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions); // if(status1)printf(" Store file : %s status = %d\n",shuttleFile.Data(),status1);
534 // Transferring pedestal file to FES (be sure that env variable DAQDALIB_PATH is set)
536 status1 = daqDA_FES_storeFile(shuttleFile.Data(),"GAINS");
537 if (status1) { detail=Form("%s: !!! ERROR: Failed to export calibration file : %s to FES \n",prefixLDC,shuttleFile.Data());
538 printf("%s",detail); filcout << detail ; status= -1; }
539 // else { detail=Form("%s : ---- STORE calibration FILE in FES : OK ---- \n",prefixLDC); printf("%s",detail); filcout << detail ;}
544 std::ifstream in(shuttleFile.Data());
545 ostringstream stringout;
547 while ( in.getline(line,1024) )
548 stringout << line << "\n";
551 amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
552 TObjString gaindata(stringout.str().c_str());
553 Int_t amoreStatus = amoreDA.Send("Gains",&gaindata);
555 cout << prefixLDC << " : !!! ERROR: Failed to write Gains in the AMORE database : " << amoreStatus << endl ; status=-1 ;
557 cout << prefixLDC << " : amoreDA.Send(Gains) ok" << endl;
559 cout << prefixLDC << " : Warning: MCH DA not compiled with AMORE support" << endl;
563 if(!status)printf("\n%s : -------- End execution : %s -------- (status= %d) \n",prefixLDC,prefixDA,status);
564 else { printf("\n%s : -------- %s ending in ERROR !!!! -------- (status= %d) \n",prefixLDC,prefixDA,status);}
567 printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());