]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/icepack/iceconvert/IceF2k.cxx
25-sep-2006 NvE AliSample extended with memberfunction GetSpread() and also
[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  // Storage of the used parameters in the IceF2k device
392  AliSignal params;
393  params.SetNameTitle("IceF2k","IceF2k processor parameters");
394  params.SetSlotName("Toffset",1);
395  params.SetSlotName("Mctoffset",2);
396  params.SetSlotName("Mctracks",3);
397
398  // Initialise the job working environment
399  SetMainObject(evt);
400  if (fOutfile)
401  {
402   AddObject(fOutfile);
403   AddObject(otree);
404  }
405
406  TString inputfile;
407
408  cout << " ***" << endl;
409  cout << " *** Start processing of job " << GetName() << " ***" << endl;
410  cout << " ***" << endl;
411  for (Int_t i=0; i<ninfiles; i++)
412  {
413   TObjString* sx=(TObjString*)fInfiles->At(i);
414   if (!sx) continue;
415   inputfile=sx->GetString(); 
416   cout << " F2K input file : " << inputfile.Data() << endl;
417  }
418  cout << " Maximum number of events to be processed : " << fMaxevt << endl;
419  cout << " Print frequency : " << fPrintfreq << endl;
420  if (fOutfile)
421  {
422   cout << " ROOT output file : " << fOutfile->GetName() << endl;
423   cout << " Output characteristics : splitlevel = " << fSplit << " buffersize = " << fBsize << endl;
424  }
425
426  ListEnvironment();
427
428  Int_t nevt=0;
429  for (Int_t ifile=0; ifile<ninfiles; ifile++)
430  {
431   TObjString* sx=(TObjString*)fInfiles->At(ifile);
432   if (!sx) continue;
433
434   inputfile=sx->GetString(); 
435   if (inputfile=="") continue;
436
437   // Open the input file in the default ascii format (autodetection) for reading 
438   fInput=rdmc_mcopen(inputfile.Data(),"r",RDMC_DEFAULT_ASCII_F);
439
440   if (!fInput)
441   {
442    cout << " *IceF2k Exec* No input file found with name : " << inputfile.Data() << endl;
443    continue;
444   }
445
446   // Initialise the event structure 
447   rdmc_init_mevt(&fEvent);
448
449   // Read the file header information
450   rdmc_rarr(fInput,&fHeader);
451
452   // Fill the database with geometry, calib. etc... parameters
453   // for all the devices
454   FillOMdbase();
455
456   // Set the fit definitions according to the F2000 header info
457   SetFitdefs();
458
459   // Set the trigger definitions according to the F2000 header info
460   SetTrigdefs();
461  
462   while (!rdmc_revt(fInput,&fHeader,&fEvent))
463   {
464    if (fMaxevt>-1 && nevt>=fMaxevt) break;
465
466    // Reset the complete Event structure
467    evt->Reset();
468
469    evt->SetRunNumber(fEvent.nrun);
470    evt->SetEventNumber(fEvent.enr);
471    evt->SetMJD(fEvent.mjd,fEvent.secs,fEvent.nsecs);
472
473    // Take trigger offset into account which might have been
474    // introduced during the filtering process.
475    // For simulated data this will be treated separately in PutMcTracks().
476    fToffset=fEvent.t_offset;
477
478    PutTrigger();
479
480    PutMcTracks();
481
482    PutRecoTracks();
483
484    PutHits();
485
486    // Enter the IceF2k processor parameters into the event structure 
487    params.SetSignal(fToffset,1);
488    params.SetSignal(fMctoffset,2);
489    params.SetSignal(fMctracks,3);
490    evt->AddDevice(params);
491
492    // Invoke all available sub-tasks (if any)
493    CleanTasks();
494    ExecuteTasks(opt);
495
496    if (fPrintfreq)
497    {
498     if (!(nevt%fPrintfreq)) evt->HeaderData();
499    }
500
501    // Write the complete structure to the output Tree
502    if (otree) otree->Fill();
503
504    // Update event counter
505    nevt++;
506   }
507   if (fMaxevt>-1 && nevt>=fMaxevt) break;
508  }
509
510  // Flush possible memory resident data to the output file
511  if (fOutfile) fOutfile->Write();
512
513  // Remove the IceEvent object from the environment
514  // and delete it as well
515  if (evt)
516  {
517   RemoveObject(evt);
518   delete evt;
519  }
520 }
521 ///////////////////////////////////////////////////////////////////////////
522 void IceF2k::FillOMdbase()
523 {
524 // Fill the database with geometry, calib. etc... parameters 
525 // for all the devices.
526
527  if (fHeader.nch<=0)
528  {
529   if (fOmdb)
530   {
531    delete fOmdb;
532    fOmdb=0;
533   }
534   return;
535  }
536
537  Int_t adccal=fHeader.is_calib.adc;
538  Int_t tdccal=fHeader.is_calib.tdc;
539  Int_t totcal=fHeader.is_calib.tot;
540
541  TF1 fadccal("fadccal","(x-[1])*[0]");
542  TF1 fadcdecal("fadcdecal","(x/[0])+[1]");
543  fadccal.SetParName(0,"BETA-ADC");
544  fadccal.SetParName(1,"PED-ADC");
545  fadcdecal.SetParName(0,"BETA-ADC");
546  fadcdecal.SetParName(1,"PED-ADC");
547
548  TF1 ftdccal("ftdccal","(x*[0])-[1]-([0]-1.)*32767.-[2]/sqrt([3])");
549  TF1 ftdcdecal("ftdcdecal","(x+([0]-1.)*32767.+[1]+[2]/sqrt([3]))/[0]");
550  ftdccal.SetParName(0,"BETA-TDC");
551  ftdccal.SetParName(1,"T0");
552  ftdccal.SetParName(2,"ALPHA-TDC");
553  ftdccal.SetParName(3,"ADC-SLEW");
554  ftdcdecal.SetParName(0,"BETA-TDC");
555  ftdcdecal.SetParName(1,"T0");
556  ftdcdecal.SetParName(2,"ALPHA-TDC");
557  ftdcdecal.SetParName(3,"ADC-SLEW");
558
559  TF1 ftotcal("ftotcal","x*[0]");
560  TF1 ftotdecal("ftotdecal","x/[0]");
561  ftotcal.SetParName(0,"BETA-TOT");
562  ftotdecal.SetParName(0,"BETA-TOT");
563
564  if (fOmdb)
565  {
566   fOmdb->Reset();
567  }
568  else
569  {
570   fOmdb=new AliObjMatrix();
571   fOmdb->SetNameTitle("OMDBASE","The OM geometry, calib. etc... database");
572   fOmdb->SetOwner();
573  }
574
575  IceAOM* dev=0;
576  Double_t pos[3]={0,0,0};
577  for (Int_t i=0; i<fHeader.nch; i++)
578  {
579   dev=new IceAOM();
580   dev->SetUniqueID(i+1);
581   // Slots to hold the various (de)calibration functions  
582   dev->SetSlotName("ADC",1);
583   dev->SetSlotName("LE",2);
584   dev->SetSlotName("TOT",3);
585   // Slots to hold hardware parameters
586   dev->SetSlotName("TYPE",4);
587   dev->SetSlotName("ORIENT",5);
588   dev->SetSlotName("THRESH",6);
589   dev->SetSlotName("SENSIT",7);
590
591   pos[0]=fHeader.x[i];
592   pos[1]=fHeader.y[i];
593   pos[2]=fHeader.z[i];
594   dev->SetPosition(pos,"car");
595
596   fadccal.SetParameter(0,fHeader.cal[i].beta_a);
597   fadccal.SetParameter(1,fHeader.cal[i].ped);
598   fadcdecal.SetParameter(0,fHeader.cal[i].beta_a);
599   if (!fHeader.cal[i].beta_a) fadcdecal.SetParameter(0,1);
600   fadcdecal.SetParameter(1,fHeader.cal[i].ped);
601
602   ftdccal.SetParameter(0,fHeader.cal[i].beta_t);
603   ftdccal.SetParameter(1,fHeader.cal[i].t_0);
604   ftdccal.SetParameter(2,fHeader.cal[i].alpha_t);
605   ftdccal.SetParameter(3,1.e20);
606   ftdcdecal.SetParameter(0,fHeader.cal[i].beta_t);
607   if (!fHeader.cal[i].beta_t) ftdcdecal.SetParameter(0,1);
608   ftdcdecal.SetParameter(1,fHeader.cal[i].t_0);
609   ftdcdecal.SetParameter(2,fHeader.cal[i].alpha_t);
610   ftdcdecal.SetParameter(3,1.e20);
611
612   ftotcal.SetParameter(0,fHeader.cal[i].beta_tot);
613   ftotdecal.SetParameter(0,fHeader.cal[i].beta_tot);
614   if (!fHeader.cal[i].beta_tot) ftotdecal.SetParameter(0,1);
615
616   if (adccal)
617   {
618    dev->SetDecalFunction(&fadcdecal,1);
619   }
620   else
621   {
622    dev->SetCalFunction(&fadccal,1);
623   }
624
625   if (tdccal)
626   {
627    dev->SetDecalFunction(&ftdcdecal,2);
628   }
629   else
630   {
631    dev->SetCalFunction(&ftdccal,2);
632   }
633
634   if (totcal)
635   {
636    dev->SetDecalFunction(&ftotdecal,3);
637   }
638   else
639   {
640    dev->SetCalFunction(&ftotcal,3);
641   }
642
643   dev->SetSignal(fHeader.type[i],4);
644   dev->SetSignal((Float_t)fHeader.costh[i],5);
645   dev->SetSignal(fHeader.thresh[i],6);
646   dev->SetSignal(fHeader.sensit[i],7);
647
648   fOmdb->EnterObject(i+1,1,dev);
649  }
650 }
651 ///////////////////////////////////////////////////////////////////////////
652 void IceF2k::SetFitdefs()
653 {
654 // Obtain the names of the variables for each fit procedure from the
655 // f2000 header. Each different fit procedure is then stored as a separate
656 // "hit" of an AliDevice object and the various fit variables are stored
657 // as separate signal slots of the corresponding "hit".
658 // Via the GetFitdefs() memberfunction this AliDevice object can be
659 // retrieved and stored in the ROOT output file if wanted.
660 // The name of the object is FitDefinitions and the stored data can be
661 // inspected via the AliDevice::Data() memberfunction and looks as follows :
662 //
663 //  *AliDevice::Data* Id :0 Name : FitDefinitions
664 //    Position Vector in car coordinates : 0 0 0
665 //    Err. in car coordinates : 0 0 0
666 //  The following 8 hits are registered :
667 //  *AliSignal::Data* Id :0
668 //    Position Vector in car coordinates : 0 0 0
669 //    Err. in car coordinates : 0 0 0
670 //    Owned by device : AliDevice Name : FitDefinitions
671 //    Slot : 1 Signal value : 0 name : id
672 //    Slot : 2 Signal value : 0 name : rchi2
673 //    Slot : 3 Signal value : 0 name : prob
674 //    Slot : 4 Signal value : 0 name : sigth
675 //    Slot : 5 Signal value : 0 name : covmin
676 //    Slot : 6 Signal value : 0 name : covmax
677 //    Slot : 7 Signal value : 0 name : cutflag
678 //    Slot : 8 Signal value : 0 name : chi2
679 //  *AliSignal::Data* Id :1
680 //    Position Vector in car coordinates : 0 0 0
681 //    Err. in car coordinates : 0 0 0
682 //    Owned by device : AliDevice Name : FitDefinitions
683 //    Slot : 1 Signal value : 0 name : id
684 //    Slot : 2 Signal value : 0 name : rchi2
685 //    Slot : 3 Signal value : 0 name : prob
686 // etc....  
687 //
688 // This memberfunction is based on the original idea/code by Adam Bouchta.
689
690  if (fHeader.n_fit<=0)
691  {
692   if (fFitdefs)
693   {
694    delete fFitdefs;
695    fFitdefs=0;
696   }
697   return;
698  }
699
700  if (fFitdefs)
701  {
702   fFitdefs->Reset(1);
703  }
704  else
705  {
706   fFitdefs=new AliDevice();
707  }
708
709  fFitdefs->SetName("FitDefinitions");
710  fFitdefs->SetHitCopy (1);
711
712  AliSignal s;
713  s.Reset();
714
715  for (Int_t i=0; i<fHeader.n_fit; i++)
716  {
717   s.SetUniqueID(fHeader.def_fit[i].id);
718   s.SetName(TString(fHeader.def_fit[i].tag));
719
720   for (Int_t j=0; j<fHeader.def_fit[i].nwords; j++)
721   {
722    s.SetSlotName(TString(fHeader.def_fit[i].words[j]),j+1);
723    s.SetSignal(0,j+1);
724   }
725
726   fFitdefs->AddHit(s);
727   s.Reset(1);
728  }
729 }
730 ///////////////////////////////////////////////////////////////////////////
731 void IceF2k::SetTrigdefs()
732 {
733 // Obtain the names of the variables for each trigger procedure from the
734 // f2000 header. Each different trigger procedure is then stored as a separate
735 // "hit" of an AliDevice object and the various trigger variables are stored
736 // as separate signal slots of the corresponding "hit".
737 // Via the GetFitdefs() memberfunction this AliDevice object can be
738 // retrieved and stored in the ROOT output file if wanted.
739 // The name of the object is TrigDefinitions and the stored data can be
740 // inspected via the AliDevice::Data() memberfunction and looks as follows :
741 //
742 //  *AliDevice::Data* Id : 0 Name : TrigDefinitions
743 //    Position Vector in car (rad) coordinates : 0 0 0
744 //    Err. in car (rad) coordinates : 0 0 0
745 //  The following 9 hits are registered : 
746 //  *AliSignal::Data* Id : 1 Name : main
747 //    Position Vector in car (rad) coordinates : 0 0 0
748 //    Err. in car (rad) coordinates : 0 0 0
749 //    Owned by device : AliDevice Id : 0 Name : TrigDefinitions
750 //     Slot : 1 Signal value : 0 name : trig_pulse_le
751 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
752 //     Slot : 3 Signal value : 0 name : regi_flag
753 //   *AliSignal::Data* Id : 2 Name : amaa
754 //     Position Vector in car (rad) coordinates : 0 0 0
755 //     Err. in car (rad) coordinates : 0 0 0
756 //     Owned by device : AliDevice Id : 0 Name : TrigDefinitions
757 //     Slot : 1 Signal value : 0 name : trig_pulse_le
758 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
759 //     Slot : 3 Signal value : 0 name : regi_flag
760 //   *AliSignal::Data* Id : 3 Name : amab10
761 //     Position Vector in car (rad) coordinates : 0 0 0
762 //     Err. in car (rad) coordinates : 0 0 0
763 //     Owned by device : AliDevice Id : 0 Name : TrigDefinitions
764 //     Slot : 1 Signal value : 0 name : trig_pulse_le
765 //     Slot : 2 Signal value : 0 name : trig_pulse_tot
766 //     Slot : 3 Signal value : 0 name : regi_flag
767 // etc....  
768
769  if (fHeader.n_trigger<=0)
770  {
771   if (fTrigdefs)
772   {
773    delete fTrigdefs;
774    fTrigdefs=0;
775   }
776   return;
777  }
778
779  if (fTrigdefs)
780  {
781   fTrigdefs->Reset(1);
782  }
783  else
784  {
785   fTrigdefs=new AliDevice();
786  }
787
788  fTrigdefs->SetName("TrigDefinitions");
789  fTrigdefs->SetHitCopy (1);
790
791  AliSignal s;
792  s.Reset();
793
794  for (Int_t i=0; i<fHeader.n_trigger; i++)
795  {
796   s.SetUniqueID(fHeader.def_trig[i].id);
797   s.SetName(TString(fHeader.def_trig[i].tag));
798
799   for (Int_t j=0; j<fHeader.def_trig[i].nwords; j++)
800   {
801    s.SetSlotName(TString(fHeader.def_trig[i].words[j]),j+1);
802    s.SetSignal(0,j+1);
803   }
804
805   fTrigdefs->AddHit(s);
806   s.Reset(1);
807  }
808 }
809 ///////////////////////////////////////////////////////////////////////////
810 void IceF2k::PutMcTracks()
811 {
812 // Get the MC tracks from the F2000 file into the IcePack structure.
813 // Note : MC tracks are given negative track id's in the event structure.
814 // This memberfunction is based on the original code by Adam Bouchta.
815
816  IceEvent* evt=(IceEvent*)GetMainObject();
817  if (!evt || fEvent.ntrack<=0) return;
818
819  // User defined trigger offset in case of simulated data.
820  // The offset in the F2K file is meant to put the primary interaction
821  // well ahead of the detector trigger.
822  // See the introductory docs of this IceF2k class for further details.
823  fToffset=fMctoffset;
824
825  if (!fMctracks) return;
826
827  // Loop over all the tracks and add them to the current event
828  AliTrack t;
829  Double_t vec[3];
830  AliPosition r;
831  Ali3Vector p;
832  Int_t tid=0;
833  Int_t idpdg=0;
834  Int_t idf2k=0;
835  for (Int_t i=0; i<fEvent.ntrack; i++)
836  {
837   t.Reset ();
838
839   // Beginpoint of the track
840   vec[0]=fEvent.gen[i].x;
841   vec[1]=fEvent.gen[i].y;
842   vec[2]=fEvent.gen[i].z;
843   r.SetPosition(vec,"car");
844   t.SetBeginPoint(r);
845
846   // Endpoint of the track
847   vec[0]+=fEvent.gen[i].length*fEvent.gen[i].px;
848   vec[1]+=fEvent.gen[i].length*fEvent.gen[i].py;
849   vec[2]+=fEvent.gen[i].length*fEvent.gen[i].pz;
850   r.SetPosition(vec,"car");
851   t.SetEndPoint(r);
852
853   // Momentum in GeV/c
854   vec[0]=fEvent.gen[i].e*fEvent.gen[i].px*1e-3;
855   vec[1]=fEvent.gen[i].e*fEvent.gen[i].py*1e-3;
856   vec[2]=fEvent.gen[i].e*fEvent.gen[i].pz*1e-3;
857   p.SetVector (vec,"car");
858   t.Set3Momentum(p);
859
860   // MC tracks are indicated by negative track id's
861   tid=fEvent.gen[i].tag;
862   t.SetId(-abs(tid));
863
864   idf2k=fEvent.gen[i].id;
865   idpdg=0;
866   if (idf2k>1000)
867   {
868    idpdg=idf2k+10000000;
869   }
870   else if (idf2k <= 48)
871   {
872    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
873   }
874   else
875   {
876    if (idf2k==201) idpdg=12;
877    if (idf2k==202) idpdg=14;
878    if (idf2k==203) idpdg=16;
879    if (idf2k==204) idpdg=-12;
880    if (idf2k==205) idpdg=-14;
881    if (idf2k==206) idpdg=-16;
882   }
883
884   // Check for the user selected MC track storage
885   if (fMctracks==1) // Store only muon and muon-neutrino tracks
886   {
887    if (abs(idpdg)!=13 && abs(idpdg)!=14) continue;
888   }
889   else if (fMctracks==2) // Store all lepton tracks
890   {
891    if (abs(idpdg)<11 || abs(idpdg)>16) continue;
892   }
893
894   t.SetParticleCode(idpdg);
895   t.SetName(fPdg->GetParticle(idpdg)->GetName());
896   t.SetTitle("MC track");
897   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
898   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
899
900   evt->AddTrack(t);
901  }
902
903  // Create the pointers to each particle's parent particle.
904  Int_t txid=0;
905  Int_t parid=0;
906  for (Int_t itk=1; itk<=evt->GetNtracks (); itk++)
907  {
908   AliTrack* tx=evt->GetTrack(itk);
909
910   if (!tx) continue;
911
912   txid=tx->GetId();
913
914   parid=-1;
915   for (Int_t j=0; j<fEvent.ntrack; j++)
916   {
917    tid=fEvent.gen[j].tag;
918    if (-abs(tid) == txid) parid=fEvent.gen[j].parent;
919   }
920
921   if (parid<0) continue;
922
923   AliTrack* tpar=evt->GetIdTrack(-abs(parid));
924
925   if (!tpar) continue;
926
927   tx->SetParentTrack(tpar);
928  }
929 }
930 ///////////////////////////////////////////////////////////////////////////
931 void IceF2k::PutRecoTracks()
932 {
933 // Get the reconstructed tracks from the F2000 file into the IcePack structure.
934 // Note : Reco tracks are given positive track id's in the event structure.
935 // This memberfunction is based on the original code by Adam Bouchta.
936
937  IceEvent* evt=(IceEvent*)GetMainObject();
938  if (!evt || fEvent.nfit<=0) return;
939
940  // Loop over all the tracks and add them to the current event
941  AliTrack t;
942  Double_t vec[3];
943  AliPosition r;
944  Ali3Vector p;
945  Int_t tid=0;
946  Int_t idpdg=0;
947  Int_t idf2k=0;
948  for (Int_t i=0; i<fEvent.nfit; i++)
949  {
950   t.Reset ();
951
952   // Beginpoint of the track
953   vec[0]=fEvent.rec[i].x;
954   vec[1]=fEvent.rec[i].y;
955   vec[2]=fEvent.rec[i].z;
956   r.SetPosition(vec,"car");
957   t.SetBeginPoint(r);
958
959   // Endpoint of the track
960   vec[0]+=fEvent.rec[i].length*fEvent.rec[i].px;
961   vec[1]+=fEvent.rec[i].length*fEvent.rec[i].py;
962   vec[2]+=fEvent.rec[i].length*fEvent.rec[i].pz;
963   r.SetPosition(vec,"car");
964   t.SetEndPoint(r);
965
966   // Momentum in GeV/c
967   if (fEvent.rec[i].e > 0)
968   {
969    vec[0]=fEvent.rec[i].e*fEvent.rec[i].px*1e-3;
970    vec[1]=fEvent.rec[i].e*fEvent.rec[i].py*1e-3;
971    vec[2]=fEvent.rec[i].e*fEvent.rec[i].pz*1e-3;
972   }
973   else // Give the track a nominal momentum of 1 GeV/c
974   {
975    vec[0]=fEvent.rec[i].px;
976    vec[1]=fEvent.rec[i].py;
977    vec[2]=fEvent.rec[i].pz;
978   }
979   p.SetVector (vec,"car");
980   t.Set3Momentum(p);
981
982   // Use the fit number as track id
983   tid=fEvent.rec[i].tag;
984   t.SetId(abs(tid));
985
986   idf2k=fEvent.rec[i].id;
987   idpdg=0;
988   if (idf2k>1000)
989   {
990    idpdg=idf2k+10000000;
991   }
992   else if (idf2k <= 48)
993   {
994    idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
995   }
996   else
997   {
998    if (idf2k==201) idpdg=12;
999    if (idf2k==202) idpdg=14;
1000    if (idf2k==203) idpdg=16;
1001    if (idf2k==204) idpdg=-12;
1002    if (idf2k==205) idpdg=-14;
1003    if (idf2k==206) idpdg=-16;
1004   }
1005
1006   t.SetParticleCode(idpdg);
1007   t.SetNameTitle("Sieglinde","RECO track");
1008   t.SetMass(fPdg->GetParticle(idpdg)->Mass());
1009   t.SetCharge(fPdg->GetParticle(idpdg)->Charge()/3.);
1010
1011   // Retrieve the various fit parameters for this track
1012   AliSignal* fitdata=fFitdefs->GetIdHit(i);
1013   for (Int_t jval=0; jval<fEvent.fresult[i].nval; jval++)
1014   {
1015    fitdata->SetSignal(fEvent.fresult[i].val[jval],jval+1);
1016   }
1017
1018   // Store the various fit parameters for this track
1019   t.SetFitDetails(fitdata);
1020
1021   // Store the various reco tracks as track hypotheses.
1022   // A copy of the first reco track is entered as a new track instance
1023   // into the event and all reco tracks (incl. the first one) are
1024   // stored as hypotheses linked to this new reco track.
1025   if (i==0)
1026   {
1027    evt->AddTrack(t);
1028    AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1029    Int_t nrec=evt->GetNtracks(1);
1030    tx->SetId(nrec+1);
1031   }
1032   AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1033   if (tx) tx->AddTrackHypothesis(t);
1034  }
1035 }
1036 ///////////////////////////////////////////////////////////////////////////
1037 void IceF2k::PutHits()
1038 {
1039 // Get the hit and waveform info from the F2000 file into the IcePack structure.
1040 // This memberfunction is based on the original code by Adam Bouchta.
1041
1042  IceEvent* evt=(IceEvent*)GetMainObject();
1043  if (!evt) return;
1044
1045  // Loop over all the hits and add them to the current event
1046  IceAOM om;
1047  AliSignal s;
1048  s.SetSlotName("ADC",1);
1049  s.SetSlotName("LE",2);
1050  s.SetSlotName("TOT",3);
1051  Int_t chan=0;
1052  Int_t maxchan=800;
1053  if (fOmdb) maxchan=fHeader.nch;
1054  IceAOM* omx=0;
1055  AliSignal* sx=0;
1056  Int_t tid=0;
1057  AliTrack* tx=0;
1058  Float_t adc=0;
1059  Float_t adcfirst=0; // Adc value of the first hit of an OM
1060  for (Int_t i=0; i<fEvent.nhits; i++)
1061  {
1062   chan=fEvent.h[i].ch+1;
1063   if (chan>maxchan) continue; // Channels 9001, 9002 etc are trigger channels
1064
1065   // Get corresponding device from the current event structure  
1066   omx=(IceAOM*)evt->GetIdDevice(chan);
1067   if (!omx)
1068   {
1069    if (fOmdb)
1070    {
1071     omx=(IceAOM*)fOmdb->GetObject(chan,1);
1072     evt->AddDevice(omx);
1073    }
1074    else
1075    {
1076     om.Reset(1);
1077     om.SetUniqueID(chan);
1078     evt->AddDevice(om);
1079    }
1080    omx=(IceAOM*)evt->GetIdDevice(chan);
1081   }
1082
1083   if (!omx) continue;
1084
1085   adc=fEvent.h[i].amp;
1086
1087   // Multiple hits in the same OM with the same ADC value
1088   // are indicated by "*" in the F2K file.
1089   // This corresponds to a value of -2 in the data structure.
1090   if (int(adc) == -2)
1091   {
1092    adc=adcfirst;
1093   }
1094   else
1095   {
1096    adcfirst=adc;
1097   }
1098   s.Reset();
1099   s.SetUniqueID(fEvent.h[i].id);
1100   s.SetSignal(adc,1);
1101   s.SetSignal((fEvent.h[i].t-fToffset),2);
1102   s.SetSignal(fEvent.h[i].tot,3);
1103
1104   omx->AddHit(s);
1105
1106   sx=omx->GetHit(omx->GetNhits());
1107   if (!sx) continue;
1108
1109   // ADC dependent TDC (de)calibration function for this hit
1110   TF1* fcal=omx->GetCalFunction("LE");
1111   TF1* fdecal=omx->GetDecalFunction("LE");
1112   if (fcal) sx->SetCalFunction(fcal,2);
1113   if (fdecal) sx->SetDecalFunction(fdecal,2);
1114   fcal=sx->GetCalFunction(2);
1115   fdecal=sx->GetDecalFunction(2);
1116   adc=sx->GetSignal(1,-4);
1117   if (adc>0)
1118   {
1119    if (fcal) fcal->SetParameter(3,adc);
1120    if (fdecal) fdecal->SetParameter(3,adc);
1121   }
1122   else
1123   {
1124    if (fcal) fcal->SetParameter(3,1.e20);
1125    if (fdecal) fdecal->SetParameter(3,1.e20);
1126   }
1127
1128   // Bi-directional link between this hit and the track that caused the ADC value.
1129   // This F2K info is probably only present for MC tracks.
1130   tid=fEvent.h[i].ma;
1131   if (tid > 0)
1132   {
1133    tx=evt->GetIdTrack(tid); // Reco tracks
1134    if (!tx) tx=evt->GetIdTrack(-tid); // MC tracks
1135    if (tx) sx->AddTrack(*tx);
1136   }
1137   else
1138   {
1139    if (tid == -2) sx->SetNameTitle("N","Noise");
1140    if (tid == -3) sx->SetNameTitle("A","Afterpulse");
1141   }
1142  }
1143
1144  // Loop over all the waveforms and add the histo(s) to the corresponding OM's
1145  TH1F histo;
1146  Int_t nbins=0;
1147  Float_t xlow=0;
1148  Float_t xup=0;
1149  TString hname;
1150  for (Int_t iwf=0; iwf<fEvent.nwf; iwf++)
1151  {
1152   chan=fEvent.wf[iwf].om;
1153   if (chan<=0 || chan>maxchan) continue; // Skip trigger channels
1154
1155   // Get corresponding device from the current event structure  
1156   omx=(IceAOM*)evt->GetIdDevice(chan);
1157   if (!omx)
1158   {
1159    if (fOmdb)
1160    {
1161     omx=(IceAOM*)fOmdb->GetObject(chan,1);
1162     evt->AddDevice(omx);
1163    }
1164    else
1165    {
1166     om.Reset(1);
1167     om.SetUniqueID(chan);
1168     evt->AddDevice(om);
1169    }
1170    omx=(IceAOM*)evt->GetIdDevice(chan);
1171   }
1172
1173   if (!omx) continue;
1174
1175   hname="BASELINE-WF";
1176   hname+=omx->GetNwaveforms()+1;
1177   omx->AddNamedSlot(hname);
1178   omx->SetSignal(fEvent.wf[iwf].baseline,hname);
1179
1180   // Fill the waveform histogram
1181   hname="OM";
1182   hname+=chan;
1183   hname+="-WF";
1184   hname+=omx->GetNwaveforms()+1;
1185
1186   histo.Reset();
1187   histo.SetName(hname.Data());
1188   nbins=fEvent.wf[iwf].ndigi;
1189   xlow=fEvent.wf[iwf].t_start;
1190   xup=xlow+float(nbins)*fEvent.wf[iwf].t_bin;
1191   histo.SetBins(nbins,xlow,xup);
1192
1193   for (Int_t jbin=1; jbin<=fEvent.wf[iwf].ndigi; jbin++)
1194   {
1195    histo.SetBinContent(jbin,fEvent.wf[iwf].baseline-fEvent.wf[iwf].digi[jbin-1]);
1196   }
1197
1198   omx->SetWaveform(&histo,omx->GetNwaveforms()+1);
1199  }
1200
1201  // Set bi-directional links between hits and reco track hypotheses.
1202  // Note : Reco tracks are recognised by their positive id.
1203  Int_t hid=0;
1204  TObjArray* rectracks=evt->GetTracks(1);
1205  for (Int_t jtk=0; jtk<rectracks->GetEntries(); jtk++)
1206  {
1207   tx=(AliTrack*)rectracks->At(jtk);
1208   if (!tx) continue;
1209   
1210   for (Int_t jhyp=1; jhyp<=tx->GetNhypotheses(); jhyp++)
1211   {
1212    AliTrack* hypx=tx->GetTrackHypothesis(jhyp);
1213    if (!hypx) continue;
1214
1215    // Loop over all combinations of F2K fits and used OM hits
1216    for (Int_t k=0; k<fEvent.nfit_uses; k++)
1217    {
1218     if (fEvent.fit_uses[k].useid != hypx->GetId()) continue;
1219     hid=fEvent.fit_uses[k].hitid;
1220     sx=evt->GetIdHit(hid,"IceAOM");
1221     if (sx) sx->AddTrack(*hypx);
1222    }
1223   }
1224  }
1225 }
1226 ///////////////////////////////////////////////////////////////////////////
1227 void IceF2k::PutTrigger()
1228 {
1229 // Get the trigger info from the F2000 file into the IcePack structure.
1230
1231  if (!fTrigdefs) return;
1232
1233  IceEvent* evt=(IceEvent*)GetMainObject();
1234  if (!evt || fEvent.ntrig<=0) return;
1235
1236  AliDevice trig;
1237  trig.SetNameTitle("Trigger","Amanda/IceCube event triggers");
1238  AliSignal s;
1239  TString trigname;
1240  TString slotname;
1241  Int_t id=0;
1242  Int_t nval=0;
1243  for (Int_t i=0; i<fEvent.ntrig; i++)
1244  {
1245   id=fEvent.ptrig[i].id;
1246   nval=fEvent.ptrig[i].nval;
1247   if (!nval) continue;
1248   AliSignal* tdef=fTrigdefs->GetIdHit(id+1);
1249   if (!tdef) continue;
1250   trigname=tdef->GetName();
1251   s.Reset(1);
1252   s.SetName(trigname);
1253   s.SetUniqueID(id+1);
1254   for (Int_t jval=0; jval<fEvent.ptrig[i].nval; jval++)
1255   {
1256    slotname=tdef->GetSlotName(jval+1);
1257    s.SetSlotName(slotname,jval+1);
1258    s.SetSignal(fEvent.ptrig[i].val[jval],jval+1);
1259   }
1260   trig.AddHit(s);
1261  }
1262
1263  // Store the trigger data into the IceEvent structure
1264  evt->AddDevice(trig);
1265 }
1266 ///////////////////////////////////////////////////////////////////////////