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