]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/ITSSPDVertexDiamondda.cxx
During simulation: fill STU region w/ non null time sums
[u/mrichter/AliRoot.git] / ITS / ITSSPDVertexDiamondda.cxx
1 /*
2   Contact: cvetan.cheshkov@cern.ch
3   Link: http://alisoft.cern.ch/viewvc/trunk/ITS/ITSSPDVertexDiamondda.cxx?root=AliRoot&view=log , /afs/cern.ch/user/c/cheshkov/public/ITS/VD_da_test.date , /afs/cern.ch/user/c/cheshkov/public/08000058338016.30.root.date.gz
4   Reference Run: 58338
5   Run Type: PHYSICS
6   DA Type: MON
7   Number of events needed: 100
8   Input Files: GRP/Geometry/Data , ITS/Align/Data , spd_noisy_ocdb , spd_dead_ocdb , TRIGGER/SPD/PITConditions (all the files are taken from SPD daqDetDB)
9   Output Files: SPDVertexDiamondDA.root
10   Trigger types used: PHYSICS, SPD-F0 
11 */
12
13 #define OUTPUT_FILE "SPDVertexDiamondDA.root"
14
15 extern "C" {
16 #include "daqDA.h"
17 }
18
19 #include "event.h"
20 #include "monitor.h"
21
22 #ifdef ALI_AMORE
23 #include <AmoreDA.h>
24 //int amore::da::Updated(char const*) {}
25 #endif
26
27 #include <TPluginManager.h>
28 #include <TROOT.h>
29 #include <TH1.h>
30 #include <TH2.h>
31 #include <TPaveText.h>
32 #include <TSystem.h>
33 #include <TGeoGlobalMagField.h>
34
35 #include "AliLog.h"
36 #include "AliMagF.h"
37 #include "AliRawReaderDate.h"
38 #include "AliCDBManager.h"
39 #include "AliITSMeanVertexer.h"
40
41 int main(int argc, char **argv) {
42
43   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
44                                         "*",
45                                         "TStreamerInfo",
46                                         "RIO",
47                                         "TStreamerInfo()"); 
48
49   int status;
50   if (argc<2) {
51     printf("Wrong number of arguments\n");
52     return -1;
53   }
54   // get mean vertex quality cuts
55   status = daqDA_DB_getFile("ITSSPD_VertexQualityTuning_DA.config","./ITSSPD_VertexQualityTuning_DA.config");
56   if (status) {
57     printf("Failed to get config file (ITSSPD_VertexQualityTuning_DA.config) from DAQ DB, status=%d\n", status);
58     return -1;
59   }
60   /* open the config file and retrieve running parameters */
61   Float_t errX  = -1;
62   Float_t r     = -1;
63   UInt_t minClInner = 999 ;
64   UInt_t maxClInner = 999 ;
65   Int_t nEvFirstLoop = 0;
66   Int_t nEvAUTOSAVE = 0; 
67   char name[10][10];
68
69   FILE *fpConfig = fopen("ITSSPD_VertexQualityTuning_DA.config","r");
70   fscanf(fpConfig,"%s %f\n %s %f\n %s %d\n %s %d \n %s %d \n %s %d ",&name[0], &errX, &name[1], &r, &name[2], &minClInner,&name[3], &maxClInner, &name[4],&nEvFirstLoop,&name[5],&nEvAUTOSAVE);
71   fclose(fpConfig);
72
73   printf("\n\n Mean Vertex quality cuts : \n- errX = %f\n- r = %f\n- minSPD0 = %d maxSPD0 = %d\n- nEventsFirstLoop = %d nEventsAUTOSAVE = %d \n\n\n",errX,r,minClInner,maxClInner,nEvFirstLoop,nEvAUTOSAVE);
74
75   /* define data source : this is argument 1 */  
76   status=monitorSetDataSource( argv[1] );
77   if (status!=0) {
78     printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
79     return -1;
80   }
81
82   /* declare monitoring program */
83   status=monitorDeclareMp( __FILE__ );
84   if (status!=0) {
85     printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
86     return -1;
87   }
88
89   /* define wait event timeout - 1s max */
90   monitorSetNowait();
91   monitorSetNoWaitNetworkTimeout(1000);
92   
93   /* log start of process */
94   printf("Vertex-Diamond SPD DA started\n");  
95
96   /* init some counters */
97   int nevents_with_vertex = 0;
98   int nevents_physics=0;
99   int nevents_total=0;
100
101   struct eventHeaderStruct *event;
102   eventTypeType eventT;
103
104   // Get run number
105   if (getenv("DATE_RUN_NUMBER")==0) {
106     printf("DATE_RUN_NUMBER not properly set.\n");
107     return -1;
108   }
109   int runNr = atoi(getenv("DATE_RUN_NUMBER"));
110
111   // Get the necessary OCDB files from the DAQ detector DB
112   if (gSystem->AccessPathName("localOCDB/GRP/Geometry/Data",kFileExists)) {
113     if (gSystem->mkdir("localOCDB/GRP/Geometry/Data",kTRUE) != 0) {
114       printf("Failed to create directory: localOCDB/GRP/Geometry/Data");
115       return -1;
116     }
117   }
118   status = daqDA_DB_getFile("GRP/Geometry/Data","localOCDB/GRP/Geometry/Data/Run0_999999999_v0_s0.root");
119   if (status) {
120     printf("Failed to get geometry file (GRP/Geometry/Data) from DAQdetDB, status=%d\n", status);
121     return -1;
122   }
123
124   if (gSystem->AccessPathName("localOCDB/ITS/Align/Data",kFileExists)) {
125     if (gSystem->mkdir("localOCDB/ITS/Align/Data",kTRUE) != 0) {
126       printf("Failed to create directory: localOCDB/ITS/Align/Data");
127       return -1;
128     }
129   }
130   status = daqDA_DB_getFile("ITS/Align/Data","localOCDB/ITS/Align/Data/Run0_999999999_v0_s0.root");
131   if (status) {
132     printf("Failed to get its-alignment file (ITS/Align/Data) from DAQdetDB, status=%d\n", status);
133     return -1;
134   }
135
136   if (gSystem->AccessPathName("localOCDB/ITS/Calib/SPDNoisy",kFileExists)) {
137     if (gSystem->mkdir("localOCDB/ITS/Calib/SPDNoisy",kTRUE) != 0) {
138       printf("Failed to create directory: localOCDB/ITS/Calib/SPDNoisy");
139       return -1;
140     }
141   }
142   status = daqDA_DB_getFile("spd_noisy_ocdb","localOCDB/ITS/Calib/SPDNoisy/Run0_999999999_v0_s0.root");
143   if (status) {
144     printf("Failed to get spd file (spd_noisy_ocdb) from DAQdetDB, status=%d\n", status);
145     return -1;
146   }
147
148   if (gSystem->AccessPathName("localOCDB/ITS/Calib/SPDDead",kFileExists)) {
149     if (gSystem->mkdir("localOCDB/ITS/Calib/SPDDead",kTRUE) != 0) {
150       printf("Failed to create directory: localOCDB/ITS/Calib/SPDDead");
151       return -1;
152     }
153   }
154   status = daqDA_DB_getFile("spd_dead_ocdb","localOCDB/ITS/Calib/SPDDead/Run0_999999999_v0_s0.root");
155   if (status) {
156     printf("Failed to get spd file (spd_dead_ocdb) from DAQdetDB, status=%d\n", status);
157     return -1;
158   }
159
160   if (gSystem->AccessPathName("localOCDB/ITS/Calib/SPDSparseDead",kFileExists)) {
161     if (gSystem->mkdir("localOCDB/ITS/Calib/SPDSparseDead",kTRUE) != 0) {
162       printf("Failed to create directory: localOCDB/ITS/Calib/SPDSparseDead");
163       return -1;
164     }
165   }
166
167   status = daqDA_DB_getFile("spd_sparsedead_ocdb","localOCDB/ITS/Calib/SPDSparseDead/Run0_999999999_v0_s0.root");
168   if (status) {
169     printf("Failed to get spd file (spd_sparsedead_ocdb) from DAQdetDB, status=%d\n", status);
170     return -1;
171   }
172
173   if (gSystem->AccessPathName("localOCDB/TRIGGER/SPD/PITConditions",kFileExists)) {
174     if (gSystem->mkdir("localOCDB/TRIGGER/SPD/PITConditions",kTRUE) != 0) {
175       printf("Failed to create directory: localOCDB/TRIGGER/SPD/PITConditions");
176       return -1;
177     }
178   }
179   status = daqDA_DB_getFile("TRIGGER/SPD/PITConditions","localOCDB/TRIGGER/SPD/PITConditions/Run0_999999999_v0_s0.root");
180   if (status) {
181     printf("Failed to get spd trigger file (TRIGGER/SPD/PITConditions) from DAQdetDB, status=%d\n", status);
182     return -1;
183   }
184  
185   status = daqDA_DB_getFile("mfchebKGI_sym.root","localOCDB/mfchebKGI_sym.root");
186   if (status) {
187     printf("Failed to get spd file (mfchebKGI_sym.root) from DAQdetDB, status=%d\n", status);
188     return -1;
189   }
190
191   // Global initializations
192
193   // The B filed is required in AliITSClusterFinderV2SPD
194   // for the Lorentz angle correction. B set to 0.      
195   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 0., 0., AliMagF::k5kGUniform,AliMagF::kBeamTypeAA,-1,2,15,"localOCDB/mfchebKGI_sym.root"));  
196   AliLog::SetGlobalLogLevel(AliLog::kError);
197   AliCDBManager *man = AliCDBManager::Instance();
198   man->SetDefaultStorage("local://localOCDB");
199   man->SetRun(runNr);
200
201   // Init mean vertexer
202   AliITSMeanVertexer *mv = new AliITSMeanVertexer();
203   if (!mv->Init()) {
204     printf("Initialization of mean vertexer object failed ! Check the log for details");
205     return -1;
206   }
207
208   mv->SetCutOnErrX(errX);
209   mv->SetCutOnR(r);
210   mv->SetCutOnCls(minClInner,maxClInner);
211
212   // Initialization of AMORE sender
213 #ifdef ALI_AMORE
214   amore::da::AmoreDA vtxAmore(amore::da::AmoreDA::kSender);
215 #endif
216   /* main loop (infinite) */
217   for(;;) {
218     
219     /* check shutdown condition */
220     if (daqDA_checkShutdown()) {break;}
221     
222     /* get next event (blocking call until timeout) */
223     status=monitorGetEventDynamic((void **)&event);
224     if (status==MON_ERR_EOF) {
225       printf ("End of File detected\n");
226       break; /* end of monitoring file has been reached */
227     }
228     
229     if (status!=0) {
230       printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
231       break;
232     }
233     
234     /* retry if got no event */
235     if (event==NULL) {
236       continue;
237     }
238
239     nevents_total++;
240     eventT=event->eventType;
241     switch (event->eventType){
242       
243       /* START OF RUN */
244     case START_OF_RUN:
245       break;
246       /* END START OF RUN */
247       
248       /* END OF RUN */
249     case END_OF_RUN:
250       break;
251       
252     case PHYSICS_EVENT:
253       nevents_physics++;
254       AliRawReader *rawReader = new AliRawReaderDate((void*)event);
255
256       // Run mean-vertexer reco
257       if (mv->Reconstruct(rawReader)) nevents_with_vertex++;
258       if(nevents_physics < nEvFirstLoop) continue;
259       // Auto save
260       if ((nevents_physics%nEvAUTOSAVE) == 0) {
261         TH2F *histo = ((TH2F*)mv->GetVertexXY());
262         histo->SetStats(kFALSE);
263         histo->SetTitle("");
264         histo->GetListOfFunctions()->SetOwner(kTRUE);
265         if(histo->GetListOfFunctions()->GetEntries()<1) histo->GetListOfFunctions()->Add(new TPaveText(-5,4.5,5.,6.2,"br"));
266         for(Int_t i=0; i<histo->GetListOfFunctions()->GetEntries(); i++){
267           TString funcName = histo->GetListOfFunctions()->At(i)->ClassName();
268           if(funcName.Contains("TPaveText")){
269             TPaveText *p = (TPaveText*)histo->GetListOfFunctions()->At(i);
270             p->Clear();
271             p->AddText(Form("%f events with vertex (%i out of %i processed events)",((Double_t)nevents_with_vertex)/((Double_t)nevents_physics),nevents_with_vertex,nevents_physics));
272             p->AddText(Form("%f events with good vertex (%i out of %i events with vertex)",histo->GetEntries()/((Double_t)nevents_with_vertex),(Int_t)histo->GetEntries(),nevents_with_vertex));
273           }
274         }
275         mv->WriteVertices(OUTPUT_FILE);
276
277 #ifdef ALI_AMORE
278         // send the histos to AMORE pool
279         printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexXY()->GetName(),mv->GetVertexXY()));
280         printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexZ()->GetName(),mv->GetVertexZ()));
281 #endif
282       }
283
284       delete rawReader;
285     }
286     
287     /* free resources */
288     free(event);
289     
290     /* exit when last event received, no need to wait for TERM signal */
291     if (eventT==END_OF_RUN) {
292       printf("EOR event detected\n");
293       break;
294     }
295   }
296
297   TH2F *histo = ((TH2F*)mv->GetVertexXY());
298   histo->SetStats(kFALSE);
299   histo->SetTitle("");
300   histo->GetListOfFunctions()->SetOwner(kTRUE);
301   if(histo->GetListOfFunctions()->GetEntries()<1) histo->GetListOfFunctions()->Add(new TPaveText(-5,4.5,5.,6.2,"br"));
302   for(Int_t i=0; i<histo->GetListOfFunctions()->GetEntries(); i++){
303     TString funcName = histo->GetListOfFunctions()->At(i)->ClassName();
304     if(funcName.Contains("TPaveText")){
305       TPaveText *p = (TPaveText*)histo->GetListOfFunctions()->At(i);
306       p->Clear(); 
307       p->AddText(Form("%f events with vertex (%i out of %i processed events)",((Double_t)nevents_with_vertex)/((Double_t)nevents_physics),nevents_with_vertex,nevents_physics));
308       p->AddText(Form("%f events with good vertex (%i out of %i events with vertex)",histo->GetEntries()/((Double_t)nevents_with_vertex),(Int_t)histo->GetEntries(),nevents_with_vertex));
309     }
310   }
311   mv->WriteVertices(OUTPUT_FILE);
312
313 #ifdef ALI_AMORE
314   // send the histos to AMORE pool
315   printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexXY()->GetName(),mv->GetVertexXY()));
316   printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexZ()->GetName(),mv->GetVertexZ()));
317 #endif
318
319   delete mv;
320
321   /* write report */
322   printf("Run #%s, received %d events with vertex, out of %d physics and out of %d total events\n",getenv("DATE_RUN_NUMBER"),nevents_with_vertex,nevents_physics,nevents_total);
323
324   status=0;
325
326   /* export file to FXS */
327   if (daqDA_FES_storeFile(OUTPUT_FILE, "VertexDiamond")) {
328     status=-2;
329   }
330   
331   return status;
332 }