]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/icepack/iceconvert/IceF2k.cxx
24-jun-2005 NvE Selection/overruling of specific job folder creation introduced in...
[u/mrichter/AliRoot.git] / RALICE / icepack / iceconvert / IceF2k.cxx
1 /*******************************************************************************
2  * Copyright(c) 2003, IceCube Experiment at the South Pole. All rights reserved.
3  *
4  * Author: The IceCube RALICE-based Offline Project.
5  * Contributors are mentioned in the code where appropriate.
6  *
7  * Permission to use, copy, modify and distribute this software and its
8  * documentation strictly for non-commercial purposes is hereby granted
9  * without fee, provided that the above copyright notice appears in all
10  * copies and that both the copyright notice and this permission notice
11  * appear in the supporting documentation.
12  * The authors make no claims about the suitability of this software for
13  * any purpose. It is provided "as is" without express or implied warranty.
14  *******************************************************************************/
15
16 // $Id$
17
18 ///////////////////////////////////////////////////////////////////////////
19 // Class IceF2k
20 // Conversion of Amanda F2K data into IceEvent physics event structures.
21 // This class is derived from AliJob providing a task-based processing
22 // structure on an event-by-event basis.
23 // The main object in the job environment is an IceEvent* pointer.
24 // In case the user has provided sub-tasks, these will be executed
25 // on an event-by-event basis after the IceEvent structure has been filled
26 // with the F2K data and before the final structures are written out.
27 // Note that the data structures are only written out if an outputfile has
28 // been specified via the SetOutputFile memberfunction.
29 // In case no outputfile has been specified, this class provides a facility
30 // to investigate/analyse F2K data using the Ralice/IcePack analysis tools. 
31 //
32 // Usage example :
33 // ---------------
34 //
35 // gSystem->Load("ralice");
36 // gSystem->Load("icepack");
37 // gSystem->Load("iceconvert");
38 //
39 // IceF2k q("IceF2k","F2K to IcePack data structure conversion");
40 //
41 // // Limit the number of entries for testing
42 // q.SetMaxEvents(10);
43 //
44 // // Print frequency to produce a short summary print every printfreq events
45 // q.SetPrintFreq(1);
46 //
47 // // Split level for the output structures
48 // q.SetSplitLevel(2);
49 //
50 // // Buffer size for the output structures
51 // q.SetBufferSize(32000);
52 //
53 // // The F2K input filename
54 // q.SetInputFile("run7825.f2k");
55 //
56 // // Output file for the event structures
57 // TFile* ofile=new TFile("events.root","RECREATE","F2K data in IceEvent structure");
58 // q.SetOutputFile(ofile);
59 //
60 // ///////////////////////////////////////////////////////////////////
61 // // Here the user can specify his/her sub-tasks to be executed
62 // // on an event-by-event basis after the IceEvent structure
63 // // has been filled and before the data is written out.
64 // // Sub-tasks (i.e. a user classes derived from TTask) are entered
65 // // as follows :
66 // //
67 // //    MyXtalk task1("task1","Cross talk correction");
68 // //    MyClean task2("task2","Hit cleaning");
69 // //    q.Add(&task1);
70 // //    q.Add(&task2);
71 // //
72 // // The sub-tasks will be executed in the order as they are entered.
73 // ///////////////////////////////////////////////////////////////////
74 //
75 // // Perform the conversion and execute subtasks (if any)
76 // // on an event-by-event basis
77 // q.ExecuteJob();
78 //
79 // // Select various objects to be added to the output file
80 //
81 // ofile->cd(); // Switch to the output file directory
82 //
83 // AliObjMatrix* omdb=q.GetOMdbase();
84 // if (omdb) omdb->Write();
85 //
86 // AliDevice* fitdefs=q.GetFitdefs();
87 // if (fitdefs) fitdefs->Write();
88 //
89 // TDatabasePDG* pdg=q.GetPDG();
90 // if (pdg) pdg->Write();
91 //
92 // // Flush the output file.
93 // // The output file is not explicitly closed here
94 // // to allow ineractive investigation of the data tree
95 // // when this macro is run in an interactive ROOT/CINT session.
96 // ofile->Write();
97 //
98 //--- Author: Nick van Eijndhoven 11-mar-2005 Utrecht University
99 //- Modified: NvE $Date$ Utrecht University
100 ///////////////////////////////////////////////////////////////////////////
101  
102 #include "IceF2k.h"
103 #include "Riostream.h"
104
105 ClassImp(IceF2k) // Class implementation to enable ROOT I/O
106
107 IceF2k::IceF2k(const char* name,const char* title) : AliJob(name,title)
108 {
109 // Default constructor.
110 // By default maxevent=-1, split=99, bsize=32000, printfreq=1.
111
112  fSplit=99;
113  fBsize=32000;
114  fMaxevt=-1;
115  fPrintfreq=1;
116  fInfile="";
117  fOutfile=0;
118
119  fPdg=0;
120  fOmdb=0;
121  fFitdefs=0;
122 }
123 ///////////////////////////////////////////////////////////////////////////
124 IceF2k::~IceF2k()
125 {
126 // Default destructor.
127
128  if (fPdg)
129  {
130   delete fPdg;
131   fPdg=0;
132  }
133
134  if (fOmdb)
135  {
136   delete fOmdb;
137   fOmdb=0;
138  }
139
140  if (fFitdefs)
141  {
142   delete fFitdefs;
143   fFitdefs=0;
144  }
145 }
146 ///////////////////////////////////////////////////////////////////////////
147 void IceF2k::SetMaxEvents(Int_t n)
148 {
149 // Set the maximum number of events to be processed.
150 // n=-1 implies processing of the complete input file, which is the default
151 // initialisation in the constructor.
152  fMaxevt=n;
153 }
154 ///////////////////////////////////////////////////////////////////////////
155 void IceF2k::SetPrintFreq(Int_t f)
156 {
157 // Set the printfrequency to produce info every f events.
158 // f=1 is the default initialisation in the constructor.
159  if (f>0) fPrintfreq=f;
160 }
161 ///////////////////////////////////////////////////////////////////////////
162 void IceF2k::SetSplitLevel(Int_t split)
163 {
164 // Set the split level for the ROOT data file.
165 // split=99 is the default initialisation in the constructor.
166  if (split>=0) fSplit=split;
167 }
168 ///////////////////////////////////////////////////////////////////////////
169 void IceF2k::SetBufferSize(Int_t bsize)
170 {
171 // Set the buffer size for the ROOT data file.
172 // bsize=32000 is the default initialisation in the constructor.
173  if (bsize>=0) fBsize=bsize;
174 }
175 ///////////////////////////////////////////////////////////////////////////
176 void IceF2k::SetInputFile(TString name)
177 {
178 // Set the name of the F2K input file.
179  fInfile=name;
180 }
181 ///////////////////////////////////////////////////////////////////////////
182 void IceF2k::SetOutputFile(TFile* ofile)
183 {
184 // Set the output file for the ROOT data.
185  fOutfile=ofile;
186 }
187 ///////////////////////////////////////////////////////////////////////////
188 TDatabasePDG* IceF2k::GetPDG()
189 {
190 // Provide pointer to the PDG database
191  return fPdg;
192 }
193 ///////////////////////////////////////////////////////////////////////////
194 AliObjMatrix* IceF2k::GetOMdbase()
195 {
196 // Provide pointer to the OM geometry, calib. etc... database
197  return fOmdb;
198 }
199 ///////////////////////////////////////////////////////////////////////////
200 AliDevice* IceF2k::GetFitdefs()
201 {
202 // Provide pointer to the fit definitions
203  return fFitdefs;
204 }
205 ///////////////////////////////////////////////////////////////////////////
206 void IceF2k::Exec(Option_t* opt)
207 {
208 // Job to loop over the specified number of events and convert the 
209 // F2K data into the IceEvent structure.
210 // If maxevents<0 (default) all the entries of the input file
211 // will be processed.
212 // Every "printfreq" events a short event summary will be printed.
213 // The default value is printfreq=1.
214 // The output will be written on a standard output tree named "T".
215 //
216 // Notes :
217 // -------
218 // 1) This class is derived from AliJob, allowing a task based processing.
219 //    After the conversion of an F2K event into an IceEvent structure,
220 //    the processing of all available sub-tasks (if any) is invoked.
221 //    This provides an event-by-event (sub)task processing before the
222 //    final data structures are written out.
223 // 2) The main object in this job environment is an IceEvent* pointer.
224
225  if (fInfile=="")
226  {
227   cout << " *IceF2k Exec* No data input file specified." << endl;
228   return;
229  }
230
231  // Open the input file in the default ascii format (autodetection) for reading 
232  fInput=rdmc_mcopen(fInfile.Data(),"r",RDMC_DEFAULT_ASCII_F);
233
234  if (!fInput)
235  {
236   cout << " *IceF2k Exec* No input file found with name : " << fInfile.Data() << endl;
237   return;
238  }
239
240  // Initialise the event structure 
241  rdmc_init_mevt(&fEvent);
242
243  // Read the file header information
244  rdmc_rarr(fInput,&fHeader);
245
246  TTree* otree=0;
247  if (fOutfile)
248  {
249   otree=new TTree("T","F2K Data converted to IceEvent structures");
250   otree->SetDirectory(fOutfile);
251  }
252
253  IceEvent* evt=new IceEvent();
254  evt->SetTrackCopy(1);
255  evt->SetDevCopy(1);
256
257  // Branch in the tree for the event structure
258  if (otree) otree->Branch("IceEvent","IceEvent",&evt,fBsize,fSplit); 
259
260  // Create the particle database and extend it with some F2000 specific definitions
261  if (!fPdg) fPdg=new TDatabasePDG();
262  Double_t me=fPdg->GetParticle(11)->Mass();
263  fPdg->AddParticle("brems"   ,"brems"   ,0,1,0,0,"none",10001001,0,0);
264  fPdg->AddParticle("deltae"  ,"deltae"  ,me,1,0,-3,"Lepton",10001002,0,0);
265  fPdg->AddParticle("pairprod","pairprod",0,1,0,0,"none",10001003,0,0);
266  fPdg->AddParticle("nucl_int","nucl_Int",0,1,0,0,"none",10001004,0,0);
267  fPdg->AddParticle("mu_pair" ,"mu_pair" ,0,1,0,0,"none",10001005,0,0);
268  fPdg->AddParticle("hadrons" ,"hadrons" ,0,1,0,0,"none",10001006,0,0);
269  fPdg->AddParticle("fiberlaser","fiberlaser",0,1,0,0,"none",10002100,0,0);
270  fPdg->AddParticle("n2laser"   ,"n2laser"   ,0,1,0,0,"none",10002101,0,0);
271  fPdg->AddParticle("yaglaser"  ,"yaglaser"  ,0,1,0,0,"none",10002201,0,0);
272  fPdg->AddParticle("z_primary","z_primary",0,1,0,0,"none",10003000,0,0);
273  fPdg->AddParticle("a_primary","a_primary",0,1,0,0,"none",10003500,0,0);
274
275  // Fill the database with geometry, calib. etc... parameters
276  // for all the devices
277  FillOMdbase();
278
279  // Set the fit definitions according to the F2000 header info
280  SetFitdefs();
281
282  // Initialise the job working environment
283  SetMainObject(evt);
284  if (fOutfile)
285  {
286   AddObject(fOutfile);
287   AddObject(otree);
288  }
289
290  cout << " ***" << endl;
291  cout << " *** Start processing of job " << GetName() << " ***" << endl;
292  cout << " ***" << endl;
293  cout << " F2K input file : " << fInfile.Data() << endl;
294  cout << " Maximum number of events to be processed : " << fMaxevt << endl;
295  cout << " Print frequency : " << fPrintfreq << endl;
296  if (fOutfile)
297  {
298   cout << " ROOT output file : " << fOutfile->GetName() << endl;
299   cout << " Output characteristics : splitlevel = " << fSplit << " buffersize = " << fBsize << endl;
300  }
301
302  ListEnvironment();
303  
304  Int_t nevt=0;
305  while (!rdmc_revt(fInput,&fHeader,&fEvent))
306  {
307   if (fMaxevt>-1 && nevt>=fMaxevt) break;
308
309   // Reset the complete Event structure
310   evt->Reset();
311
312   evt->SetRunNumber(fEvent.nrun);
313   evt->SetEventNumber(fEvent.enr);
314   evt->SetMJD(fEvent.mjd,fEvent.secs,fEvent.nsecs);
315
316   PutMcTracks();
317
318   PutRecoTracks();
319
320   PutHits();
321
322   // Invoke all available sub-tasks (if any)
323   ExecuteTasks(opt);
324
325   if (!(nevt%fPrintfreq)) evt->HeaderData();
326
327   // Write the complete structure to the output Tree
328   if (otree) otree->Fill();
329
330   // Update event counter
331   nevt++;
332  }
333
334  // Remove the IceEvent object from the environment
335  // and delete it as well
336  if (evt)
337  {
338   RemoveObject(evt);
339   delete evt;
340  }
341 }
342 ///////////////////////////////////////////////////////////////////////////
343 void IceF2k::FillOMdbase()
344 {
345 // Fill the database with geometry, calib. etc... parameters 
346 // for all the devices.
347
348  if (fHeader.nch<=0) return;
349
350  if (fOmdb)
351  {
352   fOmdb->Reset();
353  }
354  else
355  {
356   fOmdb=new AliObjMatrix();
357   fOmdb->SetNameTitle("OMDBASE","The OM geometry, calib. etc... database");
358   fOmdb->SetOwner();
359  }
360
361  IceAOM* dev=0;
362  Double_t pos[3]={0,0,0};
363  for (Int_t i=0; i<fHeader.nch; i++)
364  {
365   dev=new IceAOM();
366   dev->SetUniqueID(i+1);
367   dev->SetSlotName("TYPE",1);
368   dev->SetSlotName("ORIENT",2);
369   dev->SetSlotName("T0",3);
370   dev->SetSlotName("ALPHA",4);
371   dev->SetSlotName("KADC",5);
372   dev->SetSlotName("KTOT",6);
373   dev->SetSlotName("KTDC",7);
374
375   pos[0]=fHeader.x[i];
376   pos[1]=fHeader.y[i];
377   pos[2]=fHeader.z[i];
378   dev->SetPosition(pos,"car");
379   dev->SetSignal(fHeader.type[i],1);
380   dev->SetSignal((Float_t)fHeader.costh[i],2);
381   dev->SetSignal(fHeader.cal[i].t_0,3);
382   dev->SetSignal(fHeader.cal[i].alpha_t,4);
383   dev->SetSignal(fHeader.cal[i].beta_a,5);
384   dev->SetSignal(fHeader.cal[i].beta_tot,6);
385   dev->SetSignal(fHeader.cal[i].beta_t,7);
386   fOmdb->EnterObject(i+1,1,dev);
387  }
388 }
389 ///////////////////////////////////////////////////////////////////////////
390 void IceF2k::SetFitdefs()
391 {
392 // Obtain the names of the variables for each fit procedure from the
393 // f2000 header. Each different fit procedure is then stored as a separate
394 // hit of an AliDevice object and the various fit variables are stored
395 // as separate signal slots of the corresponding hit.
396 // Via the GetFitdefs() memberfunction this AliDevice object can be
397 // retrieved and stored in the ROOT output file if wanted.
398 // The name of the object is FitDefinitions and the stored data can be
399 // inspected via the AliDevice::Data() memberfunction and looks as follows :
400 //
401 //  *AliDevice::Data* Id :0 Name : FitDefinitions
402 //    Position Vector in car coordinates : 0 0 0
403 //    Err. in car coordinates : 0 0 0
404 //  The following 8 hits are registered :
405 //  *AliSignal::Data* Id :0
406 //    Position Vector in car coordinates : 0 0 0
407 //    Err. in car coordinates : 0 0 0
408 //    Owned by device : AliDevice Name : FitDefinitions
409 //    Slot : 1 Signal value : 1 name : id
410 //    Slot : 2 Signal value : 2 name : rchi2
411 //    Slot : 3 Signal value : 3 name : prob
412 //    Slot : 4 Signal value : 4 name : sigth
413 //    Slot : 5 Signal value : 5 name : covmin
414 //    Slot : 6 Signal value : 6 name : covmax
415 //    Slot : 7 Signal value : 7 name : cutflag
416 //    Slot : 8 Signal value : 8 name : chi2
417 //  *AliSignal::Data* Id :1
418 //    Position Vector in car coordinates : 0 0 0
419 //    Err. in car coordinates : 0 0 0
420 //    Owned by device : AliDevice Name : FitDefinitions
421 //    Slot : 1 Signal value : 1 name : id
422 //    Slot : 2 Signal value : 2 name : rchi2
423 //    Slot : 3 Signal value : 3 name : prob
424 // etc....  
425 //
426 // This memberfunction is based on the original idea/code by Adam Bouchta.
427
428  if (fHeader.n_fit<=0) return;
429
430  if (fFitdefs)
431  {
432   fFitdefs->Reset(1);
433  }
434  else
435  {
436   fFitdefs=new AliDevice();
437  }
438
439  fFitdefs->SetName("FitDefinitions");
440  fFitdefs->SetHitCopy (1);
441
442  AliSignal s;
443  s.Reset();
444
445  for (Int_t i=0; i<fHeader.n_fit; i++)
446  {
447   s.SetUniqueID(fHeader.def_fit[i].id);
448
449   for (Int_t j=0; j<fHeader.def_fit[i].nwords; j++)
450   {
451    s.SetSlotName(TString(fHeader.def_fit[i].words[j]),j+1);
452    s.SetSignal(j+1,j+1);
453   }
454
455   fFitdefs->AddHit(s);
456   s.Reset(1);
457  }
458 }
459 ///////////////////////////////////////////////////////////////////////////
460 void IceF2k::PutMcTracks()
461 {
462 // Get the MC tracks from the F2000 file into the IcePack structure.
463 // Note : MC tracks are given negative track id's in the event structure.
464 // This memberfunction is based on the original code by Adam Bouchta.
465
466  IceEvent* evt=(IceEvent*)GetMainObject();
467  if (!evt || fEvent.ntrack<=0) return;
468
469  // Loop over all the tracks and add them to the current event
470  AliTrack t;
471  Double_t vec[3];
472  AliPosition r;
473  Ali3Vector p;
474  Int_t tid=0;
475  Int_t idpdg=0;
476  Int_t idf2k=0;
477  for (Int_t i=0; i<fEvent.ntrack; i++)
478  {
479   t.Reset ();
480
481   // Beginpoint of the track
482   vec[0]=fEvent.gen[i].x;
483   vec[1]=fEvent.gen[i].y;
484   vec[2]=fEvent.gen[i].z;
485   r.SetPosition(vec,"car");
486   t.SetBeginPoint(r);
487
488   // Endpoint of the track
489   vec[0]+=fEvent.gen[i].length*fEvent.gen[i].px;
490   vec[1]+=fEvent.gen[i].length*fEvent.gen[i].py;
491   vec[2]+=fEvent.gen[i].length*fEvent.gen[i].pz;
492   r.SetPosition(vec,"car");
493   t.SetEndPoint(r);
494
495   // Momentum in GeV/c
496   vec[0]=fEvent.gen[i].e*fEvent.gen[i].px*1e-3;
497   vec[1]=fEvent.gen[i].e*fEvent.gen[i].py*1e-3;
498   vec[2]=fEvent.gen[i].e*fEvent.gen[i].pz*1e-3;
499   p.SetVector (vec,"car");
500   t.Set3Momentum(p);
501
502   // MC tracks are indicated by negative track id's
503   tid=fEvent.gen[i].tag;
504   t.SetId(-abs(tid));
505
506   idf2k=fEvent.gen[i].id;
507   idpdg=0;
508   if (idf2k>1000)
509   {
510    idpdg=idf2k+10000000;
511   }
512   else if (idf2k <= 48)
513   {
514    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
515   }
516   else
517   {
518    if (idf2k==201) idpdg=12;
519    if (idf2k==202) idpdg=14;
520    if (idf2k==203) idpdg=16;
521    if (idf2k==204) idpdg=-12;
522    if (idf2k==205) idpdg=-14;
523    if (idf2k==206) idpdg=-16;
524   }
525
526   t.SetParticleCode(idpdg);
527   t.SetName(fPdg->GetParticle(idpdg)->GetName());
528   t.SetTitle("MC track");
529   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
530   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
531
532   evt->AddTrack(t);
533  }
534
535  // Create the pointers to each particle's parent particle.
536  Int_t txid=0;
537  Int_t parid=0;
538  for (Int_t itk=1; itk<=evt->GetNtracks (); itk++)
539  {
540   AliTrack* tx=evt->GetTrack(itk);
541
542   if (!tx) continue;
543
544   txid=tx->GetId();
545
546   parid=-1;
547   for (Int_t j=0; j<fEvent.ntrack; j++)
548   {
549    tid=fEvent.gen[j].tag;
550    if (-abs(tid) == txid) parid=fEvent.gen[j].parent;
551   }
552
553   if (parid<0) continue;
554
555   AliTrack* tpar=evt->GetIdTrack(-abs(parid));
556
557   if (!tpar) continue;
558
559   tx->SetParentTrack(tpar);
560  }
561 }
562 ///////////////////////////////////////////////////////////////////////////
563 void IceF2k::PutRecoTracks()
564 {
565 // Get the reconstructed tracks from the F2000 file into the IcePack structure.
566 // Note : Reco tracks are given positive track id's in the event structure.
567 // This memberfunction is based on the original code by Adam Bouchta.
568
569  IceEvent* evt=(IceEvent*)GetMainObject();
570  if (!evt || fEvent.nfit<=0) return;
571
572  // Loop over all the tracks and add them to the current event
573  AliTrack t;
574  Double_t vec[3];
575  AliPosition r;
576  Ali3Vector p;
577  Int_t tid=0;
578  Int_t idpdg=0;
579  Int_t idf2k=0;
580  for (Int_t i=0; i<fEvent.nfit; i++)
581  {
582   t.Reset ();
583
584   // Beginpoint of the track
585   vec[0]=fEvent.rec[i].x;
586   vec[1]=fEvent.rec[i].y;
587   vec[2]=fEvent.rec[i].z;
588   r.SetPosition(vec,"car");
589   t.SetBeginPoint(r);
590
591   // Endpoint of the track
592   vec[0]+=fEvent.rec[i].length*fEvent.rec[i].px;
593   vec[1]+=fEvent.rec[i].length*fEvent.rec[i].py;
594   vec[2]+=fEvent.rec[i].length*fEvent.rec[i].pz;
595   r.SetPosition(vec,"car");
596   t.SetEndPoint(r);
597
598   // Momentum in GeV/c
599   if (fEvent.rec[i].e > 0)
600   {
601    vec[0]=fEvent.rec[i].e*fEvent.rec[i].px*1e-3;
602    vec[1]=fEvent.rec[i].e*fEvent.rec[i].py*1e-3;
603    vec[2]=fEvent.rec[i].e*fEvent.rec[i].pz*1e-3;
604   }
605   else // Give the track a nominal momentum of 1 GeV/c
606   {
607    vec[0]=fEvent.rec[i].px;
608    vec[1]=fEvent.rec[i].py;
609    vec[2]=fEvent.rec[i].pz;
610   }
611   p.SetVector (vec,"car");
612   t.Set3Momentum(p);
613
614   // Use the fit number as track id
615   tid=fEvent.rec[i].tag;
616   t.SetId(abs(tid));
617
618   idf2k=fEvent.rec[i].id;
619   idpdg=0;
620   if (idf2k>1000)
621   {
622    idpdg=idf2k+10000000;
623   }
624   else if (idf2k <= 48)
625   {
626    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
627   }
628   else
629   {
630    if (idf2k==201) idpdg=12;
631    if (idf2k==202) idpdg=14;
632    if (idf2k==203) idpdg=16;
633    if (idf2k==204) idpdg=-12;
634    if (idf2k==205) idpdg=-14;
635    if (idf2k==206) idpdg=-16;
636   }
637
638   t.SetParticleCode(idpdg);
639   t.SetName(fPdg->GetParticle(idpdg)->GetName());
640   t.SetTitle("RECO track");
641   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
642   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
643
644   // Retrieve the various fit parameters for this track
645   AliSignal* fitdata=fFitdefs->GetIdHit(i);
646   for (Int_t jval=0; jval<fEvent.fresult[i].nval; jval++)
647   {
648    fitdata->SetSignal(fEvent.fresult[i].val[jval],jval+1);
649   }
650
651   // Store the various fit parameters for this track
652   t.SetFitDetails(fitdata);
653
654   // Store the various reco tracks as track hypotheses.
655   // A copy of the first reco track is entered as a new track instance
656   // into the event and all reco tracks (incl. the first one) are
657   // stored as hypotheses linked to this new reco track.
658   if (i==0)
659   {
660    evt->AddTrack(t);
661    AliTrack* tx=evt->GetTrack(evt->GetNtracks());
662    Int_t nrec=evt->GetNtracks(1);
663    tx->SetId(nrec+1);
664   }
665   AliTrack* tx=evt->GetTrack(evt->GetNtracks());
666   if (tx) tx->AddTrackHypothesis(t);
667  }
668 }
669 ///////////////////////////////////////////////////////////////////////////
670 void IceF2k::PutHits()
671 {
672 // Get the hit and waveform info from the F2000 file into the IcePack structure.
673 // This memberfunction is based on the original code by Adam Bouchta.
674
675  IceEvent* evt=(IceEvent*)GetMainObject();
676  if (!evt) return;
677
678  // Loop over all the hits and add them to the current event
679  IceAOM om;
680  AliSignal s;
681  s.SetSlotName("ADC",1);
682  s.SetSlotName("LE",2);
683  s.SetSlotName("TOT",3);
684  Int_t chan=0;
685  Int_t maxchan=800;
686  if (fOmdb) maxchan=fHeader.nch;
687  IceAOM* omx=0;
688  AliSignal* sx=0;
689  Int_t tid=0;
690  AliTrack* tx=0;
691  for (Int_t i=0; i<fEvent.nhits; i++)
692  {
693   chan=fEvent.h[i].ch+1;
694   if (chan>maxchan) continue; // Channels 9001, 9002 etc are trigger channels
695
696   // Get corresponding device from the current event structure  
697   omx=(IceAOM*)evt->GetIdDevice(chan);
698   if (!omx)
699   {
700    if (fOmdb)
701    {
702     omx=(IceAOM*)fOmdb->GetObject(chan,1);
703     evt->AddDevice(omx);
704    }
705    else
706    {
707     om.Reset(1);
708     om.SetUniqueID(chan);
709     evt->AddDevice(om);
710    }
711    omx=(IceAOM*)evt->GetIdDevice(chan);
712   }
713
714   if (!omx) continue;
715
716   s.Reset();
717   s.SetUniqueID(fEvent.h[i].id);
718   s.SetSignal(fEvent.h[i].amp,1);
719   s.SetSignal(fEvent.h[i].t,2);
720   s.SetSignal(fEvent.h[i].tot,3);
721
722   omx->AddHit(s);
723
724   sx=omx->GetHit(omx->GetNhits());
725   if (!sx) continue;
726
727   // Bi-directional link between this hit and the track that caused the ADC value.
728   // This F2K info is probably only present for MC tracks.
729   tid=fEvent.h[i].ma;
730   if (tid > 0)
731   {
732    tx=evt->GetIdTrack(tid); // Reco tracks
733    if (!tx) tx=evt->GetIdTrack(-tid); // MC tracks
734    if (tx) sx->AddLink(tx);
735   }
736   else
737   {
738    if (tid == -2) s.SetNameTitle("N","Noise");
739    if (tid == -3) s.SetNameTitle("A","Afterpulse");
740   }
741  }
742
743  // Loop over all the waveforms and add the histo(s) to the corresponding OM's
744  TH1F histo;
745  Int_t nbins=0;
746  Float_t xlow=0;
747  Float_t xup=0;
748  TString hname;
749  for (Int_t iwf=0; iwf<fEvent.nwf; iwf++)
750  {
751   chan=fEvent.wf[iwf].om;
752   if (chan<=0 || chan>maxchan) continue; // Skip trigger channels
753
754   // Get corresponding device from the current event structure  
755   omx=(IceAOM*)evt->GetIdDevice(chan);
756   if (!omx)
757   {
758    if (fOmdb)
759    {
760     omx=(IceAOM*)fOmdb->GetObject(chan,1);
761     evt->AddDevice(omx);
762    }
763    else
764    {
765     om.Reset(1);
766     om.SetUniqueID(chan);
767     evt->AddDevice(om);
768    }
769    omx=(IceAOM*)evt->GetIdDevice(chan);
770   }
771
772   if (!omx) continue;
773
774   omx->SetSlotName("BASELINE",omx->GetNnames()+1);
775   omx->SetSignal(fEvent.wf[iwf].baseline,"BASELINE");
776
777   // Fill the waveform histogram
778   hname="OM";
779   hname+=chan;
780   hname+="-WF";
781   hname+=omx->GetNwaveforms()+1;
782
783   histo.Reset();
784   histo.SetName(hname.Data());
785   nbins=fEvent.wf[iwf].ndigi;
786   xlow=fEvent.wf[iwf].t_start;
787   xup=xlow+float(nbins)*fEvent.wf[iwf].t_bin;
788   histo.SetBins(nbins,xlow,xup);
789
790   for (Int_t jbin=1; jbin<=fEvent.wf[iwf].ndigi; jbin++)
791   {
792    histo.SetBinContent(jbin,fEvent.wf[iwf].digi[jbin]);
793   }
794
795   omx->SetWaveform(&histo,omx->GetNwaveforms()+1);
796  }
797
798  // Set bi-directional links between hits and reco track hypotheses.
799  // Note : Reco tracks are recognised by their positive id.
800  Int_t hid=0;
801  TObjArray* rectracks=evt->GetTracks(1);
802  for (Int_t jtk=0; jtk<rectracks->GetEntries(); jtk++)
803  {
804   tx=(AliTrack*)rectracks->At(jtk);
805   if (!tx) continue;
806   
807   for (Int_t jhyp=1; jhyp<=tx->GetNhypotheses(); jhyp++)
808   {
809    AliTrack* hypx=tx->GetTrackHypothesis(jhyp);
810    if (!hypx) continue;
811
812    // Loop over all combinations of F2K fits and used OM hits
813    for (Int_t k=0; k<fEvent.nfit_uses; k++)
814    {
815     if (fEvent.fit_uses[k].useid != hypx->GetId()) continue;
816     hid=fEvent.fit_uses[k].hitid;
817     sx=evt->GetIdHit(hid,"IceAOM");
818     if (sx) sx->AddLink(hypx);
819    }
820   }
821  }
822 }
823 ///////////////////////////////////////////////////////////////////////////