]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSGAINda.cxx
fatal error in case OCDB entry with alignment objects is not found for one or more...
[u/mrichter/AliRoot.git] / PHOS / PHOSGAINda.cxx
1 /*
2 contact: Boris.Polishchuk@cern.ch
3 link: see comments in the $ALICE_ROOT/PHOS/AliPHOSRcuDA1.cxx
4 reference run: /alice/data/2009/LHC09b_PHOS/000075883/raw/09000075883017.20.root
5 run type: PHYSICS
6 DA type: MON 
7 number of events needed: 1000
8 input files: RCU0.data  RCU1.data  RCU2.data  RCU3.data  zs.txt
9 Output files: PHOS_ModuleN_Calib.root, where N is the module number (0-5).
10 Trigger types used: PHYSICS
11 */
12
13
14 #include "event.h"
15 #include "monitor.h"
16 extern "C" {
17 #include "daqDA.h"
18 }
19
20 #include <stdio.h>
21 #include <stdlib.h>
22
23 #include <TSystem.h>
24 #include <TROOT.h>
25 #include <TPluginManager.h>
26
27 #include "AliRawReader.h"
28 #include "AliRawReaderDate.h"
29 #include "AliPHOSRcuDA1.h"
30 #include "AliPHOSRawFitterv0.h"
31 #include "AliCaloAltroMapping.h"
32 #include "AliCaloRawStreamV3.h"
33 #include "AliLog.h"
34
35 /* Main routine
36       Arguments: 
37       1- monitoring data source
38 */
39 int main(int argc, char **argv) {
40
41   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
42                                         "*",
43                                         "TStreamerInfo",
44                                         "RIO",
45                                         "TStreamerInfo()");
46
47   AliLog::SetGlobalDebugLevel(0) ;
48   AliLog::SetGlobalLogLevel(AliLog::kFatal);
49   
50   int status;
51   
52   if (argc!=2) {
53     printf("Wrong number of arguments\n");
54     return -1;
55   }
56   
57   /* Retrieve ZS parameters from DAQ DB */
58   const char* zsfile = "zs.txt";
59   int failZS = daqDA_DB_getFile(zsfile, zsfile);
60   
61   Int_t offset,threshold;
62   
63   if(!failZS) {
64     FILE *f = fopen(zsfile,"r");
65     int scan = fscanf(f,"%d %d",&offset,&threshold);
66   }
67   
68   /* Retrieve mapping files from DAQ DB */ 
69   const char* mapFiles[4] = {"RCU0.data","RCU1.data","RCU2.data","RCU3.data"};
70   
71   for(Int_t iFile=0; iFile<4; iFile++) {
72     int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
73     if(failed) { 
74       printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);
75       return -1;
76     }
77   }
78   
79   /* Open mapping files */
80   AliAltroMapping *mapping[4];
81   TString path = "./";
82   path += "RCU";
83   TString path2;
84   for(Int_t i = 0; i < 4; i++) {
85     path2 = path;
86     path2 += i;
87     path2 += ".data";
88     mapping[i] = new AliCaloAltroMapping(path2.Data());
89   }
90   
91   /* define data source : this is argument 1 */  
92   status=monitorSetDataSource( argv[1] );
93   if (status!=0) {
94     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
95     return -1;
96   }
97
98
99   /* declare monitoring program */
100   status=monitorDeclareMp( __FILE__ );
101   if (status!=0) {
102     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
103     return -1;
104   }
105
106
107   /* define wait event timeout - 1s max */
108   monitorSetNowait();
109   monitorSetNoWaitNetworkTimeout(1000);
110   
111    /* init some counters */
112   int nevents_physics=0;
113   int nevents_total=0;
114
115   AliRawReader *rawReader = NULL;
116   AliPHOSRcuDA1* dAs[5];
117
118   for(Int_t iMod=0; iMod<5; iMod++) {
119     dAs[iMod] = 0;
120   }
121   
122   Float_t e[64][56][2];
123   Float_t t[64][56][2];
124
125   for(Int_t iX=0; iX<64; iX++) {
126     for(Int_t iZ=0; iZ<56; iZ++) {
127       for(Int_t iGain=0; iGain<2; iGain++) {
128         e[iX][iZ][iGain] = 0.;
129         t[iX][iZ][iGain] = 0.;
130       }
131     }
132   }
133   
134   Int_t cellX    = -1;
135   Int_t cellZ    = -1;
136   Int_t nBunches =  0;
137   Int_t sigStart, sigLength;
138   Int_t caloFlag;
139
140   /* main loop (infinite) */
141   for(;;) {
142     struct eventHeaderStruct *event;
143     eventTypeType eventT;
144   
145     /* check shutdown condition */
146     if (daqDA_checkShutdown()) {break;}
147     
148     /* get next event (blocking call until timeout) */
149     status=monitorGetEventDynamic((void **)&event);
150     if (status==MON_ERR_EOF) {
151       printf ("End of File detected\n");
152       break; /* end of monitoring file has been reached */
153     }
154     
155     if (status!=0) {
156       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
157       break;
158     }
159
160     /* retry if got no event */
161     if (event==NULL) {
162       continue;
163     }
164
165
166     /* use event - here, just write event id to result file */
167     eventT=event->eventType;
168     
169     if (eventT==PHYSICS_EVENT) {
170       
171       rawReader = new AliRawReaderDate((void*)event);
172       AliCaloRawStreamV3 stream(rawReader,"PHOS",mapping);
173       AliPHOSRawFitterv0 fitter;
174       fitter.SubtractPedestals(kTRUE); // assume that data is non-ZS
175       
176       if(!failZS) {
177         fitter.SubtractPedestals(kFALSE);
178         fitter.SetAmpOffset(offset);
179         fitter.SetAmpThreshold(threshold);
180       }
181       
182       while (stream.NextDDL()) {
183         while (stream.NextChannel()) {
184
185           cellX    = stream.GetCellX();
186           cellZ    = stream.GetCellZ();
187           caloFlag = stream.GetCaloFlag();  // 0=LG, 1=HG, 2=TRU
188           
189           if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
190           
191           // In case of oscillating signals with ZS, 
192           //a channel can have several bunches.
193           
194           nBunches = 0;
195           while (stream.NextBunch()) {
196             nBunches++;
197             if (nBunches > 1) continue;
198             sigStart  = stream.GetStartTimeBin();
199             sigLength = stream.GetBunchLength();
200             fitter.SetChannelGeo(stream.GetModule(),cellX,cellZ,caloFlag);
201             fitter.Eval(stream.GetSignals(),sigStart,sigLength);
202           } // End of NextBunch()
203           
204           if (nBunches>1) continue;
205           
206           e[cellX][cellZ][caloFlag] = fitter.GetEnergy();
207           t[cellX][cellZ][caloFlag] = fitter.GetTime();
208         }
209         
210         if(stream.GetModule()<0 || stream.GetModule()>4) continue;
211         
212         if(dAs[stream.GetModule()])
213           dAs[stream.GetModule()]->FillHistograms(e,t);
214         else {
215           dAs[stream.GetModule()] = new AliPHOSRcuDA1(stream.GetModule(),-1);
216           dAs[stream.GetModule()]->FillHistograms(e,t);
217         }
218         
219         for(Int_t iX=0; iX<64; iX++) {
220           for(Int_t iZ=0; iZ<56; iZ++) {
221             for(Int_t iGain=0; iGain<2; iGain++) {
222               e[iX][iZ][iGain] = 0.;
223               t[iX][iZ][iGain] = 0.;
224             }
225           }
226         }
227
228       }
229
230 //       da1.FillHistograms(e,t);
231 //     //da1.UpdateHistoFile();
232       
233         delete rawReader;     
234         nevents_physics++;
235     }
236     
237     nevents_total++;
238     
239     /* free resources */
240     free(event);
241     
242     /* exit when last event received, no need to wait for TERM signal */
243     if (eventT==END_OF_RUN) {
244       printf("EOR event detected\n");
245       break;
246     }
247   }
248   
249   for(Int_t i = 0; i < 4; i++) delete mapping[i];  
250   
251   /* Be sure that all histograms are saved */
252   
253   char localfile[128];
254   
255   for(Int_t iMod=0; iMod<5; iMod++) {
256     if(!dAs[iMod]) continue;
257     
258     dAs[iMod]->UpdateHistoFile();
259     dAs[iMod]->SetWriteToFile(kFALSE);    
260     
261     /* Store output files to the File Exchange Server */
262     sprintf(localfile,"PHOS_Module%d_Calib.root",iMod);
263     daqDA_FES_storeFile(localfile,"AMPLITUDES");
264   }
265   
266   return status;
267 }