]>
Commit | Line | Data |
---|---|---|
eb35e591 | 1 | /* |
2 | Contact: cvetan.cheshkov@cern.ch | |
049f9cac | 3 | Link: http://alisoft.cern.ch/viewvc/trunk/ITS/ITSSPDVertexDiamondda.cxx?root=AliRoot&view=log , /afs/cern.ch/user/c/cheshkov/public/08000058338016.30.root.date.gz |
6bf5b296 | 4 | Reference Run: 58338 |
eb35e591 | 5 | Run Type: PHYSICS |
6 | DA Type: MON | |
049f9cac | 7 | Number of events needed: 100 |
6bf5b296 | 8 | Input Files: GRP/Geometry/Data , ITS/Align/Data , spd_noisy_ocdb , spd_dead_ocdb |
9 | Output Files: SPDVertexDiamondDA.root | |
eb35e591 | 10 | Trigger types used: PHYSICS |
11 | */ | |
12 | ||
eb35e591 | 13 | #define OUTPUT_FILE "SPDVertexDiamondDA.root" |
eb35e591 | 14 | #define N_EVENTS_AUTOSAVE 50 |
eb35e591 | 15 | |
16 | extern "C" { | |
17 | #include "daqDA.h" | |
18 | } | |
19 | ||
20 | #include "event.h" | |
21 | #include "monitor.h" | |
22 | ||
55c5e86d | 23 | #ifdef ALI_AMORE |
24 | #include <AmoreDA.h> | |
308c2f7c | 25 | //int amore::da::Updated(char const*) {} |
55c5e86d | 26 | #endif |
27 | ||
eb35e591 | 28 | #include <TPluginManager.h> |
29 | #include <TROOT.h> | |
cf0b1cea | 30 | #include <TH1.h> |
31 | #include <TH2.h> | |
6bf5b296 | 32 | #include <TSystem.h> |
eb35e591 | 33 | |
34 | #include "AliRawReaderDate.h" | |
eb35e591 | 35 | #include "AliCDBManager.h" |
cf0b1cea | 36 | #include "AliITSMeanVertexer.h" |
eb35e591 | 37 | |
38 | int main(int argc, char **argv) { | |
39 | ||
40 | gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo", | |
41 | "*", | |
42 | "TStreamerInfo", | |
43 | "RIO", | |
44 | "TStreamerInfo()"); | |
eb35e591 | 45 | |
46 | int status; | |
47 | if (argc<2) { | |
48 | printf("Wrong number of arguments\n"); | |
49 | return -1; | |
50 | } | |
51 | ||
52 | /* define data source : this is argument 1 */ | |
53 | status=monitorSetDataSource( argv[1] ); | |
54 | if (status!=0) { | |
55 | printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status)); | |
56 | return -1; | |
57 | } | |
58 | ||
59 | /* declare monitoring program */ | |
60 | status=monitorDeclareMp( __FILE__ ); | |
61 | if (status!=0) { | |
62 | printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status)); | |
63 | return -1; | |
64 | } | |
65 | ||
66 | /* define wait event timeout - 1s max */ | |
67 | monitorSetNowait(); | |
68 | monitorSetNoWaitNetworkTimeout(1000); | |
69 | ||
70 | /* log start of process */ | |
71 | printf("Vertex-Diamond SPD DA started\n"); | |
72 | ||
73 | /* init some counters */ | |
74 | int nevents_with_vertex = 0; | |
75 | int nevents_physics=0; | |
76 | int nevents_total=0; | |
77 | ||
78 | struct eventHeaderStruct *event; | |
79 | eventTypeType eventT; | |
80 | ||
81 | // Get run number | |
82 | if (getenv("DATE_RUN_NUMBER")==0) { | |
83 | printf("DATE_RUN_NUMBER not properly set.\n"); | |
84 | return -1; | |
85 | } | |
86 | int runNr = atoi(getenv("DATE_RUN_NUMBER")); | |
87 | ||
6bf5b296 | 88 | // Get the necessary OCDB files from the DAQ detector DB |
89 | if (gSystem->AccessPathName("localOCDB/GRP/Geometry/Data",kFileExists)) { | |
90 | if (gSystem->mkdir("localOCDB/GRP/Geometry/Data",kTRUE) != 0) { | |
91 | printf("Failed to create directory: localOCDB/GRP/Geometry/Data"); | |
92 | return -1; | |
93 | } | |
94 | } | |
95 | status = daqDA_DB_getFile("GRP/Geometry/Data","localOCDB/GRP/Geometry/Data/Run0_999999999_v0_s0.root"); | |
96 | if (status) { | |
97 | printf("Failed to get geometry file (GRP/Geometry/Data) from DAQdetDB, status=%d\n", status); | |
98 | return -1; | |
99 | } | |
100 | ||
101 | if (gSystem->AccessPathName("localOCDB/ITS/Align/Data",kFileExists)) { | |
102 | if (gSystem->mkdir("localOCDB/ITS/Align/Data",kTRUE) != 0) { | |
103 | printf("Failed to create directory: localOCDB/ITS/Align/Data"); | |
104 | return -1; | |
105 | } | |
106 | } | |
107 | status = daqDA_DB_getFile("ITS/Align/Data","localOCDB/ITS/Align/Data/Run0_999999999_v0_s0.root"); | |
108 | if (status) { | |
109 | printf("Failed to get its-alignment file (ITS/Align/Data) from DAQdetDB, status=%d\n", status); | |
110 | return -1; | |
111 | } | |
112 | ||
113 | if (gSystem->AccessPathName("localOCDB/ITS/Calib/SPDNoisy",kFileExists)) { | |
114 | if (gSystem->mkdir("localOCDB/ITS/Calib/SPDNoisy",kTRUE) != 0) { | |
115 | printf("Failed to create directory: localOCDB/ITS/Calib/SPDNoisy"); | |
116 | return -1; | |
117 | } | |
118 | } | |
119 | status = daqDA_DB_getFile("spd_noisy_ocdb","localOCDB/ITS/Calib/SPDNoisy/Run0_999999999_v0_s0.root"); | |
120 | if (status) { | |
121 | printf("Failed to get spd file (spd_noisy_ocdb) from DAQdetDB, status=%d\n", status); | |
122 | return -1; | |
123 | } | |
124 | ||
125 | if (gSystem->AccessPathName("localOCDB/ITS/Calib/SPDDead",kFileExists)) { | |
126 | if (gSystem->mkdir("localOCDB/ITS/Calib/SPDDead",kTRUE) != 0) { | |
127 | printf("Failed to create directory: localOCDB/ITS/Calib/SPDDead"); | |
128 | return -1; | |
129 | } | |
130 | } | |
131 | status = daqDA_DB_getFile("spd_dead_ocdb","localOCDB/ITS/Calib/SPDDead/Run0_999999999_v0_s0.root"); | |
132 | if (status) { | |
133 | printf("Failed to get spd file (spd_dead_ocdb) from DAQdetDB, status=%d\n", status); | |
134 | return -1; | |
135 | } | |
136 | ||
eb35e591 | 137 | // Global initializations |
138 | AliCDBManager *man = AliCDBManager::Instance(); | |
6bf5b296 | 139 | man->SetDefaultStorage("local://localOCDB"); |
eb35e591 | 140 | man->SetRun(runNr); |
eb35e591 | 141 | |
37add1d1 | 142 | // Init mean vertexer |
143 | AliITSMeanVertexer *mv = new AliITSMeanVertexer(); | |
144 | if (!mv->Init()) { | |
145 | printf("Initialization of mean vertexer object failed ! Check the log for details"); | |
146 | return -1; | |
147 | } | |
eb35e591 | 148 | |
55c5e86d | 149 | // Initialization of AMORE sender |
150 | #ifdef ALI_AMORE | |
151 | amore::da::AmoreDA vtxAmore(amore::da::AmoreDA::kSender); | |
152 | #endif | |
eb35e591 | 153 | /* main loop (infinite) */ |
154 | for(;;) { | |
155 | ||
156 | /* check shutdown condition */ | |
157 | if (daqDA_checkShutdown()) {break;} | |
158 | ||
159 | /* get next event (blocking call until timeout) */ | |
160 | status=monitorGetEventDynamic((void **)&event); | |
161 | if (status==MON_ERR_EOF) { | |
162 | printf ("End of File detected\n"); | |
163 | break; /* end of monitoring file has been reached */ | |
164 | } | |
165 | ||
166 | if (status!=0) { | |
167 | printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status)); | |
168 | break; | |
169 | } | |
170 | ||
171 | /* retry if got no event */ | |
172 | if (event==NULL) { | |
173 | continue; | |
174 | } | |
175 | ||
176 | nevents_total++; | |
177 | eventT=event->eventType; | |
178 | switch (event->eventType){ | |
179 | ||
180 | /* START OF RUN */ | |
181 | case START_OF_RUN: | |
182 | break; | |
183 | /* END START OF RUN */ | |
184 | ||
185 | /* END OF RUN */ | |
186 | case END_OF_RUN: | |
187 | break; | |
188 | ||
189 | case PHYSICS_EVENT: | |
190 | nevents_physics++; | |
191 | AliRawReader *rawReader = new AliRawReaderDate((void*)event); | |
192 | ||
37add1d1 | 193 | // Run mean-vertexer reco |
194 | if (mv->Reconstruct(rawReader)) nevents_with_vertex++; | |
eb35e591 | 195 | |
37add1d1 | 196 | // Auto save |
6bf5b296 | 197 | if ((nevents_physics%N_EVENTS_AUTOSAVE) == 0) { |
37add1d1 | 198 | mv->WriteVertices(OUTPUT_FILE); |
eb35e591 | 199 | |
6bf5b296 | 200 | #ifdef ALI_AMORE |
201 | // send the histos to AMORE pool | |
202 | printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexXY()->GetName(),mv->GetVertexXY())); | |
203 | printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexZ()->GetName(),mv->GetVertexZ())); | |
204 | #endif | |
205 | } | |
206 | ||
eb35e591 | 207 | delete rawReader; |
eb35e591 | 208 | } |
209 | ||
210 | /* free resources */ | |
211 | free(event); | |
212 | ||
213 | /* exit when last event received, no need to wait for TERM signal */ | |
214 | if (eventT==END_OF_RUN) { | |
215 | printf("EOR event detected\n"); | |
216 | break; | |
217 | } | |
218 | } | |
eb35e591 | 219 | |
37add1d1 | 220 | mv->WriteVertices(OUTPUT_FILE); |
eb35e591 | 221 | |
55c5e86d | 222 | #ifdef ALI_AMORE |
223 | // send the histos to AMORE pool | |
cf0b1cea | 224 | printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexXY()->GetName(),mv->GetVertexXY())); |
225 | printf("AMORE send status: %d\n",vtxAmore.Send(mv->GetVertexZ()->GetName(),mv->GetVertexZ())); | |
55c5e86d | 226 | #endif |
227 | ||
37add1d1 | 228 | delete mv; |
eb35e591 | 229 | |
230 | /* write report */ | |
231 | 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); | |
232 | ||
233 | status=0; | |
234 | ||
235 | /* export file to FXS */ | |
236 | if (daqDA_FES_storeFile(OUTPUT_FILE, "VertexDiamond")) { | |
237 | status=-2; | |
238 | } | |
239 | ||
240 | return status; | |
241 | } |