]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTRKGAINda.cxx
In AliMUONTriggerQADataMakerRec:
[u/mrichter/AliRoot.git] / MUON / MUONTRKGAINda.cxx
1 /*
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: 92796-92798-92800-92801-92803-92804-92807-92807-92809-92810-92811
5  Run Type: CALIBRATION
6  DA Type: LDC
7  Number of events needed: 400 events for each calibration run (11)
8  Input Files: mutrkcalibvalues and config_ldc-MTRK-S3-0
9  Output Files: local dir (not persistent) -> MUONTRKGAINda.par   FXS -> run<#>_MCH_<ldc>_GAINS
10  Trigger types used:
11  */
12
13 /**************************************************************************
14  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
15  *                                                                        *
16  * Author: The ALICE Off-line Project.                                    *
17  * Contributors are mentioned in the code where appropriate.              *
18  *                                                                        *
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  **************************************************************************/
27
28 /* $Id$ */
29
30 /*
31  -------------------------------------------------------------------------
32  2009-11-05 New version: MUONTRKGAINda.cxx,v 1.4
33  -------------------------------------------------------------------------
34  
35  Version for MUONTRKGAINda MUON tracking
36  (A. Baldisseri, J.-L. Charvet)
37  
38  
39  Rem:  AliMUON2DMap stores all channels, even those which are not connected
40  if pedMean == -1, channel not connected to a pad  
41  
42  
43  */
44 extern "C" {
45 #include <daqDA.h>
46 }
47
48 #include "event.h"
49 #include "monitor.h"
50
51 #include <Riostream.h>
52 #include <stdio.h>
53 #include <stdlib.h>
54 #include <sstream>
55 #include <math.h> 
56
57 //AliRoot
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"
67
68 //ROOT
69 #include "TFile.h"
70 #include "TSystem.h"
71 #include "TTree.h"
72 #include "TH1F.h"
73 #include "TString.h"
74 #include "TStopwatch.h"
75 #include "TMath.h"
76 #include "TTimeStamp.h"
77 #include "TGraphErrors.h"
78 #include "TF1.h"
79 #include "TROOT.h"
80 #include "TPluginManager.h"
81 #include "TFitter.h"
82 #include "TObjString.h"
83 #include "THashTable.h"
84 #include <THashList.h>
85 //
86 //AMORE
87 //
88 #ifdef ALI_AMORE
89 #include <AmoreDA.h>
90 #endif
91
92 #include "AliMUONPedestal.h"
93 #include "AliMUONGain.h"
94 #include "AliMUONErrorCounter.h"
95
96 // main routine
97 int main(Int_t argc, const char** argv) 
98 {
99   Int_t status=0;
100   TStopwatch timers;
101   timers.Start(kTRUE); 
102   
103   const char* prefixDA = "MUONTRKGAINda"; // program prefix
104   printf(" ######## Begin execution : %s ######## \n\n",prefixDA); 
105   
106   TString inputFile;
107   // decode the input line
108   if (argc!=2) { printf("Wrong number of arguments\n");  return -1; }
109   
110   inputFile = argv[1];
111   
112   // needed for streamer application
113   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
114                                         "*",
115                                         "TStreamerInfo",
116                                         "RIO",
117                                         "TStreamerInfo()"); 
118   // needed for Minuit plugin
119   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer",
120                                         "Minuit",
121                                         "TMinuitMinimizer",
122                                         "Minuit",
123                                         "TMinuitMinimizer(const char*)");
124   Int_t skipEvents = 0;
125   Int_t maxEvents  = 1000000;
126   Int_t maxDateEvents  = 1000000;
127   
128   Int_t nDateEvents = 0;
129   Int_t nGlitchErrors= 0;
130   Int_t nParityErrors= 0;
131   Int_t nPaddingErrors= 0;
132   Int_t nEventsRecovered = 0;
133   Int_t nEvents = 0;
134   
135   TString logOutputFile;
136   Char_t flatFile[256]="";
137   TString shuttleFile;
138   
139   UInt_t runNumber   = 0;
140   ofstream filcout;
141   Int_t nIndex = -1; 
142   
143   // For DA Gain
144   Int_t nEntries = daqDA_ECS_getTotalIteration(); // usually = 11 = Nb of calibration runs
145   Int_t nInit=1;  // = 0 all DAC values ; = 1 DAC=0 excluded (default=1)
146   Int_t nbpf1=6;  // nb of points for linear fit (default=6) 
147   Int_t printLevel  = 0;  // printout (default=0, =1 =>.ped , => .peak & .param)
148   Int_t plotLevel  = 1;  // plotout (default=1 => tree , =2 tree+Tgraph+fit)
149   Int_t nbev=0; 
150   Int_t injCharge; 
151   
152   // Reading current iteration
153   nIndex = daqDA_ECS_getCurrentIteration();
154   if(nIndex<0 || nIndex>nEntries) {printf("\n Failed: nIndex = %d\n",nIndex); return -1 ;}
155   
156   //Gain object
157   AliMUONGain* muonGain = new AliMUONGain();
158   muonGain->SetprefixDA(prefixDA);
159   muonGain->SetAliRootDataFileName(); // MUONTRKGAINda_data.root
160   
161   UShort_t manuId;  
162   UChar_t channelId;
163   UShort_t charge;
164   
165   // DAC values stored in array vDAC (reading dbfile in DETDB)
166   //   Int_t vDAC[11] = {0,200,400,800,1200,1600,2000,2500,3000,3500,4000}
167   Int_t nConfig = 1; // flag to read or not configuration ascii file in detDB
168   Int_t vDAC[11]; // DAC values
169   Char_t dbfile[256]="";
170   sprintf(dbfile,"mutrkcalibvalues");
171   status=daqDA_DB_getFile(dbfile,dbfile);
172   if(status) {printf(" Failed  : input file %s is missing, status = %d\n",dbfile,status); return -1; } 
173   
174   ifstream filein(dbfile,ios::in); Int_t k=0, kk;
175   while (k<nEntries ) { filein >> kk >> vDAC[k] ; k++; }
176   injCharge=vDAC[nIndex-1];
177   
178   filein >> nInit; // = 0 all DAC values fitted ; = 1 DAC=0 excluded (default=1)
179   filein >> nbpf1; // nb of points for linear fit (default=6) 
180   filein >> printLevel;  // printout (default=0, =1 =>.ped /run, =2 => .peak & .param)
181   filein >> plotLevel;  // plotout (default=1 => tree , =2 tree+Tgraph+fit)
182   filein >> nConfig;  //nConfig (default=1 => read config in DetDB, otherwise =0)
183   filein >> nbev;  // Nb of events to read  (default = 0 => reading all events)
184   if(nbev>0)maxEvents=nbev;
185   
186   printf(" *** Copy: %s from DetDB to working directory  ***      Config= %d\n",dbfile,nConfig);
187   printf(" Input parameters:  nInit= %d   Nb linear pts= %d   Print level= %d   Plot Level= %d",nInit,nbpf1,printLevel,plotLevel);
188   if(nbev==0)printf("\n");
189   else printf("  Nb_max evt = %d\n",maxEvents);
190   
191   muonGain->SetAliPrintLevel(printLevel);
192   muonGain->SetAliPlotLevel(plotLevel);
193   muonGain->SetconfigDA(nConfig);
194   
195   if(nConfig)
196   {
197     // Laod configuration ascii file from DetDB
198     sprintf(dbfile,"config_%s",getenv("DATE_ROLE_NAME"));
199     status=daqDA_DB_getFile(dbfile,dbfile);
200     if(status) {printf(" Failed  : Configuration file %s is missing, status = %d\n",dbfile,status); return -1; }
201     else printf(" *** Copy ascii config file: %s from DetDB to working directory and reading ...*** \n",dbfile);
202     muonGain->LoadConfig(dbfile);  
203   } 
204   
205   // Rawdeader, RawStreamHP
206   AliRawReader* rawReader = AliRawReader::Create(inputFile.Data());
207   AliMUONRawStreamTrackerHP* rawStream  = new AliMUONRawStreamTrackerHP(rawReader);    
208   rawStream->DisableWarnings();
209   rawStream->EnabbleErrorLogger();
210   
211   cout << "\n" << prefixDA << " : Reading data from file " << inputFile.Data()  << endl;
212   
213   while (rawReader->NextEvent())
214   {
215     if (nDateEvents >= maxDateEvents) break;
216     if (nEvents >= maxEvents) break;
217     if (nDateEvents>0 &&  nDateEvents % 100 == 0)       
218       cout<<"Cumulated:  DATE events = " << nDateEvents << "   Used events = " << nEvents << endl;
219     
220     // check shutdown condition 
221     if (daqDA_checkShutdown()) 
222       break;
223     //Skip events
224     while (skipEvents)
225     {
226       rawReader->NextEvent();
227       skipEvents--;
228     }
229     
230     Int_t eventType = rawReader->GetType();
231     runNumber = rawReader->GetRunNumber();
232     
233     // Output log file initialisations
234     if(nDateEvents==0)
235     {
236       sprintf(flatFile,"%s.log",prefixDA);
237       logOutputFile=flatFile;
238       
239       filcout.open(logOutputFile.Data());
240       filcout<<"//=================================================" << endl;
241       filcout<<"//       " << prefixDA << " for run = " << runNumber << "  (DAC=" << injCharge << ")" << endl;
242       filcout<<"//=================================================" << endl;
243       filcout<<"//   * Date          : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
244       
245       cout<<"\n ********  " << prefixDA << " for run = " << runNumber << "  (Index= " << nIndex << "/" << nEntries << "  DAC=" << injCharge << ") ********\n" << endl;
246       cout<<" * Date : " << muonGain->GetDate()->AsString("l") << "\n" << endl;
247     }
248     
249     muonGain->SetAlifilcout(&filcout);
250     
251     nDateEvents++;
252     if (eventType != PHYSICS_EVENT)
253       continue; // for the moment
254     
255     // First lopp over DDL's to find good events
256     // Error counters per event (counters in the decoding lib are for each DDL)
257     Bool_t eventIsErrorMessage = kFALSE;
258     int eventGlitchErrors = 0;
259     int eventParityErrors = 0;
260     int eventPaddingErrors = 0;
261     rawStream->First();
262     do
263     {
264       if (rawStream->IsErrorMessage()) eventIsErrorMessage = kTRUE;
265       eventGlitchErrors += rawStream->GetGlitchErrors();
266       eventParityErrors += rawStream->GetParityErrors();
267       eventPaddingErrors += rawStream->GetPaddingErrors();
268     } while(rawStream->NextDDL()); 
269     
270     AliMUONRawStreamTrackerHP::AliBusPatch* busPatch;
271     if (!eventIsErrorMessage) 
272     {
273       // Good events (no error) -> compute pedestal for all channels
274       rawStream->First(); 
275       while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
276             {
277               for(int i = 0; i < busPatch->GetLength(); ++i)
278         {
279           busPatch->GetData(i, manuId, channelId, charge);
280           muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
281         }
282             }
283       nEvents++;
284     }
285     else
286     {
287       // Events with errors
288       if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
289             {
290               // Recover parity errors -> compute pedestal for all good buspatches
291               if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
292                                    ATTR_ORBIT_BC )) 
293         {
294           filcout <<"Event recovered -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
295           <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
296           <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                            
297         } 
298               else 
299         {
300           filcout <<"Event recovered -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
301           <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
302           <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
303         }
304               rawStream->First();
305               while( (busPatch = (AliMUONRawStreamTrackerHP::AliBusPatch*) rawStream->Next())) 
306         {
307           // Check the buspatch -> if error not use it in the pedestal calculation
308           int errorCount = 0;
309           for(int i = 0; i < busPatch->GetLength(); ++i)
310           {
311             if (!busPatch->IsParityOk(i)) errorCount++;
312           }
313           if (!errorCount) 
314           {
315             // Good buspatch
316             for(int i = 0; i < busPatch->GetLength(); ++i)
317             {
318               busPatch->GetData(i, manuId, channelId, charge);
319               muonGain->MakePed(busPatch->GetBusPatchId(), (Int_t)manuId, (Int_t)channelId, (Int_t)charge);
320             }
321           }
322           else
323           {
324             char bpname[256];
325             AliMUONErrorCounter* errorCounter;
326             // Bad buspatch -> not used (just print)
327             filcout<<"bpId "<<busPatch->GetBusPatchId()<<" words "<<busPatch->GetLength()
328             <<" parity errors "<<errorCount<<endl;
329             // Number of events where this buspatch is missing
330             sprintf(bpname,"bp%d",busPatch->GetBusPatchId());                                           
331             if (!(errorCounter = (AliMUONErrorCounter*) (muonGain->GetErrorBuspatchTable()->FindObject(bpname))))
332             {
333               // New buspatch
334               errorCounter = new AliMUONErrorCounter(busPatch->GetBusPatchId());
335               errorCounter->SetName(bpname);
336               muonGain->GetErrorBuspatchTable()->Add(errorCounter);
337             }
338             else
339             {
340               // Existing buspatch
341               errorCounter->Increment();
342             }   
343             // errorCounter->Print();                                           
344           } // end of if (!errorCount)
345         } // end of while( (busPatch = (AliMUONRawStreamTrackerHP ...
346               nEvents++;
347               nEventsRecovered++;
348             } //end of if (eventParityErrors && !eventGlitchErrors&& !eventPaddingErrors)
349       else
350             {
351               // Fatal errors reject the event
352               if ( TEST_SYSTEM_ATTRIBUTE( rawReader->GetAttributes(),
353                                    ATTR_ORBIT_BC )) 
354         {
355           filcout <<"Event rejected -> Period:"<<EVENT_ID_GET_PERIOD( rawReader->GetEventId() )
356           <<" Orbit:"<<EVENT_ID_GET_ORBIT( rawReader->GetEventId() )
357           <<" BunchCrossing:"<<EVENT_ID_GET_BUNCH_CROSSING( rawReader->GetEventId() )<<endl;                            
358         } 
359               else 
360         {
361           filcout <<"Event rejected -> nbInRun:"<<EVENT_ID_GET_NB_IN_RUN( rawReader->GetEventId() )
362           <<" burstNb:"<<EVENT_ID_GET_BURST_NB( rawReader->GetEventId() )
363           <<" nbInBurst:"<<EVENT_ID_GET_NB_IN_BURST( rawReader->GetEventId() )<<endl;
364           
365         }
366             } // end of if (!rawStream->GetGlitchErrors() && !rawStream->GetPaddingErrors() ...
367       filcout<<"Number of errors : Glitch "<<eventGlitchErrors
368       <<" Parity "<<eventParityErrors
369       <<" Padding "<<eventPaddingErrors<<endl;
370       filcout<<endl;                    
371     } // end of if (!rawStream->IsErrorMessage())
372     
373     if (eventGlitchErrors)  nGlitchErrors++;
374     if (eventParityErrors)  nParityErrors++;
375     if (eventPaddingErrors) nPaddingErrors++;
376     
377   } // while (rawReader->NextEvent())
378   delete rawReader;
379   delete rawStream;
380   
381   // process and store mean peak values in .root file
382   sprintf(flatFile,"%s.par",prefixDA);
383   if(shuttleFile.IsNull())shuttleFile=flatFile;
384   injCharge=vDAC[nIndex-1];
385   muonGain->SetAliIndex(nIndex); // fIndex 
386   muonGain->SetAliInjCharge(injCharge);
387   muonGain->SetAliNEvents(nEvents);
388   muonGain->SetAliRunNumber(runNumber);
389   muonGain->MakePedStoreForGain(shuttleFile);
390   
391   
392   // writing some counters
393   cout << endl;
394   cout << prefixDA << " : Nb of DATE events           = " << nDateEvents    << endl;
395   cout << prefixDA << " : Nb of Glitch errors         = "   << nGlitchErrors  << endl;
396   cout << prefixDA << " : Nb of Parity errors         = "   << nParityErrors  << endl;
397   cout << prefixDA << " : Nb of Padding errors        = "   << nPaddingErrors << endl;          
398   cout << prefixDA << " : Nb of events recovered      = "   << nEventsRecovered<< endl;
399   cout << prefixDA << " : Nb of events without errors = "   << nEvents-nEventsRecovered<< endl;
400   cout << prefixDA << " : Nb of events used           = "   << nEvents        << endl;
401   
402   filcout << endl;
403   filcout << prefixDA << " : Nb of DATE events           = " << nDateEvents    << endl;
404   filcout << prefixDA << " : Nb of Glitch errors         = "   << nGlitchErrors << endl;
405   filcout << prefixDA << " : Nb of Parity errors         = "   << nParityErrors << endl;
406   filcout << prefixDA << " : Nb of Padding errors        = "   << nPaddingErrors << endl;
407   filcout << prefixDA << " : Nb of events recovered      = "   << nEventsRecovered<< endl;      
408   filcout << prefixDA << " : Nb of events without errors = "   << nEvents-nEventsRecovered<< endl;
409   filcout << prefixDA << " : Nb of events used           = "   << nEvents        << endl;
410   
411   // Computing gain 
412   if(nIndex==nEntries)
413   {
414     muonGain->SetAliInit(nInit); // fnInit
415     muonGain->SetAliEntries(nEntries); // fnEntries
416     muonGain->SetAliNbpf1(nbpf1); // fnbpf1
417     muonGain->MakeGainStore(shuttleFile);
418 #ifdef ALI_AMORE  
419     std::ifstream in(shuttleFile.Data());
420     ostringstream stringout;
421     char line[1024];
422     while ( in.getline(line,1024) )
423       stringout << line << "\n";  
424     in.close();
425           
426     amore::da::AmoreDA amoreDA(amore::da::AmoreDA::kSender);
427     TObjString gaindata(stringout.str().c_str());
428     status = amoreDA.Send("Gains",&gaindata);
429     if ( status )
430       cout << "Warning: Failed to write Pedestals in the AMORE database : " << status << endl;
431     else 
432       cout << "amoreDA.Send(Gains) ok" << endl;  
433 #else
434     cout << "Warning: MCH DA not compiled with AMORE support" << endl;
435 #endif
436   }
437   
438   // ouput files
439   filcout << endl;
440   filcout << prefixDA << " : Root data file         : " << muonGain->GetRootDataFileName() << endl;
441   filcout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
442   filcout << prefixDA << " : Gain Histo file        : " << muonGain->GetHistoFileName() << endl;
443   filcout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;
444   
445   //     Copying files to local DB folder defined by DAQ_DETDB_LOCAL
446   Char_t *dir;
447   dir= getenv("DAQ_DETDB_LOCAL");
448   unsigned int nLastVersions = 90;
449   printf("\n ***  Local DataBase: %s  (Max= %d) ***\n",dir,nLastVersions);
450   status = daqDA_localDB_storeFile(logOutputFile.Data(),nLastVersions);
451   if(status)printf(" Store file : %s   status = %d\n",logOutputFile.Data(),status);
452   if(nIndex==nEntries)
453   {
454     status = daqDA_localDB_storeFile(muonGain->GetRootDataFileName(),nLastVersions);
455     if(status)printf(" Store file : %s   status = %d\n",muonGain->GetRootDataFileName(),status);
456     status = daqDA_localDB_storeFile(muonGain->GetHistoFileName(),nLastVersions);
457     if(status)printf(" Store file : %s   status = %d\n",muonGain->GetHistoFileName(),status);
458     status = daqDA_localDB_storeFile(shuttleFile.Data(),nLastVersions);
459     if(status)printf(" Store file : %s   status = %d\n",shuttleFile.Data(),status);
460   }      
461   
462   
463   // ouput files
464   cout << endl;
465   cout << prefixDA << " : Root data file         : " << muonGain->GetRootDataFileName() << endl;
466   cout << prefixDA << " : Output logfile         : " << logOutputFile  << endl;
467   cout << prefixDA << " : Gain Histo file        : " << muonGain->GetHistoFileName() << endl;
468   cout << prefixDA << " : Gain file (to SHUTTLE) : " << shuttleFile << endl;   
469   
470   filcout.close();
471   
472   // Transferring to OCDB via the SHUTTLE
473   // be sure that env variable DAQDALIB_PATH is set in script file
474   //       gSystem->Setenv("DAQDALIB_PATH", "$DATE_SITE/infoLogger");
475   printf("\n *****  STORE FILE in FES ****** \n");
476   status = daqDA_FES_storeFile(shuttleFile.Data(),"GAINS");
477   if (status) { printf(" Failed to export file : %s , status = %d\n",shuttleFile.Data(),status); return -1; }
478   else printf(" %s successfully exported to FES  \n",shuttleFile.Data());
479   
480   printf("\n ######## End execution : %s ######## \n",prefixDA); 
481   timers.Stop();
482   printf("\nExecution time : R:%7.2fs C:%7.2fs\n", timers.RealTime(), timers.CpuTime());
483   return status;
484 }
485