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