]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/HMPIDphysda.cxx
attempting to add support for par arhcives in forward
[u/mrichter/AliRoot.git] / HMPID / HMPIDphysda.cxx
1 /*
2
3 HMPID PHYSICS DA for noise and DDL monitoring
4
5 Contact:                  Levente.Molnar@cern.ch, Giacomo.Volpe@ba.infn.it
6 Link:                     https://twiki.cern.ch/twiki/bin/view/ALICE/DAInstructions 
7 Run Type:                 PHSYICS and STANDALONE
8 DA Type:                  MON
9 Reference Run:            104157
10 Number of events needed:  ~ 2000 events
11 Input Files:              none
12 Output Files:             HmpPhysicsDaNoiseMap.root
13 Trigger types used:       PHSYICS and STANDALONE
14
15 */
16
17 #define NUM_PHYSICS_EVENTS 2000
18 #define HMP_OUTPUT_FILE "HmpPhysicsDaNoiseMap.root"
19 #include "event.h"
20 #include "monitor.h"
21 //DAQ
22 extern "C" {
23 #include <daqDA.h>
24 }
25 //
26 #include <stdio.h>
27 #include <stdlib.h>
28 //ROOT
29 #include <TPluginManager.h>
30 #include <TROOT.h>
31 #include <TSystem.h>
32 #include <TFile.h>
33 #include <TH2.h>
34 //AliRoot
35 #include "AliRawReaderDate.h"
36 #include "AliHMPIDRawStream.h"
37 #include "AliHMPIDParam.h"
38
39
40 int main(int argc, char **argv) {
41
42   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
43                                         "*",
44                                         "TStreamerInfo",
45                                         "RIO",
46                                         "TStreamerInfo()"); 
47   int status;  
48   if (argc<2) { printf("HMP PhysicsDa: Wrong number of arguments\n"); return -1; }
49
50   /* define data source : this is argument 1 */  
51   status=monitorSetDataSource( argv[1] );
52   if (status!=0) {
53     printf("HMP PhysicsDa: monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
54     return -1;
55   }
56
57
58   /* declare monitoring program */
59   status=monitorDeclareMp( __FILE__ );
60   if (status!=0) {
61     printf("HMP PhysicsDa: monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
62     return -1;
63   }
64
65
66   /* define wait event timeout - 1s max */
67   monitorSetNowait();
68   monitorSetNoWaitNetworkTimeout(1000);
69   
70
71   /* log start of process */
72   printf("HMP PhysicsDa: HMPID PHYSICS Monitoring DA program has started\n");  
73
74
75   /* init some counters */
76   int nevents_physics=0;
77   int nevents_total=0;
78   
79   struct eventHeaderStruct *event;
80   eventTypeType eventT;
81   
82   Int_t runNumber=0;
83   AliHMPIDDigit dig;
84   TH2F *hHmpNoiseMaps=new TH2F("hHmpNoiseMaps","Noise Maps Ch: 0-7;Ch 0-7: pad X;Ch0, Ch1, Ch2, Ch3, Ch4, Ch5, Ch6 pad Y ;Noise level (\%)",160,0,160,144*7,0,144*7); //In y we have all 7 chambers
85  
86   /* main loop (infinite) */
87   for(;;) {
88     /* check shutdown condition */
89     if (daqDA_checkShutdown()) {break;}
90     
91     /* get next event (blocking call until timeout) */
92     status=monitorGetEventDynamic((void **)&event);
93     if (status==MON_ERR_EOF) {
94       printf ("HMP PhysicsDa: End of File detected\n");
95       break; /* end of monitoring file has been reached */
96     }
97     
98     if (status!=0) {
99       printf("HMP PhysicsDa: monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
100       break;
101     }
102     
103     /* retry if got no event */
104     if (event==NULL) {
105       continue;
106     }
107
108     /* use event - here, just write event id to result file */
109     eventT=event->eventType;
110     
111     if (eventT==PHYSICS_EVENT) {
112       runNumber=(unsigned long)event->eventRunNb;   
113       AliRawReader *reader = new AliRawReaderDate((void*)event);
114       AliHMPIDRawStream stream(reader);stream.SetTurbo(kTRUE);                  //raw data decoding without error checks SetTurbo(kTRUE)
115       while(stream.Next())
116          {
117             for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
118               dig.SetPad(stream.GetPadArray()[iPad]);
119               hHmpNoiseMaps->Fill( dig.PadChX(), dig.Ch()*144+dig.PadChY() );
120             } //pads            
121           }//while -- loop on det load in one event
122           
123       nevents_physics++;
124     }//physics event
125     nevents_total++;
126   
127     /* free resources */
128     free(event);
129     
130     /* exit when last event received, no need to wait for TERM signal */
131     if (eventT==END_OF_RUN) {
132       printf("HMP PhysicsDa: EOR event detected\n");
133       break;
134     }
135     /* Exit after the NUM_PHYSICS_EVENTS physics events */
136     if( nevents_physics == NUM_PHYSICS_EVENTS ) 
137     {
138      printf("HMP PhysicsDa: The number of required physics events (%d) is reached!\n",NUM_PHYSICS_EVENTS);
139      break;  
140     }    
141   }//main loop
142
143
144   /* write report */
145   printf("HMP PhysicsDa: Run #%s, received %d physics events out of %d\n",getenv("DATE_RUN_NUMBER"),nevents_physics,nevents_total);
146
147   /* set histogram properties */
148
149   hHmpNoiseMaps->SetTitle(Form("Run number: %d Tested events: %d",runNumber,nevents_physics)); 
150   if(nevents_physics!=0) hHmpNoiseMaps->Scale(1.0/nevents_physics);
151
152
153   /* write output file */
154   TFile *fout = new TFile(HMP_OUTPUT_FILE,"recreate");
155   hHmpNoiseMaps->Write();
156   fout->Close();
157   
158   /* send file to DAQ FXS */
159   status=0;
160   status=daqDA_FES_storeFile(HMP_OUTPUT_FILE,HMP_OUTPUT_FILE);
161   if(status) return -1;
162     
163   return status;
164 }