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