2 PMD DA for online calibration
4 contact: basanta@phy.iitb.ac.in
6 Reference run:/afs/cern.ch/user/s/sjena/public/run83496.raw
9 Number of events needed: 1 million for PB+PB, 200 milion for p+p
10 Input Files: Run0_999999999_v0_s0.root,PMD_PED.root, PMD_GAIN_CONFIGFILE, pmd_gain_tempfile.dat
11 Output Files: PMDGAINS.root, to be exported to the DAQ FES
12 Trigger types used: PHYSICS_EVENT
22 #include <Riostream.h>
27 #include "AliRawReaderDate.h"
28 #include "AliPMDCalibPedestal.h"
29 #include "AliPMDCalibGain.h"
31 #include "AliCDBManager.h"
36 #include "TBenchmark.h"
39 #include "TPluginManager.h"
45 1- monitoring data source
47 int main(int argc, char **argv) {
49 /* magic line from Rene */
50 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
60 Int_t filestatus = -1, xvar = 5;
61 Int_t totevt = -1, maxevt = -1;
62 Int_t hotevtsize = -1;
63 Bool_t hotfilestatus = false;
65 // Reads the pedestal file and keep the values in memory for subtraction
67 AliPMDCalibGain *calibgain = new AliPMDCalibGain();
69 // Fetch the pedestal file - PMD_PED.root
71 status = daqDA_DB_getFile("PMD_PED.root","PMD_PED.root");
75 printf("*** Pedestal file retrieved from DB *** \n");
79 printf("*** Pedestal file NOT retrieved from DB *** \n");
83 Int_t pstatus = calibgain->ExtractPedestal("PMD_PED.root");
85 if(pstatus == -3) return -3;
90 // Retrieve the PMD_GAIN_CONFIGFILE
91 status = daqDA_DB_getFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
95 fp1 = fopen("PMD_GAIN_CONFIGFILE","r");
99 printf("*** PMD GAIN Configfile doesn't exist,Provide one ***\n");
104 fscanf(fp1,"%d %d %d %d %d\n",&filestatus, &xvar, &totevt, &maxevt, &hotevtsize);
105 //printf("%d %d %d %d %d\n",filestatus, xvar, totevt, maxevt, hotevtsize);
112 // Retrieve the Temporray ascii file from DB
113 status = daqDA_DB_getFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
116 calibgain->ReadTempFile("pmd_gain_tempfile.dat");
120 printf("--- pmd_gain_tempfile.dat: not retrieved from DB ---\n");
122 // Retrieve the hot cell file from DB - PMD_HOT.root
123 status = daqDA_DB_getFile("PMD_HOT.root","PMD_HOT.root");
126 calibgain->ExtractHotChannel("PMD_HOT.root");
130 printf("--- pmd_gain_tempfile.dat: not retrieved from DB ---\n");
135 // decoding the events
139 printf("Wrong number of arguments\n");
144 /* define data source : this is argument 1 */
145 status=monitorSetDataSource( argv[1] );
147 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
151 /* declare monitoring program */
152 status=monitorDeclareMp( __FILE__ );
154 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
158 /* define wait event timeout - 1s max */
160 monitorSetNoWaitNetworkTimeout(1000);
162 /* log start of process */
163 printf("PMD GAIN DA - strted generating the gain of a cell\n");
165 /* init some counters */
166 int nevents_physics=0;
169 struct eventHeaderStruct *event;
170 eventTypeType eventT = 0;
176 if (getenv("DATE_RUN_NUMBER")==0) {
177 printf("DATE_RUN_NUMBER not properly set.\n");
180 int runNr = atoi(getenv("DATE_RUN_NUMBER"));
184 if (gSystem->AccessPathName("localOCDB/PMD/Calib/Mapping",kFileExists))
186 if (gSystem->mkdir("localOCDB/PMD/Calib/Mapping",kTRUE) != 0)
188 printf("Failed to create directory: localOCDB/PMD/Calib/Mapping");
192 status = daqDA_DB_getFile("PMD/Calib/Mapping","localOCDB/PMD/Calib/Mapping/Run0_999999999_v0_s0.root");
195 printf("Failed to get PMD-Mapping file (PMD/Calib/Mapping) from DAQdetDB, status=%d\n", status);
199 // Global initializations
200 AliLog::SetGlobalLogLevel(AliLog::kError);
201 AliCDBManager *man = AliCDBManager::Instance();
202 man->SetDefaultStorage("local://localOCDB");
206 /* main loop (infinite) */
209 /* check shutdown condition */
210 if (daqDA_checkShutdown()) {break;}
212 /* get next event (blocking call until timeout) */
213 status=monitorGetEventDynamic((void **)&event);
214 if (status==MON_ERR_EOF) {
215 printf ("End of File detected\n");
216 break; /* end of monitoring file has been reached */
220 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
224 /* retry if got no event */
231 /* use event - here, just write event id to result file */
233 eventT=event->eventType;
234 switch (event->eventType){
239 /* END START OF RUN */
248 //if(nevents_physics%100 == 0)printf("Physis Events = %d\n",nevents_physics);
249 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
250 TObjArray *pmdddlcont = new TObjArray();
251 calibgain->ProcessEvent(rawReader, pmdddlcont);
253 if (totevt%hotevtsize == 0) hotfilestatus = true;
266 /* exit when last event received, no need to wait for TERM signal */
269 ic = new TTree("ic","PMD Gain tree");
270 meanc = new TTree("meanc","PMD Module mean tree");
274 TFile *hotRun = new TFile ("PMD_HOT.root","RECREATE");
276 TTree *hot = new TTree("hot","PMD Hot cell tree");
278 calibgain->FindHotCell(hot,xvar);
283 // store the hot cell root file in the DB
285 status = daqDA_DB_storeFile("PMD_HOT.root","PMD_HOT.root");
287 // store the hot cell root file in the file exchange server
289 printf("root file for hot cell is created and getting exported\n");
290 status = daqDA_FES_storeFile("PMD_HOT.root","PMD_HOT.root");
295 TFile *hotRun = new TFile ("PMD_HOT.root","RECREATE");
297 TTree *hot = new TTree("hot","PMD Hot cell tree");
299 calibgain->FindHotCell(hot,xvar);
304 // store the hot cell root file in the DB
306 status = daqDA_DB_storeFile("PMD_HOT.root","PMD_HOT.root");
308 // store the hot cell root file in the file exchange server
310 printf("root file for hot cell is created and getting exported\n");
311 status = daqDA_FES_storeFile("PMD_HOT.root","PMD_HOT.root");
315 fp1 = fopen("PMD_GAIN_CONFIGFILE","w+");
319 printf("-----------------------------------------------\n");
320 printf("*** Required Number of Events not reached ***\n");
321 printf("*** Number of Events processed = %d ***\n",totevt);
322 printf("*** Writing the intermediate ASCII file ***\n");
323 printf("-----------------------------------------------\n");
325 calibgain->WriteTempFile("pmd_gain_tempfile.dat");
327 // Store the Intermediate ascii file in the DB
328 status = daqDA_DB_storeFile("pmd_gain_tempfile.dat","pmd_gain_tempfile.dat");
331 fprintf(fp1,"%d %d %d %d %d\n",filestatus,xvar,totevt,maxevt,hotevtsize);
334 // Store the configfile in the DB
335 status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
338 else if (totevt >= maxevt)
340 printf("-----------------------------------------------\n");
341 printf("*** Required Number of Events reached = %d ***\n",totevt);
342 printf("*** Writing the PMDGAINS.root file ***\n");
343 printf("-----------------------------------------------\n");
345 calibgain->Analyse(ic, meanc);
347 TFile * gainRun = new TFile ("PMDGAINS.root","RECREATE");
351 TFile * meanRun = new TFile ("PMD_MEAN_SM.root","RECREATE");
358 fprintf(fp1,"%d %d %d %d %d\n",filestatus,xvar,totevt,maxevt,hotevtsize);
361 // Store the configfile in the DB
362 status = daqDA_DB_storeFile("PMD_GAIN_CONFIGFILE","PMD_GAIN_CONFIGFILE");
373 /* store the result file on FES */
377 printf("root file for cell gain is created and getting exported\n");
378 status = daqDA_FES_storeFile("PMDGAINS.root","PMDGAINS.root");
379 printf("root file for normalised means of different modules\n");
380 status = daqDA_FES_storeFile("PMD_MEAN_SM.root","PMD_MEAN_SM.root");