]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/sim/AliHLTSimulation.cxx
fixes in standalone tracker code
[u/mrichter/AliRoot.git] / HLT / sim / AliHLTSimulation.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTSimulation.cxx
20     @author Matthias Richter
21     @date   
22     @brief  Binding class for HLT simulation in AliRoot. */
23
24 #include <cassert>
25 #include <cerrno>
26 #include "TObjArray.h"
27 #include "TObjString.h"
28 #include "AliHLTSimulation.h"
29 #include "AliSimulation.h"
30 #include "AliLog.h"
31 #include "AliRun.h"
32 #include "AliRunLoader.h"
33 #include "AliHeader.h"
34 #include "AliCDBManager.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBPath.h"
37 #include "AliCDBId.h"
38 #include "AliCDBMetaData.h"
39 #include "AliCDBStorage.h"
40 #include "AliGRPObject.h"
41 #include "AliGRPManager.h"
42 #include "AliHLTSystem.h"
43 #include "AliHLTConfigurationHandler.h"
44 #include "AliHLTPluginBase.h"
45 #include "AliRawReaderFile.h"
46 #include "AliRawReaderDate.h"
47 #include "AliRawReaderRoot.h"
48 #include "AliESDEvent.h"
49 #include "AliHLTOUTComponent.h"
50 #include "AliTracker.h"
51 #include "TGeoGlobalMagField.h"
52 #include "TSystem.h"
53 #include "TMath.h"
54 #include "TGeoGlobalMagField.h"
55
56 #if ALIHLTSIMULATION_LIBRARY_VERSION != LIBHLTSIM_VERSION
57 #error library version in header file and lib*.pkg do not match
58 #endif
59
60 /** ROOT macro for the implementation of ROOT specific class methods */
61 ClassImp(AliHLTSimulation);
62
63 AliHLTSimulation::AliHLTSimulation()
64   : fOptions()
65   , fpPluginBase(new AliHLTPluginBase)
66   , fpRawReader(NULL)
67   , fNEvents(-1)
68 {
69   // see header file for class documentation
70   // or
71   // refer to README to build package
72   // or
73   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
74 }
75
76 AliHLTSimulation::~AliHLTSimulation()
77 {
78   // see header file for function documentation
79   if (fpPluginBase) delete fpPluginBase;
80   fpPluginBase=NULL;
81
82   if (fpRawReader) {
83     delete fpRawReader;
84   }
85   fpRawReader=NULL;
86 }
87
88 AliHLTSimulation* AliHLTSimulation::CreateInstance()
89 {
90   // see header file for function documentation
91   return new AliHLTSimulation;
92 }
93
94 int AliHLTSimulation::DeleteInstance(AliHLTSimulation* pSim)
95 {
96   // see header file for function documentation
97   assert(pSim!=NULL);
98   delete pSim;
99   return 0;
100 }
101
102 int AliHLTSimulation::Init(AliRunLoader* pRunLoader, const char* options)
103 {
104   // init the simulation
105   fOptions=options;
106   TString sysOp;
107
108   if(!fpPluginBase) {
109     AliError("internal initialization failed");
110     return -EINVAL;
111   }
112
113   AliHLTSystem* pSystem=fpPluginBase->GetInstance();
114   if (!pSystem) {
115     AliError("can not get AliHLTSystem instance");
116     return -ENOMEM;
117   }
118   if (pSystem->CheckStatus(AliHLTSystem::kError)) {
119     AliError("HLT system in error state");
120     return -EFAULT;
121   }
122
123   // scan options for specific entries
124   TObjArray* pTokens=fOptions.Tokenize(" ");
125   if (pTokens) {
126     int iEntries=pTokens->GetEntries();
127     for (int i=0; i<iEntries; i++) {
128       TString token=(((TObjString*)pTokens->At(i))->GetString());
129       if (token.Contains("rawfile=")) {
130         TString param=token.ReplaceAll("rawfile=", "");
131         if (param.EndsWith("/")) {
132           AliInfo(Form("creating AliRawReaderFile (%s)", param.Data()));
133           fpRawReader = new AliRawReaderFile(param);
134         } else if (param.EndsWith(".root")) {
135           AliInfo(Form("creating AliRawReaderRoot (%s)", param.Data()));
136           fpRawReader = new AliRawReaderRoot(param);
137         } else if (!param.IsNull()) {
138           AliInfo(Form("creating AliRawReaderDate (%s)", param.Data()));
139           fpRawReader = new AliRawReaderDate(param);
140         }
141         if (fpRawReader) {
142             fpRawReader->RewindEvents();
143             int count=0;
144             for ( ; fpRawReader->NextEvent(); count++) {/* empty body */};
145             if (count!=pRunLoader->GetNumberOfEvents()) {
146               AliError(Form("mismatch in event count: runloader %d, rawreader %d; ignoring rawreader", 
147                             pRunLoader->GetNumberOfEvents(), count));
148               count=0;
149             }
150             if (count>0) {
151               fpRawReader->RewindEvents();
152               fpRawReader->NextEvent();
153             } else {
154               delete fpRawReader;
155               fpRawReader=NULL;
156             }
157         }
158       } else if (token.Contains("writerawfiles=")) {
159         if (!token.ReplaceAll("writerawfiles=", "").Contains("HLT")) {
160           if (TestBit(kOneChain) && AliHLTOUTComponent::TestGlobalOption(AliHLTOUTComponent::kWriteRawFiles)) {
161             AliWarning("empty argument 'writerawfiles=' disables HLTOUTComponent mode 'raw' which was set by argument 'hltout-mode'");
162           }
163           AliHLTOUTComponent::ClearGlobalOption(AliHLTOUTComponent::kWriteRawFiles);
164         }
165       } else if (token.BeginsWith("hltout-mode=")) {
166         // this is a legacy mode to emulate the behavior before Dec 2010 where only
167         // one chain was executed on either digits or simulated raw data and the output
168         // was controlled via global flags
169         // add to the arguments for AliHLTSystem as also there the information is needed
170         if (sysOp.Length()>0) sysOp+=" ";
171         sysOp+=token;
172         TString param=token.ReplaceAll("hltout-mode=", "");
173         SetBit(kOneChain);
174         if (param.CompareTo("raw")==0) {
175           // please note that this option
176           AliHLTOUTComponent::SetGlobalOption(AliHLTOUTComponent::kWriteRawFiles);
177           AliHLTOUTComponent::ClearGlobalOption(AliHLTOUTComponent::kWriteDigits);
178         } else if (param.CompareTo("digits")==0) {
179           // please note that this option
180           AliHLTOUTComponent::ClearGlobalOption(AliHLTOUTComponent::kWriteRawFiles);
181           AliHLTOUTComponent::SetGlobalOption(AliHLTOUTComponent::kWriteDigits);
182         } else if (param.CompareTo("legacy")==0) {
183           AliHLTOUTComponent::SetGlobalOption(AliHLTOUTComponent::kWriteRawFiles);
184           AliHLTOUTComponent::SetGlobalOption(AliHLTOUTComponent::kWriteDigits);
185         } else {
186           AliError(Form("invalid parameter for argument 'hltout-mode=' %s, allowed: raw, digits, legacy ... ignoring argument  and using the standard simulation", param.Data()));
187           ResetBit(kOneChain);
188         }
189       } else if (token.Contains("events=")) {
190         fNEvents=token.ReplaceAll("events=", "").Atoi();
191       } else {
192         if (sysOp.Length()>0) sysOp+=" ";
193         sysOp+=token;
194       }
195     }
196     delete pTokens;
197   }
198   // only store the options for AliHLTSystem
199   fOptions=sysOp;
200
201   // if no specific hltout-mode has been chosen set the split mode for
202   // running separate chains for digits and raw data
203   if (!fOptions.Contains("hltout-mode=")) fOptions+=" hltout-mode=split";
204
205   AliCDBManager* man = AliCDBManager::Instance();
206   if (man && man->IsDefaultStorageSet())
207   {
208     // init solenoid field
209     // 2009-11-07 magnetic field handling fo HLT components has been switched to the
210     // global AliMagF instance, the HLT/ConfigHLT/SolenoidBz entry is obsolete
211     // The global instance is either established by the AliRoot environment or the
212     // component external interface.
213     if (TGeoGlobalMagField::Instance()->GetField()) {
214       AliDebug(0, Form("magnetic field: %f", AliTracker::GetBz()));
215     } else {
216       // workaround for bug #51285
217       AliGRPManager grpman;
218       if (grpman.ReadGRPEntry() &&
219           grpman.SetMagField()) {
220         // nothing to do any more
221       }
222       AliError(Form("can not get the AliMagF instance, falling back to GRP entry (%f)", AliTracker::GetBz()));
223     }
224   } else if (man) {
225     AliError("OCDB default storage not yet set, can not prepare OCDB entries");    
226   } else {
227     AliError("unable to get instance of AliCDBMetaData, can not prepare OCDB entries");    
228   }
229
230   // configure the main HLTSystem instance for digit simulation (pRawReader NULL)
231   return ConfigureHLTSystem(pSystem, fOptions.Data(), pRunLoader, TestBit(kOneChain)?fpRawReader:NULL);
232 }
233
234 int AliHLTSimulation::ConfigureHLTSystem(AliHLTSystem* pSystem, const char* options, AliRunLoader* pRunLoader, AliRawReader* pRawReader) const
235 {
236   // scan options and configure AliHLTSystem
237   if (pSystem->ScanOptions(options)<0) {
238     AliError("error setting options for HLT system");
239     return -EINVAL;     
240   }
241
242   if (!pSystem->CheckStatus(AliHLTSystem::kReady)) {
243     if ((pSystem->Configure(pRawReader, pRunLoader))<0) {
244       AliError("error during HLT system configuration");
245       return -EFAULT;
246     }
247   }
248
249   return 0;
250 }
251
252 int AliHLTSimulation::Run(AliRunLoader* pRunLoader)
253 {
254   // HLT reconstruction for simulated data  
255   if(!fpPluginBase) {
256     AliError("internal initialization failed");
257     return -EINVAL;
258   }
259
260   if(!pRunLoader) {
261     AliError("Missing RunLoader! 0x0");
262     return -EINVAL;
263   }
264
265   int iResult=0;
266
267   AliHLTSystem* pSystem=fpPluginBase->GetInstance();
268   if (!pSystem) {
269     AliError("can not get AliHLTSystem instance");
270     return -ENOMEM;
271   }
272
273   if (pSystem->CheckStatus(AliHLTSystem::kError)) {
274     AliError("HLT system in error state");
275     return -EFAULT;
276   }
277
278   // run the main HLTSystem instance for digit simulation (pRawReader NULL)
279   // in legacy mode only one chain is run and the output is controlled via
280   // global flags
281   if (!TestBit(kOneChain)) AliInfo("running HLT simulation for digits");
282   iResult=RunHLTSystem(pSystem, pRunLoader, TestBit(kOneChain)?fpRawReader:NULL);
283
284   // now run once again with the raw data as input, a completely new HLT system
285   // with new configurations is used
286   if (fpRawReader && !TestBit(kOneChain)) {
287     AliInfo("running HLT simulation for raw data");
288     int iLocalResult=0;
289     AliHLTConfigurationHandler* confHandler=new AliHLTConfigurationHandler;
290     // note that the configuration handler is owned by the
291     // AliHLTSystem instance from now on
292     AliHLTSystem rawSimulation(kHLTLogDefault, "", NULL, confHandler);
293     if ((iLocalResult=ConfigureHLTSystem(&rawSimulation, fOptions.Data(), pRunLoader, fpRawReader))>=0) {
294       iLocalResult=RunHLTSystem(&rawSimulation, pRunLoader, fpRawReader);
295     }
296     if (iResult>=0) iResult=iLocalResult;
297   }
298
299   return iResult;
300 }
301
302 int AliHLTSimulation::RunHLTSystem(AliHLTSystem* pSystem, AliRunLoader* pRunLoader, AliRawReader* pRawReader) const
303 {
304   // run reconstruction cycle for AliHLTSystem
305   int nEvents = (fNEvents<0 || fNEvents>pRunLoader->GetNumberOfEvents())?pRunLoader->GetNumberOfEvents():fNEvents;
306   int iResult=0;
307
308   // Note: the rawreader is already placed at the first event
309   if ((iResult=pSystem->Reconstruct(1, pRunLoader, pRawReader))>=0) {
310     pSystem->FillESD(0, pRunLoader, NULL);
311     for (int i=1; i<nEvents; i++) {
312       if (pRawReader && !pRawReader->NextEvent()) {
313         AliError("mismatch in event count, rawreader corrupted");
314         break;
315       }
316       pSystem->Reconstruct(1, pRunLoader, pRawReader);
317       pSystem->FillESD(i, pRunLoader, NULL);
318     }
319     // send specific 'event' to execute the stop sequence
320     pSystem->Reconstruct(0, NULL, NULL);
321   }
322   return iResult;
323 }
324
325 AliHLTSimulation* AliHLTSimulationCreateInstance()
326 {
327   // see header file for function documentation
328   return AliHLTSimulation::CreateInstance();
329 }
330
331 int AliHLTSimulationDeleteInstance(AliHLTSimulation* pSim)
332 {
333   // see header file for function documentation
334   return AliHLTSimulation::DeleteInstance(pSim);
335 }
336
337 int AliHLTSimulationInit(AliHLTSimulation* pSim, AliRunLoader* pRunLoader, const char* options)
338 {
339   assert(pSim!=NULL);
340   if (pSim) {
341     return pSim->Init(pRunLoader, options);
342   }
343   return -ENODEV;
344 }
345
346 int AliHLTSimulationRun(AliHLTSimulation* pSim, AliRunLoader* pRunLoader)
347 {
348   assert(pSim!=NULL);
349   if (pSim) {
350     return pSim->Run(pRunLoader);
351   }
352   return -ENODEV;
353 }
354
355 int AliHLTSimulationGetLibraryVersion()
356 {
357   // see header file for function documentation
358   return LIBHLTSIM_VERSION;
359 }
360
361 int AliHLTSimulationSetup(AliHLTSimulation* /*pHLTSim*/, AliSimulation* pSim, const char* specificObjects)
362 {
363   // see header file for function documentation
364
365   // this is an attempt to solve issue #48360
366   // since there are many jobs running in parallel during the production,
367   // all the jobs want to put entries into the OCDB. The solution is to
368   // make them temporary, since they are only used to propagate information
369   // from the simulation to the reconstruction.
370
371   if (!pSim) return -EINVAL;
372   const char* entries[]={
373     NULL
374   };
375
376   TString specificStorage; 
377   specificStorage.Form("local://%s",gSystem->pwd());
378   for (const char** pEntry=entries; *pEntry!=NULL; pEntry++) {
379     const char* pObject=specificObjects?strstr(specificObjects, *pEntry):NULL;
380     if (pObject) {
381       // skip this entry if it is found in the list and either
382       // last one or separated by a blank
383       pObject+=strlen(*pEntry);
384       if (*pObject==0 || *pObject==' ') continue;
385     }
386     pSim->SetSpecificStorage(*pEntry, specificStorage.Data());
387   }
388
389   return 0;
390 }
391
392 #ifndef HAVE_COMPILEINFO
393 extern "C" void CompileInfo(const char*& date, const char*& time)
394 {
395   // the fall back compile info of the HLTsim library
396   // this is not up-to-date if other files have been changed and recompiled
397   date=__DATE__; time=__TIME__;
398   return;
399 }
400 #endif