2 Contact: Jean-Luc Charvet <jean-luc.charvet@cea.fr>
3 Link: http://aliceinfo.cern.ch/static/Offline/dimuon/muon_html/README_Mchda
6 Number of events needed: 400 events for each calibration run
7 Input Files: Rawdata file (DATE format)
8 Output Files: local dir (not persistent) -> MUONTRKGAINda_<run#>.par
9 FXS -> run<#>_MCH_<ldc>_GAINS
13 /**************************************************************************
14 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
16 * Author: The ALICE Off-line Project. *
17 * Contributors are mentioned in the code where appropriate. *
19 * Permission to use, copy, modify and distribute this software and its *
20 * documentation strictly for non-commercial purposes is hereby granted *
21 * without fee, provided that the above copyright notice appears in all *
22 * copies and that both the copyright notice and this permission notice *
23 * appear in the supporting documentation. The authors make no claims *
24 * about the suitability of this software for any purpose. It is *
25 * provided "as is" without express or implied warranty. *
26 **************************************************************************/
31 -------------------------------------------------------------------------
32 2009-05-18 New version: MUONTRKGAINda.cxx,v 1.1
33 -------------------------------------------------------------------------
35 Version for MUONTRKGAINda MUON tracking
36 (A. Baldisseri, J.-L. Charvet)
39 Rem: AliMUON2DMap stores all channels, even those which are not connected
40 if pedMean == -1, channel not connected to a pad
51 #include <Riostream.h>
58 #include "AliMUONRawStreamTrackerHP.h"
59 #include "AliRawReader.h"
60 #include "AliMUONVStore.h"
61 #include "AliMUON2DMap.h"
62 #include "AliMUONCalibParamND.h"
63 #include "AliMpIntPair.h"
64 #include "AliMpConstants.h"
65 #include "AliRawDataErrorLog.h"
66 #include "AliMUONTrackerIO.h"
74 #include "TStopwatch.h"
76 #include "TTimeStamp.h"
77 #include "TGraphErrors.h"
80 #include "TPluginManager.h"
82 #include "TObjString.h"
83 #include "THashTable.h"
84 #include <THashList.h>
92 #include "AliMUONPedestal.h"
93 #include "AliMUONGain.h"
94 #include "AliMUONErrorCounter.h"
97 int main(Int_t argc, Char_t **argv)
104 // needed for streamer application
105 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
111 // needed for Minuit plugin
112 gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
116 "TMinuitMinimizer(const char*)");
118 Char_t prefixDA[256]="MUONTRKGAINda"; // program prefix
120 Int_t skipEvents = 0;
121 Int_t maxEvents = 1000000;
122 Int_t maxDateEvents = 1000000;
123 Char_t inputFile[256]="";
125 Int_t nDateEvents = 0;
126 Int_t nGlitchErrors= 0;
127 Int_t nParityErrors= 0;
128 Int_t nPaddingErrors= 0;
130 TString logOutputFile;
132 Char_t flatFile[256]="";
135 Int_t nEventsRecovered = 0;
137 UInt_t runNumber = 0;
143 Int_t printLevel = 0; // printout (default=0, =1 =>.ped , => .peak & .param)
144 Int_t plotLevel = 1; // plotout (default=1 => tree , =2 tree+Tgraph+fit)
147 Int_t nEntries = daqDA_ECS_getTotalIteration(); // usually = 11 = Nb of calibration runs
148 Int_t nInit=1; // = 0 all DAC values ; = 1 DAC=0 excluded (default=1)
149 Int_t nbpf1=6; // nb of points for linear fit (default=6)
152 // decode the input line
153 for (Int_t i = 1; i < argc; i++) // argument 0 is the executable name
160 // If only one argument and no "-" => DA calling from ECS
163 sprintf(inputFile,argv[i]);
171 sprintf(inputFile,argv[i]);
175 shuttleFile = argv[i];
179 // printLevel=atoi(argv[i]);
183 // plotLevel=atoi(argv[i]);
187 // nbpf1=atoi(argv[i]);
191 // nInit=atoi(argv[i]);
195 skipEvents=atoi(argv[i]);
199 sscanf(argv[i],"%d",&maxDateEvents);
203 sscanf(argv[i],"%d",&maxEvents);
207 printf("\n******************* %s usage **********************",argv[0]);
208 printf("\nOnline (called from ECS) : %s <raw data file> (no inline options)\n",argv[0]);
209 printf("\n%s -options, the available options are :",argv[0]);
210 printf("\n-h help (this screen)");
213 printf("\n-f <raw data file> (default = %s)",inputFile);
216 printf("\n-a <Flat ASCII file> (default = %s)",shuttleFile.Data());
218 printf("\n Options");
219 // printf("\n-d <print level> (default = %d)",printLevel);
220 // printf("\n-g <plot level> (default = %d)",plotLevel);
221 // printf("\n-i <nb linear points> (default = %d)",nbpf1);
222 // printf("\n-j start point of fit (default = %d)",nInit);
223 printf("\n-m <max date events> (default = %d)",maxDateEvents);
224 printf("\n-s <skip events> (default = %d)",skipEvents);
225 printf("\n-n <max events> (default = %d)",maxEvents);
230 printf("%s : bad argument %s (please check %s -h)\n",argv[0],argv[i],argv[0]);
231 argc = 2; exit(-1); // exit if error
237 AliMUONGain* muonGain = new AliMUONGain();
238 muonGain->SetprefixDA(prefixDA);
239 muonGain->SetAliRootDataFileName(); // MUONTRKGAINda_data.root
241 // Reading current iteration
242 nIndex = daqDA_ECS_getCurrentIteration();
243 if(nIndex<0)nIndex=0; // compute gain directly from existing root data file
248 if(nIndex==0) // compute gain from existing root data file: MUONTRKGAINda_data.root
250 sprintf(flatFile,"%s_data.log",prefixDA);
251 logOutputFile=flatFile;
252 filcout.open(logOutputFile.Data());
253 muonGain->SetAlifilcout(&filcout);
255 sprintf(flatFile,"%s_data.par",prefixDA);
256 if(shuttleFile.IsNull())shuttleFile=flatFile;
258 muonGain->SetAliPrintLevel(printLevel);
259 muonGain->SetAliPlotLevel(plotLevel);
260 muonGain->SetAliIndex(nIndex); // fIndex
261 muonGain->SetAliInit(nInit); // fnInit
262 muonGain->SetAliNbpf1(nbpf1); // fnbpf1
263 muonGain->MakeGainStore(shuttleFile);
268 if(nIndex>0) // normal case: reading calibration runs before computing gains
274 // DAC values stored in array vDAC (reading dbfile in DETDB)
275 // Int_t vDAC[11] = {0,200,400,800,1200,1600,2000,2500,3000,3500,4000}; // DAC values
276 Int_t vDAC[11]; // DAC values
278 dbfile="mutrkcalibvalues";
279 cout << "\n *** Copy: " << dbfile << " from DetDB to working directory *** \n" << endl;
280 status=daqDA_DB_getFile(dbfile,dbfile);
281 ifstream filein(dbfile,ios::in); Int_t k=0, kk;
282 while (k<nEntries ) { filein >> kk >> vDAC[k] ; k++; }
283 filein >> nInit; // = 0 all DAC values fitted ; = 1 DAC=0 excluded (default=1)
284 filein >> nbpf1; // nb of points for linear fit (default=6)
285 filein >> printLevel; // printout (default=0, =1 =>.ped /run, =2 => .peak & .param)
286 filein >> plotLevel; // plotout (default=1 => tree , =2 tree+Tgraph+fit)
287 cout << "nInit=" << nInit << " Nb linear pts=" << nbpf1 << " Print level=" << printLevel << " Plot Level=" << plotLevel << endl;
288 muonGain->SetAliPrintLevel(printLevel);
289 muonGain->SetAliPlotLevel(plotLevel);
291 injCharge=vDAC[nIndex-1];
293 // Rawdeader, RawStreamHP
294 AliRawReader* rawReader = AliRawReader::Create(inputFile);
295 AliMUONRawStreamTrackerHP* rawStream = new AliMUONRawStreamTrackerHP(rawReader);
296 rawStream->DisableWarnings();
297 rawStream->EnabbleErrorLogger();
299 cout << "\n" << prefixDA << " : Reading data from file " << inputFile << endl;
301 while (rawReader->NextEvent())
303 if (nDateEvents >= maxDateEvents) break;
304 if (nEvents >= maxEvents) break;
305 if (nDateEvents>0 && nDateEvents % 100 == 0)
306 cout<<"Cumulated: DATE events = " << nDateEvents << " Used events = " << nEvents << endl;
308 // check shutdown condition
309 if (daqDA_checkShutdown())
315 rawReader->NextEvent();
319 Int_t eventType = rawReader->GetType();
320 runNumber = rawReader->GetRunNumber();
322 // Output log file initialisations
325 sprintf(flatFile,"%s.log",prefixDA);
326 logOutputFile=flatFile;
328 filcout.open(logOutputFile.Data());
329 filcout<<"//=================================================" << endl;
330 filcout<<"// " << prefixDA << " for run = " << runNumber << " (DAC=" << injCharge << ")" << endl;
331 filcout<<"//=================================================" << endl;
332 filcout<<"// * Date : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
334 cout<<"\n ******** " << prefixDA << " for run = " << runNumber << " (Index= " << nIndex << "/" << nEntries << " DAC=" << injCharge << ") ********\n" << endl;
335 cout<<" * Date : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
338 muonGain->SetAlifilcout(&filcout);
341 if (eventType != PHYSICS_EVENT)
342 continue; // for the moment
344 // First lopp over DDL's to find good events
345 // Error counters per event (counters in the decoding lib are for each DDL)
346 Bool_t eventIsErrorMessage = kFALSE;
347 int eventGlitchErrors = 0;
348 int eventParityErrors = 0;
349 int eventPaddingErrors = 0;
353 if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
354 eventGlitchErrors += rawStream->GetGlitchErrors();
355 eventParityErrors += rawStream->GetParityErrors();
356 eventPaddingErrors += rawStream->GetPaddingErrors();
357 } while(rawStream->NextDDL());
359 AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
360 if (!eventIsErrorMessage)
362 // Good events (no error) -> compute pedestal for all channels
364 while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next()))
366 for(int i = 0; i < busPatch->GetLength(); ++i)
368 if (nEvents == 0) nChannel++;
369 busPatch->GetData(i, manuId, channelId, charge);
370 muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
377 // Events with errors
378 if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
380 // Recover parity errors -> compute pedestal for all good buspatches
381 if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
384 filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
385 <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
386 <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;
390 filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
391 <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
392 <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
395 while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next()))
397 // Check the buspatch -> if error not use it in the pedestal calculation
399 for(int i = 0; i < busPatch->GetLength(); ++i)
401 if (!busPatch->IsParityOk(i)) errorCount++;
406 for(int i = 0; i < busPatch->GetLength(); ++i)
408 if (nEvents == 0) nChannel++;
409 busPatch->GetData(i, manuId, channelId, charge);
410 // if (busPatch->GetBusPatchId()==1719 && manuId == 1 && channelId == 0) cout <<"Recovered charge "<<charge<<endl;
411 muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
417 AliMUONErrorCounter* errorCounter;
418 // Bad buspatch -> not used (just print)
419 filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
420 <<" parity errors "<<errorCount<<endl;
421 // Number of events where this buspatch is missing
422 sprintf(bpname,"bp%d",busPatch->GetBusPatchId());
423 if (!(errorCounter = (AliMUONErrorCounter*) (muonGain->GetErrorBuspatchTable()->FindObject(bpname))))
426 errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
427 errorCounter->SetName(bpname);
428 muonGain->GetErrorBuspatchTable()->Add(errorCounter);
433 errorCounter->Increment();
435 // errorCounter->Print();
436 } // end of if (!errorCount)
437 } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
440 } //end of if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
443 // Fatal errors reject the event
444 if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
447 filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
448 <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
449 <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;
453 filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
454 <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
455 <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
458 } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
459 filcout<<"Number of errors : Glitch "<<eventGlitchErrors
460 <<" Parity "<<eventParityErrors
461 <<" Padding "<<eventPaddingErrors<<endl;
463 } // end of if (!rawStream->IsErrorMessage())
465 if (eventGlitchErrors) nGlitchErrors++;
466 if (eventParityErrors) nParityErrors++;
467 if (eventPaddingErrors) nPaddingErrors++;
469 } // while (rawReader->NextEvent())
473 // process and store mean peak values in .root file
474 sprintf(flatFile,"%s.par",prefixDA);
475 if(shuttleFile.IsNull())shuttleFile=flatFile;
476 injCharge=vDAC[nIndex-1];
477 muonGain->SetAliIndex(nIndex); // fIndex
478 muonGain->SetAliInjCharge(injCharge);
479 muonGain->SetAliNEvents(nEvents);
480 muonGain->SetAliRunNumber(runNumber);
481 muonGain->SetAliNChannel(nChannel);
482 muonGain->MakePedStoreForGain(shuttleFile);
484 // writing some counters
486 cout << prefixDA << " : Nb of DATE events = " << nDateEvents << endl;
487 cout << prefixDA << " : Nb of Glitch errors = " << nGlitchErrors << endl;
488 cout << prefixDA << " : Nb of Parity errors = " << nParityErrors << endl;
489 cout << prefixDA << " : Nb of Padding errors = " << nPaddingErrors << endl;
490 cout << prefixDA << " : Nb of events recovered = " << nEventsRecovered<< endl;
491 cout << prefixDA << " : Nb of events without errors = " << nEvents-nEventsRecovered<< endl;
492 cout << prefixDA << " : Nb of events used = " << nEvents << endl;
495 filcout << prefixDA << " : Nb of DATE events = " << nDateEvents << endl;
496 filcout << prefixDA << " : Nb of Glitch errors = " << nGlitchErrors << endl;
497 filcout << prefixDA << " : Nb of Parity errors = " << nParityErrors << endl;
498 filcout << prefixDA << " : Nb of Padding errors = " << nPaddingErrors << endl;
499 filcout << prefixDA << " : Nb of events recovered = " << nEventsRecovered<< endl;
500 filcout << prefixDA << " : Nb of events without errors = " << nEvents-nEventsRecovered<< endl;
501 filcout << prefixDA << " : Nb of events used = " << nEvents << endl;
506 muonGain->SetAliInit(nInit); // fnInit
507 muonGain->SetAliEntries(nEntries); // fnEntries
508 muonGain->SetAliNbpf1(nbpf1); // fnbpf1
509 muonGain->MakeGainStore(shuttleFile);
514 filcout << prefixDA << " : Root data file : " << muonGain->GetRootDataFileName() << endl;
515 filcout << prefixDA << " : Output logfile : " << logOutputFile << endl;
516 filcout << prefixDA << " : Gain Histo file : " << muonGain->GetHistoFileName() << endl;
517 filcout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;
519 // Copying files to local DB folder defined by DAQ_DETDB_LOCAL
521 dir= getenv("DAQ_DETDB_LOCAL");
522 unsigned int nLastVersions = 99;
523 cout << "\n *** Local DataBase: " << dir << " (Max= " << nLastVersions << ") ***" << endl;
524 status = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
525 printf(" Store file : %s status = %d\n",logOutputFile.Data(),status);
528 status = daqDA_localDB_storeFile(muonGain->GetRootDataFileName(),nLastVersions);
529 printf(" Store file : %s status = %d\n",muonGain->GetRootDataFileName(),status);
530 status = daqDA_localDB_storeFile(muonGain->GetHistoFileName(),nLastVersions);
531 printf(" Store file : %s status = %d\n",muonGain->GetHistoFileName(),status);
532 status = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions);
533 printf(" Store file : %s status = %d\n",shuttleFile.Data(),status);
540 cout << prefixDA << " : Root data file : " << muonGain->GetRootDataFileName() << endl;
541 cout << prefixDA << " : Output logfile : " << logOutputFile << endl;
542 cout << prefixDA << " : Gain Histo file : " << muonGain->GetHistoFileName() << endl;
543 cout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;
546 // Transferring to OCDB via the SHUTTLE
547 printf("\n ***** STORE FILE in FES ****** \n");
549 // be sure that env variable DAQDALIB_PATH is set in script file
550 // gSystem->Setenv("DAQDALIB_PATH", "$DATE_SITE/infoLogger");
552 status = daqDA_FES_storeFile(shuttleFile.Data(),"GAINS");
555 printf(" Failed to export file : %d\n",status);
557 else printf(" %s successfully exported to FES \n",shuttleFile.Data());
562 printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());