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