]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/icepack/iceconvert/IceF2k.cxx
22-jun-2006 NvE Linear track storage facility introduced in AliSignal to enable
[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 // Note : Sometimes the filtering/reco process which produced the F2K file
33 //        may have introduced a shift (i.e. offset) in the hit times w.r.t.
34 //        the actual trigger time. The aim of this is to obtain the hit times
35 //        centered more or less around zero.
36 //        In case of real data, this is recorded in the F2K data itself and
37 //        as such will be taken automatically into account by this IceF2k
38 //        processor such that all times will be provided again unshifted.
39 //        In other words, all times will be w.r.t. the actual trigger time
40 //        as recorded in the trigger data device named "Trigger" in the IceEvent
41 //        structure.
42 //        In case of simulated data this shift is not available in the F2K data.
43 //        The offset denoted in the F2K record is related to the time of the
44 //        primary interaction to put it well ahead of the detector trigger.
45 //        This primary interaction time, however, is irrelevant for the
46 //        reconstruction of the recorded hit patterns. 
47 //        If a user had introduced a shift in producing the MC data,
48 //        very frequently (but not always) a value of -19000 is used.
49 //        For the IceF2k processing, the user can manually introduce a
50 //        time offset in case of MC data via the memberfunction SetMcToffset().
51 //        This user defined offset value will then be used to correct all
52 //        the hit times such that they will be provided again unshifted
53 //        w.r.t. the actual trigger time as recorded in the device named 
54 //        "Trigger" in the IceEvent structure.
55 //        By default the MC time offset is set to 0 in the constructor
56 //        of this class.
57 //
58 // Usage example :
59 // ---------------
60 //
61 // Note : This example creates automatically the ROOT output file, which
62 //        is the most user friendly way of running the conversion job. 
63 //        In the subdirectory /macros the example macro icef2k.cc provides
64 //        an example of how to create a ROOT output file yourself and passing
65 //        this file via a pointer to IceF2k. 
66 //
67 // gSystem->Load("ralice");
68 // gSystem->Load("icepack");
69 // gSystem->Load("iceconvert");
70 //
71 // IceF2k q("IceF2k","F2K to IcePack data structure conversion");
72 //
73 // // Limit the number of entries for testing
74 // q.SetMaxEvents(10);
75 //
76 // // Print frequency to produce a short summary print every printfreq events
77 // q.SetPrintFreq(1);
78 //
79 // // Split level for the output structures
80 // q.SetSplitLevel(0);
81 //
82 // // Buffer size for the output structures
83 // q.SetBufferSize(32000);
84 //
85 // // The F2K input filename(s)
86 // q.AddInputFile("run7825.f2k");
87 //
88 // // Output file for the event structures
89 // q.SetOutputFile("events.root");
90 //
91 // ///////////////////////////////////////////////////////////////////
92 // // Here the user can specify his/her sub-tasks to be executed
93 // // on an event-by-event basis after the IceEvent structure
94 // // has been filled and before the data is written out.
95 // // Sub-tasks (i.e. a user classes derived from TTask) are entered
96 // // as follows :
97 // //
98 // //    MyXtalk task1("task1","Cross talk correction");
99 // //    MyClean task2("task2","Hit cleaning");
100 // //    q.Add(&task1);
101 // //    q.Add(&task2);
102 // //
103 // // The sub-tasks will be executed in the order as they are entered.
104 // ///////////////////////////////////////////////////////////////////
105 //
106 // // Perform the conversion and execute subtasks (if any)
107 // // on an event-by-event basis
108 // q.ExecuteJob();
109 //
110 // // Select various objects to be added to the output file
111 //
112 // TFile* ofile=q.GetOutputFile();
113 // 
114 // if (ofile)
115 // {
116 //  ofile->cd(); // Switch to the output file directory
117 //
118 //  AliObjMatrix* omdb=q.GetOMdbase();
119 //  if (omdb) omdb->Write();
120 //
121 //  AliDevice* fitdefs=q.GetFitdefs();
122 //  if (fitdefs) fitdefs->Write();
123 //
124 //  TDatabasePDG* pdg=q.GetPDG();
125 //  if (pdg) pdg->Write();
126 //
127 //  // Flush additional objects to the output file.
128 //  // The output file is not explicitly closed here
129 //  // to allow interactive investigation of the data tree
130 //  // when this macro is run in an interactive ROOT/CINT session.
131 //  ofile->Write();
132 // }
133 //
134 //--- Author: Nick van Eijndhoven 11-mar-2005 Utrecht University
135 //- Modified: NvE $Date$ Utrecht University
136 ///////////////////////////////////////////////////////////////////////////
137  
138 #include "IceF2k.h"
139 #include "Riostream.h"
140
141 ClassImp(IceF2k) // Class implementation to enable ROOT I/O
142
143 IceF2k::IceF2k(const char* name,const char* title) : AliJob(name,title)
144 {
145 // Default constructor.
146 // By default maxevent=-1, split=0, bsize=32000, printfreq=1.
147
148  fSplit=0;
149  fBsize=32000;
150  fMaxevt=-1;
151  fPrintfreq=1;
152  fInfiles=0;
153  fOutfile=0;
154
155  fPdg=0;
156  fOmdb=0;
157  fFitdefs=0;
158  fTrigdefs=0;
159  fToffset=0;
160  fMctoffset=0;
161  fMctracks=3;
162 }
163 ///////////////////////////////////////////////////////////////////////////
164 IceF2k::~IceF2k()
165 {
166 // Default destructor.
167
168  if (fInfiles)
169  {
170   delete fInfiles;
171   fInfiles=0;
172  }
173
174  if (fPdg)
175  {
176   delete fPdg;
177   fPdg=0;
178  }
179
180  if (fOmdb)
181  {
182   delete fOmdb;
183   fOmdb=0;
184  }
185
186  if (fFitdefs)
187  {
188   delete fFitdefs;
189   fFitdefs=0;
190  }
191
192  if (fTrigdefs)
193  {
194   delete fTrigdefs;
195   fTrigdefs=0;
196  }
197 }
198 ///////////////////////////////////////////////////////////////////////////
199 void IceF2k::SetMaxEvents(Int_t n)
200 {
201 // Set the maximum number of events to be processed.
202 // n=-1 implies processing of the complete input file, which is the default
203 // initialisation in the constructor.
204  fMaxevt=n;
205 }
206 ///////////////////////////////////////////////////////////////////////////
207 void IceF2k::SetPrintFreq(Int_t f)
208 {
209 // Set the printfrequency to produce info every f events.
210 // f=1 is the default initialisation in the constructor.
211  if (f>=0) fPrintfreq=f;
212 }
213 ///////////////////////////////////////////////////////////////////////////
214 void IceF2k::SetSplitLevel(Int_t split)
215 {
216 // Set the split level for the ROOT data file.
217 // split=0 is the default initialisation in the constructor.
218  if (split>=0) fSplit=split;
219 }
220 ///////////////////////////////////////////////////////////////////////////
221 void IceF2k::SetBufferSize(Int_t bsize)
222 {
223 // Set the buffer size for the ROOT data file.
224 // bsize=32000 is the default initialisation in the constructor.
225  if (bsize>=0) fBsize=bsize;
226 }
227 ///////////////////////////////////////////////////////////////////////////
228 void IceF2k::SetMcToffset(Float_t toffset)
229 {
230 // Set a user defined time offset for Monte Carlo data.
231 // A very frequently (but not always) used value is -19000.
232 // See the introductory docs of this class for further details.
233  fMctoffset=toffset;
234 }
235 ///////////////////////////////////////////////////////////////////////////
236 void IceF2k::SelectMcTracks(Int_t mode)
237 {
238 // User selection of MC tracks to be stored in the event structure.
239 //
240 // mode = 0 : No MC tracks are stored
241 //        1 : Only muon and muon-neutrino MC tracks are stored
242 //        2 : All lepton MC tracks are stored
243 //        3 : All MC tracks (incl. brems, pairprod etc...) are stored
244 //
245 // By default mode=3 is set in the constructor of this class.
246
247  if (mode<0 || mode >3) return;
248  fMctracks=mode;
249 }
250 ///////////////////////////////////////////////////////////////////////////
251 void IceF2k::SetInputFile(TString name)
252 {
253 // Set the name of the F2K input file.
254 // This function has become obsolete but is kept for backward compatibility.
255 // The user is advised to use AddInputFile() instead, which allows processing
256 // of multiple F2K input files.
257 // This function will reset the list of all F2K input files and put the specified
258 // filename at the first position.
259 // Additional F2K input files can be specified via AddInputFile().
260
261  if (fInfiles) delete fInfiles;
262
263  fInfiles=new TObjArray();
264  fInfiles->SetOwner();
265
266  TObjString* s=new TObjString();
267  s->SetString(name);
268  fInfiles->Add(s);
269 }
270 ///////////////////////////////////////////////////////////////////////////
271 void IceF2k::AddInputFile(TString name)
272 {
273 // Add the name of this F2K input file to the list to be processed.
274
275  if (!fInfiles)
276  {
277   fInfiles=new TObjArray();
278   fInfiles->SetOwner();
279  }
280
281  TObjString* s=new TObjString();
282  s->SetString(name);
283  fInfiles->Add(s);
284 }
285 ///////////////////////////////////////////////////////////////////////////
286 void IceF2k::SetOutputFile(TFile* ofile)
287 {
288 // Set the output file for the ROOT data.
289  if (fOutfile) delete fOutfile;
290  fOutfile=ofile;
291 }
292 ///////////////////////////////////////////////////////////////////////////
293 void IceF2k::SetOutputFile(TString name)
294 {
295 // Create the output file for the ROOT data.
296  if (fOutfile) delete fOutfile;
297  fOutfile=new TFile(name.Data(),"RECREATE","F2K data in IceEvent structure");
298 }
299 ///////////////////////////////////////////////////////////////////////////
300 TFile* IceF2k::GetOutputFile()
301 {
302 // Provide pointer to the ROOT output file.
303  return fOutfile;
304 }
305 ///////////////////////////////////////////////////////////////////////////
306 TDatabasePDG* IceF2k::GetPDG()
307 {
308 // Provide pointer to the PDG database
309  return fPdg;
310 }
311 ///////////////////////////////////////////////////////////////////////////
312 AliObjMatrix* IceF2k::GetOMdbase()
313 {
314 // Provide pointer to the OM geometry, calib. etc... database
315  return fOmdb;
316 }
317 ///////////////////////////////////////////////////////////////////////////
318 AliDevice* IceF2k::GetFitdefs()
319 {
320 // Provide pointer to the fit definitions
321  return fFitdefs;
322 }
323 ///////////////////////////////////////////////////////////////////////////
324 AliDevice* IceF2k::GetTrigdefs()
325 {
326 // Provide pointer to the trigger definitions
327  return fTrigdefs;
328 }
329 ///////////////////////////////////////////////////////////////////////////
330 void IceF2k::Exec(Option_t* opt)
331 {
332 // Job to loop over the specified number of events and convert the 
333 // F2K data into the IceEvent structure.
334 // If maxevents<0 (default) all the entries of the input file
335 // will be processed.
336 // Every "printfreq" events a short event summary will be printed.
337 // The default value is printfreq=1.
338 // The output will be written on a standard output tree named "T".
339 //
340 // Notes :
341 // -------
342 // 1) This class is derived from AliJob, allowing a task based processing.
343 //    After the conversion of an F2K event into an IceEvent structure,
344 //    the processing of all available sub-tasks (if any) is invoked.
345 //    This provides an event-by-event (sub)task processing before the
346 //    final data structures are written out.
347 // 2) The main object in this job environment is an IceEvent* pointer.
348
349  if (!fInfiles)
350  {
351   cout << " *IceF2k Exec* No data input file(s) specified." << endl;
352   return;
353  }
354
355  Int_t ninfiles=fInfiles->GetEntries();
356  if (!ninfiles)
357  {
358   cout << " *IceF2k Exec* No data input file(s) specified." << endl;
359   return;
360  }
361
362  TTree* otree=0;
363  if (fOutfile)
364  {
365   otree=new TTree("T","F2K Data converted to IceEvent structures");
366   otree->SetDirectory(fOutfile);
367  }
368
369  IceEvent* evt=new IceEvent();
370  evt->SetTrackCopy(1);
371  evt->SetDevCopy(1);
372
373  // Branch in the tree for the event structure
374  if (otree) otree->Branch("IceEvent","IceEvent",&evt,fBsize,fSplit); 
375
376  // Create the particle database and extend it with some F2000 specific definitions
377  if (!fPdg) fPdg=new TDatabasePDG();
378  Double_t me=fPdg->GetParticle(11)->Mass();
379  fPdg->AddParticle("brems"   ,"brems"   ,0,1,0,0,"none",10001001,0,0);
380  fPdg->AddParticle("deltae"  ,"deltae"  ,me,1,0,-3,"Lepton",10001002,0,0);
381  fPdg->AddParticle("pairprod","pairprod",0,1,0,0,"none",10001003,0,0);
382  fPdg->AddParticle("nucl_int","nucl_Int",0,1,0,0,"none",10001004,0,0);
383  fPdg->AddParticle("mu_pair" ,"mu_pair" ,0,1,0,0,"none",10001005,0,0);
384  fPdg->AddParticle("hadrons" ,"hadrons" ,0,1,0,0,"none",10001006,0,0);
385  fPdg->AddParticle("fiberlaser","fiberlaser",0,1,0,0,"none",10002100,0,0);
386  fPdg->AddParticle("n2laser"   ,"n2laser"   ,0,1,0,0,"none",10002101,0,0);
387  fPdg->AddParticle("yaglaser"  ,"yaglaser"  ,0,1,0,0,"none",10002201,0,0);
388  fPdg->AddParticle("z_primary","z_primary",0,1,0,0,"none",10003000,0,0);
389  fPdg->AddParticle("a_primary","a_primary",0,1,0,0,"none",10003500,0,0);
390
391  // Initialise the job working environment
392  SetMainObject(evt);
393  if (fOutfile)
394  {
395   AddObject(fOutfile);
396   AddObject(otree);
397  }
398
399  TString inputfile;
400
401  cout << " ***" << endl;
402  cout << " *** Start processing of job " << GetName() << " ***" << endl;
403  cout << " ***" << endl;
404  for (Int_t i=0; i<ninfiles; i++)
405  {
406   TObjString* sx=(TObjString*)fInfiles->At(i);
407   if (!sx) continue;
408   inputfile=sx->GetString(); 
409   cout << " F2K input file : " << inputfile.Data() << endl;
410  }
411  cout << " Maximum number of events to be processed : " << fMaxevt << endl;
412  cout << " Print frequency : " << fPrintfreq << endl;
413  if (fOutfile)
414  {
415   cout << " ROOT output file : " << fOutfile->GetName() << endl;
416   cout << " Output characteristics : splitlevel = " << fSplit << " buffersize = " << fBsize << endl;
417  }
418
419  ListEnvironment();
420
421  Int_t nevt=0;
422  for (Int_t ifile=0; ifile<ninfiles; ifile++)
423  {
424   TObjString* sx=(TObjString*)fInfiles->At(ifile);
425   if (!sx) continue;
426
427   inputfile=sx->GetString(); 
428   if (inputfile=="") continue;
429
430   // Open the input file in the default ascii format (autodetection) for reading 
431   fInput=rdmc_mcopen(inputfile.Data(),"r",RDMC_DEFAULT_ASCII_F);
432
433   if (!fInput)
434   {
435    cout << " *IceF2k Exec* No input file found with name : " << inputfile.Data() << endl;
436    continue;
437   }
438
439   // Initialise the event structure 
440   rdmc_init_mevt(&fEvent);
441
442   // Read the file header information
443   rdmc_rarr(fInput,&fHeader);
444
445   // Fill the database with geometry, calib. etc... parameters
446   // for all the devices
447   FillOMdbase();
448
449   // Set the fit definitions according to the F2000 header info
450   SetFitdefs();
451
452   // Set the trigger definitions according to the F2000 header info
453   SetTrigdefs();
454  
455   while (!rdmc_revt(fInput,&fHeader,&fEvent))
456   {
457    if (fMaxevt>-1 && nevt>=fMaxevt) break;
458
459    // Reset the complete Event structure
460    evt->Reset();
461
462    evt->SetRunNumber(fEvent.nrun);
463    evt->SetEventNumber(fEvent.enr);
464    evt->SetMJD(fEvent.mjd,fEvent.secs,fEvent.nsecs);
465
466    // Take trigger offset into account which might have been
467    // introduced during the filtering process.
468    // For simulated data this will be treated separately in PutMcTracks().
469    fToffset=fEvent.t_offset;
470
471    PutTrigger();
472
473    PutMcTracks();
474
475    PutRecoTracks();
476
477    PutHits();
478
479    // Invoke all available sub-tasks (if any)
480    CleanTasks();
481    ExecuteTasks(opt);
482
483    if (fPrintfreq)
484    {
485     if (!(nevt%fPrintfreq)) evt->HeaderData();
486    }
487
488    // Write the complete structure to the output Tree
489    if (otree) otree->Fill();
490
491    // Update event counter
492    nevt++;
493   }
494   if (fMaxevt>-1 && nevt>=fMaxevt) break;
495  }
496
497  // Flush possible memory resident data to the output file
498  if (fOutfile) fOutfile->Write();
499
500  // Remove the IceEvent object from the environment
501  // and delete it as well
502  if (evt)
503  {
504   RemoveObject(evt);
505   delete evt;
506  }
507 }
508 ///////////////////////////////////////////////////////////////////////////
509 void IceF2k::FillOMdbase()
510 {
511 // Fill the database with geometry, calib. etc... parameters 
512 // for all the devices.
513
514  if (fHeader.nch<=0)
515  {
516   if (fOmdb)
517   {
518    delete fOmdb;
519    fOmdb=0;
520   }
521   return;
522  }
523
524  Int_t adccal=fHeader.is_calib.adc;
525  Int_t tdccal=fHeader.is_calib.tdc;
526  Int_t totcal=fHeader.is_calib.tot;
527
528  TF1 fadccal("fadccal","(x-[1])*[0]");
529  TF1 fadcdecal("fadcdecal","(x/[0])+[1]");
530  fadccal.SetParName(0,"BETA-ADC");
531  fadccal.SetParName(1,"PED-ADC");
532  fadcdecal.SetParName(0,"BETA-ADC");
533  fadcdecal.SetParName(1,"PED-ADC");
534
535  TF1 ftdccal("ftdccal","(x*[0])-[1]-([0]-1.)*32767.-[2]/sqrt([3])");
536  TF1 ftdcdecal("ftdcdecal","(x+([0]-1.)*32767.+[1]+[2]/sqrt([3]))/[0]");
537  ftdccal.SetParName(0,"BETA-TDC");
538  ftdccal.SetParName(1,"T0");
539  ftdccal.SetParName(2,"ALPHA-TDC");
540  ftdccal.SetParName(3,"ADC-SLEW");
541  ftdcdecal.SetParName(0,"BETA-TDC");
542  ftdcdecal.SetParName(1,"T0");
543  ftdcdecal.SetParName(2,"ALPHA-TDC");
544  ftdcdecal.SetParName(3,"ADC-SLEW");
545
546  TF1 ftotcal("ftotcal","x*[0]");
547  TF1 ftotdecal("ftotdecal","x/[0]");
548  ftotcal.SetParName(0,"BETA-TOT");
549  ftotdecal.SetParName(0,"BETA-TOT");
550
551  if (fOmdb)
552  {
553   fOmdb->Reset();
554  }
555  else
556  {
557   fOmdb=new AliObjMatrix();
558   fOmdb->SetNameTitle("OMDBASE","The OM geometry, calib. etc... database");
559   fOmdb->SetOwner();
560  }
561
562  IceAOM* dev=0;
563  Double_t pos[3]={0,0,0};
564  for (Int_t i=0; i<fHeader.nch; i++)
565  {
566   dev=new IceAOM();
567   dev->SetUniqueID(i+1);
568
569   dev->SetSlotName("ADC",1);
570   dev->SetSlotName("LE",2);
571   dev->SetSlotName("TOT",3);
572
573   dev->SetSlotName("TYPE",4);
574   dev->SetSlotName("ORIENT",5);
575   dev->SetSlotName("THRESH",6);
576   dev->SetSlotName("SENSIT",7);
577   dev->SetSlotName("BETA-TDC",8);
578   dev->SetSlotName("T0",9);
579   dev->SetSlotName("ALPHA-TDC",10);
580   dev->SetSlotName("PED-ADC",11);
581   dev->SetSlotName("BETA-ADC",12);
582   dev->SetSlotName("KAPPA-ADC",13);
583   dev->SetSlotName("PED-TOT",14);
584   dev->SetSlotName("BETA-TOT",15);
585   dev->SetSlotName("KAPPA-TOT",16);
586
587   pos[0]=fHeader.x[i];
588   pos[1]=fHeader.y[i];
589   pos[2]=fHeader.z[i];
590   dev->SetPosition(pos,"car");
591
592   fadccal.SetParameter(0,fHeader.cal[i].beta_a);
593   fadccal.SetParameter(1,fHeader.cal[i].ped);
594   fadcdecal.SetParameter(0,fHeader.cal[i].beta_a);
595   if (!fHeader.cal[i].beta_a) fadcdecal.SetParameter(0,1);
596   fadcdecal.SetParameter(1,fHeader.cal[i].ped);
597
598   ftdccal.SetParameter(0,fHeader.cal[i].beta_t);
599   ftdccal.SetParameter(1,fHeader.cal[i].t_0);
600   ftdccal.SetParameter(2,fHeader.cal[i].alpha_t);
601   ftdccal.SetParameter(3,1.e20);
602   ftdcdecal.SetParameter(0,fHeader.cal[i].beta_t);
603   if (!fHeader.cal[i].beta_t) ftdcdecal.SetParameter(0,1);
604   ftdcdecal.SetParameter(1,fHeader.cal[i].t_0);
605   ftdcdecal.SetParameter(2,fHeader.cal[i].alpha_t);
606   ftdcdecal.SetParameter(3,1.e20);
607
608   ftotcal.SetParameter(0,fHeader.cal[i].beta_tot);
609   ftotdecal.SetParameter(0,fHeader.cal[i].beta_tot);
610   if (!fHeader.cal[i].beta_tot) ftotdecal.SetParameter(0,1);
611
612   if (adccal)
613   {
614    dev->SetDecalFunction(&fadcdecal,1);
615   }
616   else
617   {
618    dev->SetCalFunction(&fadccal,1);
619   }
620
621   if (tdccal)
622   {
623    dev->SetDecalFunction(&ftdcdecal,2);
624   }
625   else
626   {
627    dev->SetCalFunction(&ftdccal,2);
628   }
629
630   if (totcal)
631   {
632    dev->SetDecalFunction(&ftotdecal,3);
633   }
634   else
635   {
636    dev->SetCalFunction(&ftotcal,3);
637   }
638
639   dev->SetSignal(fHeader.type[i],4);
640   dev->SetSignal((Float_t)fHeader.costh[i],5);
641   dev->SetSignal(fHeader.thresh[i],6);
642   dev->SetSignal(fHeader.sensit[i],7);
643   dev->SetSignal(fHeader.cal[i].beta_t,8);
644   dev->SetSignal(fHeader.cal[i].t_0,9);
645   dev->SetSignal(fHeader.cal[i].alpha_t,10);
646   dev->SetSignal(fHeader.cal[i].ped,11);
647   dev->SetSignal(fHeader.cal[i].beta_a,12);
648   dev->SetSignal(fHeader.cal[i].kappa,13);
649   dev->SetSignal(fHeader.cal[i].ped_tot,14);
650   dev->SetSignal(fHeader.cal[i].beta_tot,15);
651   dev->SetSignal(fHeader.cal[i].kappa_tot,16);
652
653   fOmdb->EnterObject(i+1,1,dev);
654  }
655 }
656 ///////////////////////////////////////////////////////////////////////////
657 void IceF2k::SetFitdefs()
658 {
659 // Obtain the names of the variables for each fit procedure from the
660 // f2000 header. Each different fit procedure is then stored as a separate
661 // "hit" of an AliDevice object and the various fit variables are stored
662 // as separate signal slots of the corresponding "hit".
663 // Via the GetFitdefs() memberfunction this AliDevice object can be
664 // retrieved and stored in the ROOT output file if wanted.
665 // The name of the object is FitDefinitions and the stored data can be
666 // inspected via the AliDevice::Data() memberfunction and looks as follows :
667 //
668 //  *AliDevice::Data* Id :0 Name : FitDefinitions
669 //    Position Vector in car coordinates : 0 0 0
670 //    Err. in car coordinates : 0 0 0
671 //  The following 8 hits are registered :
672 //  *AliSignal::Data* Id :0
673 //    Position Vector in car coordinates : 0 0 0
674 //    Err. in car coordinates : 0 0 0
675 //    Owned by device : AliDevice Name : FitDefinitions
676 //    Slot : 1 Signal value : 0 name : id
677 //    Slot : 2 Signal value : 0 name : rchi2
678 //    Slot : 3 Signal value : 0 name : prob
679 //    Slot : 4 Signal value : 0 name : sigth
680 //    Slot : 5 Signal value : 0 name : covmin
681 //    Slot : 6 Signal value : 0 name : covmax
682 //    Slot : 7 Signal value : 0 name : cutflag
683 //    Slot : 8 Signal value : 0 name : chi2
684 //  *AliSignal::Data* Id :1
685 //    Position Vector in car coordinates : 0 0 0
686 //    Err. in car coordinates : 0 0 0
687 //    Owned by device : AliDevice Name : FitDefinitions
688 //    Slot : 1 Signal value : 0 name : id
689 //    Slot : 2 Signal value : 0 name : rchi2
690 //    Slot : 3 Signal value : 0 name : prob
691 // etc....  
692 //
693 // This memberfunction is based on the original idea/code by Adam Bouchta.
694
695  if (fHeader.n_fit<=0)
696  {
697   if (fFitdefs)
698   {
699    delete fFitdefs;
700    fFitdefs=0;
701   }
702   return;
703  }
704
705  if (fFitdefs)
706  {
707   fFitdefs->Reset(1);
708  }
709  else
710  {
711   fFitdefs=new AliDevice();
712  }
713
714  fFitdefs->SetName("FitDefinitions");
715  fFitdefs->SetHitCopy (1);
716
717  AliSignal s;
718  s.Reset();
719
720  for (Int_t i=0; i<fHeader.n_fit; i++)
721  {
722   s.SetUniqueID(fHeader.def_fit[i].id);
723   s.SetName(TString(fHeader.def_fit[i].tag));
724
725   for (Int_t j=0; j<fHeader.def_fit[i].nwords; j++)
726   {
727    s.SetSlotName(TString(fHeader.def_fit[i].words[j]),j+1);
728    s.SetSignal(0,j+1);
729   }
730
731   fFitdefs->AddHit(s);
732   s.Reset(1);
733  }
734 }
735 ///////////////////////////////////////////////////////////////////////////
736 void IceF2k::SetTrigdefs()
737 {
738 // Obtain the names of the variables for each trigger procedure from the
739 // f2000 header. Each different trigger procedure is then stored as a separate
740 // "hit" of an AliDevice object and the various trigger variables are stored
741 // as separate signal slots of the corresponding "hit".
742 // Via the GetFitdefs() memberfunction this AliDevice object can be
743 // retrieved and stored in the ROOT output file if wanted.
744 // The name of the object is TrigDefinitions and the stored data can be
745 // inspected via the AliDevice::Data() memberfunction and looks as follows :
746 //
747 //  *AliDevice::Data* Id : 0 Name : TrigDefinitions
748 //    Position Vector in car (rad) coordinates : 0 0 0
749 //    Err. in car (rad) coordinates : 0 0 0
750 //  The following 9 hits are registered : 
751 //  *AliSignal::Data* Id : 1 Name : main
752 //    Position Vector in car (rad) coordinates : 0 0 0
753 //    Err. in car (rad) coordinates : 0 0 0
754 //    Owned by device : AliDevice Id : 0 Name : TrigDefinitions
755 //     Slot : 1 Signal value : 0 name : trig_pulse_le
756 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
757 //     Slot : 3 Signal value : 0 name : regi_flag
758 //   *AliSignal::Data* Id : 2 Name : amaa
759 //     Position Vector in car (rad) coordinates : 0 0 0
760 //     Err. in car (rad) coordinates : 0 0 0
761 //     Owned by device : AliDevice Id : 0 Name : TrigDefinitions
762 //     Slot : 1 Signal value : 0 name : trig_pulse_le
763 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
764 //     Slot : 3 Signal value : 0 name : regi_flag
765 //   *AliSignal::Data* Id : 3 Name : amab10
766 //     Position Vector in car (rad) coordinates : 0 0 0
767 //     Err. in car (rad) coordinates : 0 0 0
768 //     Owned by device : AliDevice Id : 0 Name : TrigDefinitions
769 //     Slot : 1 Signal value : 0 name : trig_pulse_le
770 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
771 //     Slot : 3 Signal value : 0 name : regi_flag
772 // etc....  
773
774  if (fHeader.n_trigger<=0)
775  {
776   if (fTrigdefs)
777   {
778    delete fTrigdefs;
779    fTrigdefs=0;
780   }
781   return;
782  }
783
784  if (fTrigdefs)
785  {
786   fTrigdefs->Reset(1);
787  }
788  else
789  {
790   fTrigdefs=new AliDevice();
791  }
792
793  fTrigdefs->SetName("TrigDefinitions");
794  fTrigdefs->SetHitCopy (1);
795
796  AliSignal s;
797  s.Reset();
798
799  for (Int_t i=0; i<fHeader.n_trigger; i++)
800  {
801   s.SetUniqueID(fHeader.def_trig[i].id);
802   s.SetName(TString(fHeader.def_trig[i].tag));
803
804   for (Int_t j=0; j<fHeader.def_trig[i].nwords; j++)
805   {
806    s.SetSlotName(TString(fHeader.def_trig[i].words[j]),j+1);
807    s.SetSignal(0,j+1);
808   }
809
810   fTrigdefs->AddHit(s);
811   s.Reset(1);
812  }
813 }
814 ///////////////////////////////////////////////////////////////////////////
815 void IceF2k::PutMcTracks()
816 {
817 // Get the MC tracks from the F2000 file into the IcePack structure.
818 // Note : MC tracks are given negative track id's in the event structure.
819 // This memberfunction is based on the original code by Adam Bouchta.
820
821  IceEvent* evt=(IceEvent*)GetMainObject();
822  if (!evt || fEvent.ntrack<=0) return;
823
824  // User defined trigger offset in case of simulated data.
825  // The offset in the F2K file is meant to put the primary interaction
826  // well ahead of the detector trigger.
827  // See the introductory docs of this IceF2k class for further details.
828  fToffset=fMctoffset;
829
830  if (!fMctracks) return;
831
832  // Loop over all the tracks and add them to the current event
833  AliTrack t;
834  Double_t vec[3];
835  AliPosition r;
836  Ali3Vector p;
837  Int_t tid=0;
838  Int_t idpdg=0;
839  Int_t idf2k=0;
840  for (Int_t i=0; i<fEvent.ntrack; i++)
841  {
842   t.Reset ();
843
844   // Beginpoint of the track
845   vec[0]=fEvent.gen[i].x;
846   vec[1]=fEvent.gen[i].y;
847   vec[2]=fEvent.gen[i].z;
848   r.SetPosition(vec,"car");
849   t.SetBeginPoint(r);
850
851   // Endpoint of the track
852   vec[0]+=fEvent.gen[i].length*fEvent.gen[i].px;
853   vec[1]+=fEvent.gen[i].length*fEvent.gen[i].py;
854   vec[2]+=fEvent.gen[i].length*fEvent.gen[i].pz;
855   r.SetPosition(vec,"car");
856   t.SetEndPoint(r);
857
858   // Momentum in GeV/c
859   vec[0]=fEvent.gen[i].e*fEvent.gen[i].px*1e-3;
860   vec[1]=fEvent.gen[i].e*fEvent.gen[i].py*1e-3;
861   vec[2]=fEvent.gen[i].e*fEvent.gen[i].pz*1e-3;
862   p.SetVector (vec,"car");
863   t.Set3Momentum(p);
864
865   // MC tracks are indicated by negative track id's
866   tid=fEvent.gen[i].tag;
867   t.SetId(-abs(tid));
868
869   idf2k=fEvent.gen[i].id;
870   idpdg=0;
871   if (idf2k>1000)
872   {
873    idpdg=idf2k+10000000;
874   }
875   else if (idf2k <= 48)
876   {
877    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
878   }
879   else
880   {
881    if (idf2k==201) idpdg=12;
882    if (idf2k==202) idpdg=14;
883    if (idf2k==203) idpdg=16;
884    if (idf2k==204) idpdg=-12;
885    if (idf2k==205) idpdg=-14;
886    if (idf2k==206) idpdg=-16;
887   }
888
889   // Check for the user selected MC track storage
890   if (fMctracks==1) // Store only muon and muon-neutrino tracks
891   {
892    if (abs(idpdg)!=13 && abs(idpdg)!=14) continue;
893   }
894   else if (fMctracks==2) // Store all lepton tracks
895   {
896    if (abs(idpdg)<11 || abs(idpdg)>16) continue;
897   }
898
899   t.SetParticleCode(idpdg);
900   t.SetName(fPdg->GetParticle(idpdg)->GetName());
901   t.SetTitle("MC track");
902   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
903   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
904
905   evt->AddTrack(t);
906  }
907
908  // Create the pointers to each particle's parent particle.
909  Int_t txid=0;
910  Int_t parid=0;
911  for (Int_t itk=1; itk<=evt->GetNtracks (); itk++)
912  {
913   AliTrack* tx=evt->GetTrack(itk);
914
915   if (!tx) continue;
916
917   txid=tx->GetId();
918
919   parid=-1;
920   for (Int_t j=0; j<fEvent.ntrack; j++)
921   {
922    tid=fEvent.gen[j].tag;
923    if (-abs(tid) == txid) parid=fEvent.gen[j].parent;
924   }
925
926   if (parid<0) continue;
927
928   AliTrack* tpar=evt->GetIdTrack(-abs(parid));
929
930   if (!tpar) continue;
931
932   tx->SetParentTrack(tpar);
933  }
934 }
935 ///////////////////////////////////////////////////////////////////////////
936 void IceF2k::PutRecoTracks()
937 {
938 // Get the reconstructed tracks from the F2000 file into the IcePack structure.
939 // Note : Reco tracks are given positive track id's in the event structure.
940 // This memberfunction is based on the original code by Adam Bouchta.
941
942  IceEvent* evt=(IceEvent*)GetMainObject();
943  if (!evt || fEvent.nfit<=0) return;
944
945  // Loop over all the tracks and add them to the current event
946  AliTrack t;
947  Double_t vec[3];
948  AliPosition r;
949  Ali3Vector p;
950  Int_t tid=0;
951  Int_t idpdg=0;
952  Int_t idf2k=0;
953  for (Int_t i=0; i<fEvent.nfit; i++)
954  {
955   t.Reset ();
956
957   // Beginpoint of the track
958   vec[0]=fEvent.rec[i].x;
959   vec[1]=fEvent.rec[i].y;
960   vec[2]=fEvent.rec[i].z;
961   r.SetPosition(vec,"car");
962   t.SetBeginPoint(r);
963
964   // Endpoint of the track
965   vec[0]+=fEvent.rec[i].length*fEvent.rec[i].px;
966   vec[1]+=fEvent.rec[i].length*fEvent.rec[i].py;
967   vec[2]+=fEvent.rec[i].length*fEvent.rec[i].pz;
968   r.SetPosition(vec,"car");
969   t.SetEndPoint(r);
970
971   // Momentum in GeV/c
972   if (fEvent.rec[i].e > 0)
973   {
974    vec[0]=fEvent.rec[i].e*fEvent.rec[i].px*1e-3;
975    vec[1]=fEvent.rec[i].e*fEvent.rec[i].py*1e-3;
976    vec[2]=fEvent.rec[i].e*fEvent.rec[i].pz*1e-3;
977   }
978   else // Give the track a nominal momentum of 1 GeV/c
979   {
980    vec[0]=fEvent.rec[i].px;
981    vec[1]=fEvent.rec[i].py;
982    vec[2]=fEvent.rec[i].pz;
983   }
984   p.SetVector (vec,"car");
985   t.Set3Momentum(p);
986
987   // Use the fit number as track id
988   tid=fEvent.rec[i].tag;
989   t.SetId(abs(tid));
990
991   idf2k=fEvent.rec[i].id;
992   idpdg=0;
993   if (idf2k>1000)
994   {
995    idpdg=idf2k+10000000;
996   }
997   else if (idf2k <= 48)
998   {
999    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
1000   }
1001   else
1002   {
1003    if (idf2k==201) idpdg=12;
1004    if (idf2k==202) idpdg=14;
1005    if (idf2k==203) idpdg=16;
1006    if (idf2k==204) idpdg=-12;
1007    if (idf2k==205) idpdg=-14;
1008    if (idf2k==206) idpdg=-16;
1009   }
1010
1011   t.SetParticleCode(idpdg);
1012   t.SetNameTitle("Sieglinde","RECO track");
1013   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
1014   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
1015
1016   // Retrieve the various fit parameters for this track
1017   AliSignal* fitdata=fFitdefs->GetIdHit(i);
1018   for (Int_t jval=0; jval<fEvent.fresult[i].nval; jval++)
1019   {
1020    fitdata->SetSignal(fEvent.fresult[i].val[jval],jval+1);
1021   }
1022
1023   // Store the various fit parameters for this track
1024   t.SetFitDetails(fitdata);
1025
1026   // Store the various reco tracks as track hypotheses.
1027   // A copy of the first reco track is entered as a new track instance
1028   // into the event and all reco tracks (incl. the first one) are
1029   // stored as hypotheses linked to this new reco track.
1030   if (i==0)
1031   {
1032    evt->AddTrack(t);
1033    AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1034    Int_t nrec=evt->GetNtracks(1);
1035    tx->SetId(nrec+1);
1036   }
1037   AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1038   if (tx) tx->AddTrackHypothesis(t);
1039  }
1040 }
1041 ///////////////////////////////////////////////////////////////////////////
1042 void IceF2k::PutHits()
1043 {
1044 // Get the hit and waveform info from the F2000 file into the IcePack structure.
1045 // This memberfunction is based on the original code by Adam Bouchta.
1046
1047  IceEvent* evt=(IceEvent*)GetMainObject();
1048  if (!evt) return;
1049
1050  // Loop over all the hits and add them to the current event
1051  IceAOM om;
1052  AliSignal s;
1053  s.SetSlotName("ADC",1);
1054  s.SetSlotName("LE",2);
1055  s.SetSlotName("TOT",3);
1056  Int_t chan=0;
1057  Int_t maxchan=800;
1058  if (fOmdb) maxchan=fHeader.nch;
1059  IceAOM* omx=0;
1060  AliSignal* sx=0;
1061  Int_t tid=0;
1062  AliTrack* tx=0;
1063  Float_t adc=0;
1064  Float_t adcfirst=0; // Adc value of the first hit of an OM
1065  for (Int_t i=0; i<fEvent.nhits; i++)
1066  {
1067   chan=fEvent.h[i].ch+1;
1068   if (chan>maxchan) continue; // Channels 9001, 9002 etc are trigger channels
1069
1070   // Get corresponding device from the current event structure  
1071   omx=(IceAOM*)evt->GetIdDevice(chan);
1072   if (!omx)
1073   {
1074    if (fOmdb)
1075    {
1076     omx=(IceAOM*)fOmdb->GetObject(chan,1);
1077     evt->AddDevice(omx);
1078    }
1079    else
1080    {
1081     om.Reset(1);
1082     om.SetUniqueID(chan);
1083     evt->AddDevice(om);
1084    }
1085    omx=(IceAOM*)evt->GetIdDevice(chan);
1086   }
1087
1088   if (!omx) continue;
1089
1090   adc=fEvent.h[i].amp;
1091
1092   // Multiple hits in the same OM with the same ADC value
1093   // are indicated by "*" in the F2K file.
1094   // This corresponds to a value of -2 in the data structure.
1095   if (int(adc) == -2)
1096   {
1097    adc=adcfirst;
1098   }
1099   else
1100   {
1101    adcfirst=adc;
1102   }
1103   s.Reset();
1104   s.SetUniqueID(fEvent.h[i].id);
1105   s.SetSignal(adc,1);
1106   s.SetSignal((fEvent.h[i].t-fToffset),2);
1107   s.SetSignal(fEvent.h[i].tot,3);
1108
1109   omx->AddHit(s);
1110
1111   sx=omx->GetHit(omx->GetNhits());
1112   if (!sx) continue;
1113
1114   // ADC dependent TDC (de)calibration function for this hit
1115   TF1* fcal=omx->GetCalFunction("LE");
1116   TF1* fdecal=omx->GetDecalFunction("LE");
1117   if (fcal) sx->SetCalFunction(fcal,2);
1118   if (fdecal) sx->SetDecalFunction(fdecal,2);
1119   fcal=sx->GetCalFunction(2);
1120   fdecal=sx->GetDecalFunction(2);
1121   adc=sx->GetSignal(1,-4);
1122   if (adc>0)
1123   {
1124    if (fcal) fcal->SetParameter(3,adc);
1125    if (fdecal) fdecal->SetParameter(3,adc);
1126   }
1127   else
1128   {
1129    if (fcal) fcal->SetParameter(3,1.e20);
1130    if (fdecal) fdecal->SetParameter(3,1.e20);
1131   }
1132
1133   // Bi-directional link between this hit and the track that caused the ADC value.
1134   // This F2K info is probably only present for MC tracks.
1135   tid=fEvent.h[i].ma;
1136   if (tid > 0)
1137   {
1138    tx=evt->GetIdTrack(tid); // Reco tracks
1139    if (!tx) tx=evt->GetIdTrack(-tid); // MC tracks
1140    if (tx) sx->AddTrack(*tx);
1141   }
1142   else
1143   {
1144    if (tid == -2) sx->SetNameTitle("N","Noise");
1145    if (tid == -3) sx->SetNameTitle("A","Afterpulse");
1146   }
1147  }
1148
1149  // Loop over all the waveforms and add the histo(s) to the corresponding OM's
1150  TH1F histo;
1151  Int_t nbins=0;
1152  Float_t xlow=0;
1153  Float_t xup=0;
1154  TString hname;
1155  for (Int_t iwf=0; iwf<fEvent.nwf; iwf++)
1156  {
1157   chan=fEvent.wf[iwf].om;
1158   if (chan<=0 || chan>maxchan) continue; // Skip trigger channels
1159
1160   // Get corresponding device from the current event structure  
1161   omx=(IceAOM*)evt->GetIdDevice(chan);
1162   if (!omx)
1163   {
1164    if (fOmdb)
1165    {
1166     omx=(IceAOM*)fOmdb->GetObject(chan,1);
1167     evt->AddDevice(omx);
1168    }
1169    else
1170    {
1171     om.Reset(1);
1172     om.SetUniqueID(chan);
1173     evt->AddDevice(om);
1174    }
1175    omx=(IceAOM*)evt->GetIdDevice(chan);
1176   }
1177
1178   if (!omx) continue;
1179
1180   omx->SetSlotName("BASELINE",omx->GetNnames()+1);
1181   omx->SetSignal(-fEvent.wf[iwf].baseline,"BASELINE");
1182
1183   // Fill the waveform histogram
1184   hname="OM";
1185   hname+=chan;
1186   hname+="-WF";
1187   hname+=omx->GetNwaveforms()+1;
1188
1189   histo.Reset();
1190   histo.SetName(hname.Data());
1191   nbins=fEvent.wf[iwf].ndigi;
1192   xlow=fEvent.wf[iwf].t_start;
1193   xup=xlow+float(nbins)*fEvent.wf[iwf].t_bin;
1194   histo.SetBins(nbins,xlow,xup);
1195
1196   for (Int_t jbin=1; jbin<=fEvent.wf[iwf].ndigi; jbin++)
1197   {
1198    histo.SetBinContent(jbin,-fEvent.wf[iwf].digi[jbin-1]);
1199   }
1200
1201   omx->SetWaveform(&histo,omx->GetNwaveforms()+1);
1202  }
1203
1204  // Set bi-directional links between hits and reco track hypotheses.
1205  // Note : Reco tracks are recognised by their positive id.
1206  Int_t hid=0;
1207  TObjArray* rectracks=evt->GetTracks(1);
1208  for (Int_t jtk=0; jtk<rectracks->GetEntries(); jtk++)
1209  {
1210   tx=(AliTrack*)rectracks->At(jtk);
1211   if (!tx) continue;
1212   
1213   for (Int_t jhyp=1; jhyp<=tx->GetNhypotheses(); jhyp++)
1214   {
1215    AliTrack* hypx=tx->GetTrackHypothesis(jhyp);
1216    if (!hypx) continue;
1217
1218    // Loop over all combinations of F2K fits and used OM hits
1219    for (Int_t k=0; k<fEvent.nfit_uses; k++)
1220    {
1221     if (fEvent.fit_uses[k].useid != hypx->GetId()) continue;
1222     hid=fEvent.fit_uses[k].hitid;
1223     sx=evt->GetIdHit(hid,"IceAOM");
1224     if (sx) sx->AddTrack(*hypx);
1225    }
1226   }
1227  }
1228 }
1229 ///////////////////////////////////////////////////////////////////////////
1230 void IceF2k::PutTrigger()
1231 {
1232 // Get the trigger info from the F2000 file into the IcePack structure.
1233
1234  if (!fTrigdefs) return;
1235
1236  IceEvent* evt=(IceEvent*)GetMainObject();
1237  if (!evt || fEvent.ntrig<=0) return;
1238
1239  AliDevice trig;
1240  trig.SetNameTitle("Trigger","Amanda/IceCube event triggers");
1241  AliSignal s;
1242  TString trigname;
1243  TString slotname;
1244  Int_t id=0;
1245  Int_t nval=0;
1246  for (Int_t i=0; i<fEvent.ntrig; i++)
1247  {
1248   id=fEvent.ptrig[i].id;
1249   nval=fEvent.ptrig[i].nval;
1250   if (!nval) continue;
1251   AliSignal* tdef=fTrigdefs->GetIdHit(id+1);
1252   if (!tdef) continue;
1253   trigname=tdef->GetName();
1254   s.Reset(1);
1255   s.SetName(trigname);
1256   s.SetUniqueID(id+1);
1257   for (Int_t jval=0; jval<fEvent.ptrig[i].nval; jval++)
1258   {
1259    slotname=tdef->GetSlotName(jval+1);
1260    s.SetSlotName(slotname,jval+1);
1261    s.SetSignal(fEvent.ptrig[i].val[jval],jval+1);
1262   }
1263   trig.AddHit(s);
1264  }
1265
1266  // Store the trigger data into the IceEvent structure
1267  evt->AddDevice(trig);
1268 }
1269 ///////////////////////////////////////////////////////////////////////////