]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/icepack/iceconvert/IceF2k.cxx
05-aug-2005 NvE Typos fixed in AliAttrib::GetCalFunction and AliAttrib::GetDecalFunct...
[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 (fPrintfreq)
326   {
327    if (!(nevt%fPrintfreq)) evt->HeaderData();
328   }
329
330   // Write the complete structure to the output Tree
331   if (otree) otree->Fill();
332
333   // Update event counter
334   nevt++;
335  }
336
337  // Remove the IceEvent object from the environment
338  // and delete it as well
339  if (evt)
340  {
341   RemoveObject(evt);
342   delete evt;
343  }
344 }
345 ///////////////////////////////////////////////////////////////////////////
346 void IceF2k::FillOMdbase()
347 {
348 // Fill the database with geometry, calib. etc... parameters 
349 // for all the devices.
350
351  if (fHeader.nch<=0) return;
352
353  Int_t geocal=fHeader.is_calib.geo;
354  Int_t adccal=fHeader.is_calib.adc;
355  Int_t tdccal=fHeader.is_calib.tdc;
356  Int_t totcal=fHeader.is_calib.tot;
357  Int_t utccal=fHeader.is_calib.utc;
358
359  TF1 fadccal("fadccal","(x-[1])*[0]");
360  TF1 fadcdecal("fadcdecal","(x/[0])+[1]");
361  TF1 ftdccal("ftdccal","(x*[0])-[1]-([0]-1.)*32767.-[2]/sqrt([3])");
362  TF1 ftdcdecal("ftdcdecal","(x+([0]-1.)*32767.+[1]+[2]/sqrt([3]))/[0]");
363  TF1 ftotcal("ftotcal","x*[0]");
364  TF1 ftotdecal("ftotdecal","x/[0]");
365
366  if (fOmdb)
367  {
368   fOmdb->Reset();
369  }
370  else
371  {
372   fOmdb=new AliObjMatrix();
373   fOmdb->SetNameTitle("OMDBASE","The OM geometry, calib. etc... database");
374   fOmdb->SetOwner();
375  }
376
377  IceAOM* dev=0;
378  Double_t pos[3]={0,0,0};
379  for (Int_t i=0; i<fHeader.nch; i++)
380  {
381   dev=new IceAOM();
382   dev->SetUniqueID(i+1);
383
384   dev->SetSlotName("ADC",1);
385   dev->SetSlotName("LE",2);
386   dev->SetSlotName("TOT",3);
387
388   dev->SetSlotName("TYPE",4);
389   dev->SetSlotName("ORIENT",5);
390   dev->SetSlotName("THRESH",6);
391   dev->SetSlotName("SENSIT",7);
392   dev->SetSlotName("BETA-TDC",8);
393   dev->SetSlotName("T0",9);
394   dev->SetSlotName("ALPHA-TDC",10);
395   dev->SetSlotName("PED-ADC",11);
396   dev->SetSlotName("BETA-ADC",12);
397   dev->SetSlotName("KAPPA-ADC",13);
398   dev->SetSlotName("PED-TOT",14);
399   dev->SetSlotName("BETA-TOT",15);
400   dev->SetSlotName("KAPPA-TOT",16);
401
402   pos[0]=fHeader.x[i];
403   pos[1]=fHeader.y[i];
404   pos[2]=fHeader.z[i];
405   dev->SetPosition(pos,"car");
406
407   fadccal.SetParameter(0,fHeader.cal[i].beta_a);
408   fadccal.SetParameter(1,fHeader.cal[i].ped);
409   fadcdecal.SetParameter(0,fHeader.cal[i].beta_a);
410   if (!fHeader.cal[i].beta_a) fadcdecal.SetParameter(0,1);
411   fadcdecal.SetParameter(1,fHeader.cal[i].ped);
412
413   ftdccal.SetParameter(0,fHeader.cal[i].beta_t);
414   ftdccal.SetParameter(1,fHeader.cal[i].t_0);
415   ftdccal.SetParameter(2,fHeader.cal[i].alpha_t);
416   ftdccal.SetParameter(3,1.e20);
417   ftdcdecal.SetParameter(0,fHeader.cal[i].beta_t);
418   if (!fHeader.cal[i].beta_t) ftdcdecal.SetParameter(0,1);
419   ftdcdecal.SetParameter(1,fHeader.cal[i].t_0);
420   ftdcdecal.SetParameter(2,fHeader.cal[i].alpha_t);
421   ftdcdecal.SetParameter(3,1.e20);
422
423   ftotcal.SetParameter(0,fHeader.cal[i].beta_tot);
424   ftotdecal.SetParameter(0,fHeader.cal[i].beta_tot);
425   if (!fHeader.cal[i].beta_tot) ftotdecal.SetParameter(0,1);
426
427   if (adccal)
428   {
429    dev->SetDecalFunction(&fadcdecal,1);
430   }
431   else
432   {
433    dev->SetCalFunction(&fadccal,1);
434   }
435
436   if (tdccal)
437   {
438    dev->SetDecalFunction(&ftdcdecal,2);
439   }
440   else
441   {
442    dev->SetCalFunction(&ftdccal,2);
443   }
444
445   if (totcal)
446   {
447    dev->SetDecalFunction(&ftotdecal,3);
448   }
449   else
450   {
451    dev->SetCalFunction(&ftotcal,3);
452   }
453
454   dev->SetSignal(fHeader.type[i],4);
455   dev->SetSignal((Float_t)fHeader.costh[i],5);
456   dev->SetSignal(fHeader.thresh[i],6);
457   dev->SetSignal(fHeader.sensit[i],7);
458   dev->SetSignal(fHeader.cal[i].beta_t,8);
459   dev->SetSignal(fHeader.cal[i].t_0,9);
460   dev->SetSignal(fHeader.cal[i].alpha_t,10);
461   dev->SetSignal(fHeader.cal[i].ped,11);
462   dev->SetSignal(fHeader.cal[i].beta_a,12);
463   dev->SetSignal(fHeader.cal[i].kappa,13);
464   dev->SetSignal(fHeader.cal[i].ped_tot,14);
465   dev->SetSignal(fHeader.cal[i].beta_tot,15);
466   dev->SetSignal(fHeader.cal[i].kappa_tot,16);
467
468   fOmdb->EnterObject(i+1,1,dev);
469  }
470 }
471 ///////////////////////////////////////////////////////////////////////////
472 void IceF2k::SetFitdefs()
473 {
474 // Obtain the names of the variables for each fit procedure from the
475 // f2000 header. Each different fit procedure is then stored as a separate
476 // hit of an AliDevice object and the various fit variables are stored
477 // as separate signal slots of the corresponding hit.
478 // Via the GetFitdefs() memberfunction this AliDevice object can be
479 // retrieved and stored in the ROOT output file if wanted.
480 // The name of the object is FitDefinitions and the stored data can be
481 // inspected via the AliDevice::Data() memberfunction and looks as follows :
482 //
483 //  *AliDevice::Data* Id :0 Name : FitDefinitions
484 //    Position Vector in car coordinates : 0 0 0
485 //    Err. in car coordinates : 0 0 0
486 //  The following 8 hits are registered :
487 //  *AliSignal::Data* Id :0
488 //    Position Vector in car coordinates : 0 0 0
489 //    Err. in car coordinates : 0 0 0
490 //    Owned by device : AliDevice Name : FitDefinitions
491 //    Slot : 1 Signal value : 1 name : id
492 //    Slot : 2 Signal value : 2 name : rchi2
493 //    Slot : 3 Signal value : 3 name : prob
494 //    Slot : 4 Signal value : 4 name : sigth
495 //    Slot : 5 Signal value : 5 name : covmin
496 //    Slot : 6 Signal value : 6 name : covmax
497 //    Slot : 7 Signal value : 7 name : cutflag
498 //    Slot : 8 Signal value : 8 name : chi2
499 //  *AliSignal::Data* Id :1
500 //    Position Vector in car coordinates : 0 0 0
501 //    Err. in car coordinates : 0 0 0
502 //    Owned by device : AliDevice Name : FitDefinitions
503 //    Slot : 1 Signal value : 1 name : id
504 //    Slot : 2 Signal value : 2 name : rchi2
505 //    Slot : 3 Signal value : 3 name : prob
506 // etc....  
507 //
508 // This memberfunction is based on the original idea/code by Adam Bouchta.
509
510  if (fHeader.n_fit<=0) return;
511
512  if (fFitdefs)
513  {
514   fFitdefs->Reset(1);
515  }
516  else
517  {
518   fFitdefs=new AliDevice();
519  }
520
521  fFitdefs->SetName("FitDefinitions");
522  fFitdefs->SetHitCopy (1);
523
524  AliSignal s;
525  s.Reset();
526
527  for (Int_t i=0; i<fHeader.n_fit; i++)
528  {
529   s.SetUniqueID(fHeader.def_fit[i].id);
530
531   for (Int_t j=0; j<fHeader.def_fit[i].nwords; j++)
532   {
533    s.SetSlotName(TString(fHeader.def_fit[i].words[j]),j+1);
534    s.SetSignal(j+1,j+1);
535   }
536
537   fFitdefs->AddHit(s);
538   s.Reset(1);
539  }
540 }
541 ///////////////////////////////////////////////////////////////////////////
542 void IceF2k::PutMcTracks()
543 {
544 // Get the MC tracks from the F2000 file into the IcePack structure.
545 // Note : MC tracks are given negative track id's in the event structure.
546 // This memberfunction is based on the original code by Adam Bouchta.
547
548  IceEvent* evt=(IceEvent*)GetMainObject();
549  if (!evt || fEvent.ntrack<=0) return;
550
551  // Loop over all the tracks and add them to the current event
552  AliTrack t;
553  Double_t vec[3];
554  AliPosition r;
555  Ali3Vector p;
556  Int_t tid=0;
557  Int_t idpdg=0;
558  Int_t idf2k=0;
559  for (Int_t i=0; i<fEvent.ntrack; i++)
560  {
561   t.Reset ();
562
563   // Beginpoint of the track
564   vec[0]=fEvent.gen[i].x;
565   vec[1]=fEvent.gen[i].y;
566   vec[2]=fEvent.gen[i].z;
567   r.SetPosition(vec,"car");
568   t.SetBeginPoint(r);
569
570   // Endpoint of the track
571   vec[0]+=fEvent.gen[i].length*fEvent.gen[i].px;
572   vec[1]+=fEvent.gen[i].length*fEvent.gen[i].py;
573   vec[2]+=fEvent.gen[i].length*fEvent.gen[i].pz;
574   r.SetPosition(vec,"car");
575   t.SetEndPoint(r);
576
577   // Momentum in GeV/c
578   vec[0]=fEvent.gen[i].e*fEvent.gen[i].px*1e-3;
579   vec[1]=fEvent.gen[i].e*fEvent.gen[i].py*1e-3;
580   vec[2]=fEvent.gen[i].e*fEvent.gen[i].pz*1e-3;
581   p.SetVector (vec,"car");
582   t.Set3Momentum(p);
583
584   // MC tracks are indicated by negative track id's
585   tid=fEvent.gen[i].tag;
586   t.SetId(-abs(tid));
587
588   idf2k=fEvent.gen[i].id;
589   idpdg=0;
590   if (idf2k>1000)
591   {
592    idpdg=idf2k+10000000;
593   }
594   else if (idf2k <= 48)
595   {
596    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
597   }
598   else
599   {
600    if (idf2k==201) idpdg=12;
601    if (idf2k==202) idpdg=14;
602    if (idf2k==203) idpdg=16;
603    if (idf2k==204) idpdg=-12;
604    if (idf2k==205) idpdg=-14;
605    if (idf2k==206) idpdg=-16;
606   }
607
608   t.SetParticleCode(idpdg);
609   t.SetName(fPdg->GetParticle(idpdg)->GetName());
610   t.SetTitle("MC track");
611   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
612   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
613
614   evt->AddTrack(t);
615  }
616
617  // Create the pointers to each particle's parent particle.
618  Int_t txid=0;
619  Int_t parid=0;
620  for (Int_t itk=1; itk<=evt->GetNtracks (); itk++)
621  {
622   AliTrack* tx=evt->GetTrack(itk);
623
624   if (!tx) continue;
625
626   txid=tx->GetId();
627
628   parid=-1;
629   for (Int_t j=0; j<fEvent.ntrack; j++)
630   {
631    tid=fEvent.gen[j].tag;
632    if (-abs(tid) == txid) parid=fEvent.gen[j].parent;
633   }
634
635   if (parid<0) continue;
636
637   AliTrack* tpar=evt->GetIdTrack(-abs(parid));
638
639   if (!tpar) continue;
640
641   tx->SetParentTrack(tpar);
642  }
643 }
644 ///////////////////////////////////////////////////////////////////////////
645 void IceF2k::PutRecoTracks()
646 {
647 // Get the reconstructed tracks from the F2000 file into the IcePack structure.
648 // Note : Reco tracks are given positive track id's in the event structure.
649 // This memberfunction is based on the original code by Adam Bouchta.
650
651  IceEvent* evt=(IceEvent*)GetMainObject();
652  if (!evt || fEvent.nfit<=0) return;
653
654  // Loop over all the tracks and add them to the current event
655  AliTrack t;
656  Double_t vec[3];
657  AliPosition r;
658  Ali3Vector p;
659  Int_t tid=0;
660  Int_t idpdg=0;
661  Int_t idf2k=0;
662  for (Int_t i=0; i<fEvent.nfit; i++)
663  {
664   t.Reset ();
665
666   // Beginpoint of the track
667   vec[0]=fEvent.rec[i].x;
668   vec[1]=fEvent.rec[i].y;
669   vec[2]=fEvent.rec[i].z;
670   r.SetPosition(vec,"car");
671   t.SetBeginPoint(r);
672
673   // Endpoint of the track
674   vec[0]+=fEvent.rec[i].length*fEvent.rec[i].px;
675   vec[1]+=fEvent.rec[i].length*fEvent.rec[i].py;
676   vec[2]+=fEvent.rec[i].length*fEvent.rec[i].pz;
677   r.SetPosition(vec,"car");
678   t.SetEndPoint(r);
679
680   // Momentum in GeV/c
681   if (fEvent.rec[i].e > 0)
682   {
683    vec[0]=fEvent.rec[i].e*fEvent.rec[i].px*1e-3;
684    vec[1]=fEvent.rec[i].e*fEvent.rec[i].py*1e-3;
685    vec[2]=fEvent.rec[i].e*fEvent.rec[i].pz*1e-3;
686   }
687   else // Give the track a nominal momentum of 1 GeV/c
688   {
689    vec[0]=fEvent.rec[i].px;
690    vec[1]=fEvent.rec[i].py;
691    vec[2]=fEvent.rec[i].pz;
692   }
693   p.SetVector (vec,"car");
694   t.Set3Momentum(p);
695
696   // Use the fit number as track id
697   tid=fEvent.rec[i].tag;
698   t.SetId(abs(tid));
699
700   idf2k=fEvent.rec[i].id;
701   idpdg=0;
702   if (idf2k>1000)
703   {
704    idpdg=idf2k+10000000;
705   }
706   else if (idf2k <= 48)
707   {
708    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
709   }
710   else
711   {
712    if (idf2k==201) idpdg=12;
713    if (idf2k==202) idpdg=14;
714    if (idf2k==203) idpdg=16;
715    if (idf2k==204) idpdg=-12;
716    if (idf2k==205) idpdg=-14;
717    if (idf2k==206) idpdg=-16;
718   }
719
720   t.SetParticleCode(idpdg);
721   t.SetName(fPdg->GetParticle(idpdg)->GetName());
722   t.SetTitle("RECO track");
723   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
724   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
725
726   // Retrieve the various fit parameters for this track
727   AliSignal* fitdata=fFitdefs->GetIdHit(i);
728   for (Int_t jval=0; jval<fEvent.fresult[i].nval; jval++)
729   {
730    fitdata->SetSignal(fEvent.fresult[i].val[jval],jval+1);
731   }
732
733   // Store the various fit parameters for this track
734   t.SetFitDetails(fitdata);
735
736   // Store the various reco tracks as track hypotheses.
737   // A copy of the first reco track is entered as a new track instance
738   // into the event and all reco tracks (incl. the first one) are
739   // stored as hypotheses linked to this new reco track.
740   if (i==0)
741   {
742    evt->AddTrack(t);
743    AliTrack* tx=evt->GetTrack(evt->GetNtracks());
744    Int_t nrec=evt->GetNtracks(1);
745    tx->SetId(nrec+1);
746   }
747   AliTrack* tx=evt->GetTrack(evt->GetNtracks());
748   if (tx) tx->AddTrackHypothesis(t);
749  }
750 }
751 ///////////////////////////////////////////////////////////////////////////
752 void IceF2k::PutHits()
753 {
754 // Get the hit and waveform info from the F2000 file into the IcePack structure.
755 // This memberfunction is based on the original code by Adam Bouchta.
756
757  IceEvent* evt=(IceEvent*)GetMainObject();
758  if (!evt) return;
759
760  // Loop over all the hits and add them to the current event
761  IceAOM om;
762  AliSignal s;
763  s.SetSlotName("ADC",1);
764  s.SetSlotName("LE",2);
765  s.SetSlotName("TOT",3);
766  Int_t chan=0;
767  Int_t maxchan=800;
768  if (fOmdb) maxchan=fHeader.nch;
769  IceAOM* omx=0;
770  AliSignal* sx=0;
771  Int_t tid=0;
772  AliTrack* tx=0;
773  Float_t adc=0;
774  for (Int_t i=0; i<fEvent.nhits; i++)
775  {
776   chan=fEvent.h[i].ch+1;
777   if (chan>maxchan) continue; // Channels 9001, 9002 etc are trigger channels
778
779   // Get corresponding device from the current event structure  
780   omx=(IceAOM*)evt->GetIdDevice(chan);
781   if (!omx)
782   {
783    if (fOmdb)
784    {
785     omx=(IceAOM*)fOmdb->GetObject(chan,1);
786     evt->AddDevice(omx);
787    }
788    else
789    {
790     om.Reset(1);
791     om.SetUniqueID(chan);
792     evt->AddDevice(om);
793    }
794    omx=(IceAOM*)evt->GetIdDevice(chan);
795   }
796
797   if (!omx) continue;
798
799   s.Reset();
800   s.SetUniqueID(fEvent.h[i].id);
801   s.SetSignal(fEvent.h[i].amp,1);
802   s.SetSignal(fEvent.h[i].t,2);
803   s.SetSignal(fEvent.h[i].tot,3);
804
805   omx->AddHit(s);
806
807   sx=omx->GetHit(omx->GetNhits());
808   if (!sx) continue;
809
810   // ADC dependent TDC (de)calibration function for this hit
811   TF1* fcal=omx->GetCalFunction("LE");
812   TF1* fdecal=omx->GetDecalFunction("LE");
813   if (fcal) sx->SetCalFunction(fcal,2);
814   if (fdecal) sx->SetDecalFunction(fdecal,2);
815   fcal=sx->GetCalFunction(2);
816   fdecal=sx->GetDecalFunction(2);
817   adc=sx->GetSignal(1,-4);
818   if (adc>0)
819   {
820    if (fcal) fcal->SetParameter(3,adc);
821    if (fdecal) fdecal->SetParameter(3,adc);
822   }
823   else
824   {
825    if (fcal) fcal->SetParameter(3,0);
826    if (fdecal) fdecal->SetParameter(3,0);
827   }
828
829   // Bi-directional link between this hit and the track that caused the ADC value.
830   // This F2K info is probably only present for MC tracks.
831   tid=fEvent.h[i].ma;
832   if (tid > 0)
833   {
834    tx=evt->GetIdTrack(tid); // Reco tracks
835    if (!tx) tx=evt->GetIdTrack(-tid); // MC tracks
836    if (tx) sx->AddLink(tx);
837   }
838   else
839   {
840    if (tid == -2) s.SetNameTitle("N","Noise");
841    if (tid == -3) s.SetNameTitle("A","Afterpulse");
842   }
843  }
844
845  // Loop over all the waveforms and add the histo(s) to the corresponding OM's
846  TH1F histo;
847  Int_t nbins=0;
848  Float_t xlow=0;
849  Float_t xup=0;
850  TString hname;
851  for (Int_t iwf=0; iwf<fEvent.nwf; iwf++)
852  {
853   chan=fEvent.wf[iwf].om;
854   if (chan<=0 || chan>maxchan) continue; // Skip trigger channels
855
856   // Get corresponding device from the current event structure  
857   omx=(IceAOM*)evt->GetIdDevice(chan);
858   if (!omx)
859   {
860    if (fOmdb)
861    {
862     omx=(IceAOM*)fOmdb->GetObject(chan,1);
863     evt->AddDevice(omx);
864    }
865    else
866    {
867     om.Reset(1);
868     om.SetUniqueID(chan);
869     evt->AddDevice(om);
870    }
871    omx=(IceAOM*)evt->GetIdDevice(chan);
872   }
873
874   if (!omx) continue;
875
876   omx->SetSlotName("BASELINE",omx->GetNnames()+1);
877   omx->SetSignal(fEvent.wf[iwf].baseline,"BASELINE");
878
879   // Fill the waveform histogram
880   hname="OM";
881   hname+=chan;
882   hname+="-WF";
883   hname+=omx->GetNwaveforms()+1;
884
885   histo.Reset();
886   histo.SetName(hname.Data());
887   nbins=fEvent.wf[iwf].ndigi;
888   xlow=fEvent.wf[iwf].t_start;
889   xup=xlow+float(nbins)*fEvent.wf[iwf].t_bin;
890   histo.SetBins(nbins,xlow,xup);
891
892   for (Int_t jbin=1; jbin<=fEvent.wf[iwf].ndigi; jbin++)
893   {
894    histo.SetBinContent(jbin,fEvent.wf[iwf].digi[jbin-1]);
895   }
896
897   omx->SetWaveform(&histo,omx->GetNwaveforms()+1);
898  }
899
900  // Set bi-directional links between hits and reco track hypotheses.
901  // Note : Reco tracks are recognised by their positive id.
902  Int_t hid=0;
903  TObjArray* rectracks=evt->GetTracks(1);
904  for (Int_t jtk=0; jtk<rectracks->GetEntries(); jtk++)
905  {
906   tx=(AliTrack*)rectracks->At(jtk);
907   if (!tx) continue;
908   
909   for (Int_t jhyp=1; jhyp<=tx->GetNhypotheses(); jhyp++)
910   {
911    AliTrack* hypx=tx->GetTrackHypothesis(jhyp);
912    if (!hypx) continue;
913
914    // Loop over all combinations of F2K fits and used OM hits
915    for (Int_t k=0; k<fEvent.nfit_uses; k++)
916    {
917     if (fEvent.fit_uses[k].useid != hypx->GetId()) continue;
918     hid=fEvent.fit_uses[k].hitid;
919     sx=evt->GetIdHit(hid,"IceAOM");
920     if (sx) sx->AddLink(hypx);
921    }
922   }
923  }
924 }
925 ///////////////////////////////////////////////////////////////////////////