1 /*******************************************************************************
2 * Copyright(c) 2003, IceCube Experiment at the South Pole. All rights reserved.
4 * Author: The IceCube RALICE-based Offline Project.
5 * Contributors are mentioned in the code where appropriate.
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 *******************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
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.
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
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
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.
67 // gSystem->Load("ralice");
68 // gSystem->Load("icepack");
69 // gSystem->Load("iceconvert");
71 // IceF2k q("IceF2k","F2K to IcePack data structure conversion");
73 // // Limit the number of entries for testing
74 // q.SetMaxEvents(10);
76 // // Print frequency to produce a short summary print every printfreq events
79 // // Split level for the output structures
80 // q.SetSplitLevel(0);
82 // // Buffer size for the output structures
83 // q.SetBufferSize(32000);
85 // // The F2K input filename(s)
86 // q.AddInputFile("run7825.f2k");
88 // // Output file for the event structures
89 // q.SetOutputFile("events.root");
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
98 // // MyXtalk task1("task1","Cross talk correction");
99 // // MyClean task2("task2","Hit cleaning");
103 // // The sub-tasks will be executed in the order as they are entered.
104 // ///////////////////////////////////////////////////////////////////
106 // // Perform the conversion and execute subtasks (if any)
107 // // on an event-by-event basis
110 // // Select various objects to be added to the output file
112 // TFile* ofile=q.GetOutputFile();
116 // ofile->cd(); // Switch to the output file directory
118 // AliObjMatrix* omdb=q.GetOMdbase();
119 // if (omdb) omdb->Write();
121 // AliDevice* fitdefs=q.GetFitdefs();
122 // if (fitdefs) fitdefs->Write();
124 // TDatabasePDG* pdg=q.GetPDG();
125 // if (pdg) pdg->Write();
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.
134 //--- Author: Nick van Eijndhoven 11-mar-2005 Utrecht University
135 //- Modified: NvE $Date$ Utrecht University
136 ///////////////////////////////////////////////////////////////////////////
139 #include "Riostream.h"
141 ClassImp(IceF2k) // Class implementation to enable ROOT I/O
143 IceF2k::IceF2k(const char* name,const char* title) : AliJob(name,title)
145 // Default constructor.
146 // By default maxevent=-1, split=0, bsize=32000, printfreq=1.
163 ///////////////////////////////////////////////////////////////////////////
166 // Default destructor.
198 ///////////////////////////////////////////////////////////////////////////
199 void IceF2k::SetMaxEvents(Int_t n)
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.
206 ///////////////////////////////////////////////////////////////////////////
207 void IceF2k::SetPrintFreq(Int_t f)
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;
213 ///////////////////////////////////////////////////////////////////////////
214 void IceF2k::SetSplitLevel(Int_t split)
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;
220 ///////////////////////////////////////////////////////////////////////////
221 void IceF2k::SetBufferSize(Int_t bsize)
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;
227 ///////////////////////////////////////////////////////////////////////////
228 void IceF2k::SetMcToffset(Float_t toffset)
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.
235 ///////////////////////////////////////////////////////////////////////////
236 void IceF2k::SelectMcTracks(Int_t mode)
238 // User selection of MC tracks to be stored in the event structure.
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
245 // By default mode=3 is set in the constructor of this class.
247 if (mode<0 || mode >3) return;
250 ///////////////////////////////////////////////////////////////////////////
251 void IceF2k::SetInputFile(TString name)
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().
261 if (fInfiles) delete fInfiles;
263 fInfiles=new TObjArray();
264 fInfiles->SetOwner();
266 TObjString* s=new TObjString();
270 ///////////////////////////////////////////////////////////////////////////
271 void IceF2k::AddInputFile(TString name)
273 // Add the name of this F2K input file to the list to be processed.
277 fInfiles=new TObjArray();
278 fInfiles->SetOwner();
281 TObjString* s=new TObjString();
285 ///////////////////////////////////////////////////////////////////////////
286 void IceF2k::SetOutputFile(TFile* ofile)
288 // Set the output file for the ROOT data.
289 if (fOutfile) delete fOutfile;
292 ///////////////////////////////////////////////////////////////////////////
293 void IceF2k::SetOutputFile(TString name)
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");
299 ///////////////////////////////////////////////////////////////////////////
300 TFile* IceF2k::GetOutputFile()
302 // Provide pointer to the ROOT output file.
305 ///////////////////////////////////////////////////////////////////////////
306 TDatabasePDG* IceF2k::GetPDG()
308 // Provide pointer to the PDG database
311 ///////////////////////////////////////////////////////////////////////////
312 AliObjMatrix* IceF2k::GetOMdbase()
314 // Provide pointer to the OM geometry, calib. etc... database
317 ///////////////////////////////////////////////////////////////////////////
318 AliDevice* IceF2k::GetFitdefs()
320 // Provide pointer to the fit definitions
323 ///////////////////////////////////////////////////////////////////////////
324 AliDevice* IceF2k::GetTrigdefs()
326 // Provide pointer to the trigger definitions
329 ///////////////////////////////////////////////////////////////////////////
330 void IceF2k::Exec(Option_t* opt)
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".
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.
351 cout << " *IceF2k Exec* No data input file(s) specified." << endl;
355 Int_t ninfiles=fInfiles->GetEntries();
358 cout << " *IceF2k Exec* No data input file(s) specified." << endl;
365 otree=new TTree("T","F2K Data converted to IceEvent structures");
366 otree->SetDirectory(fOutfile);
369 IceEvent* evt=new IceEvent();
370 evt->SetTrackCopy(1);
373 // Branch in the tree for the event structure
374 if (otree) otree->Branch("IceEvent","IceEvent",&evt,fBsize,fSplit);
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);
391 // Storage of the used parameters in the IceF2k device
393 params.SetNameTitle("IceF2k","IceF2k processor parameters");
394 params.SetSlotName("Toffset",1);
395 params.SetSlotName("Mctoffset",2);
396 params.SetSlotName("Mctracks",3);
398 // Initialise the job working environment
408 cout << " ***" << endl;
409 cout << " *** Start processing of job " << GetName() << " ***" << endl;
410 cout << " ***" << endl;
411 for (Int_t i=0; i<ninfiles; i++)
413 TObjString* sx=(TObjString*)fInfiles->At(i);
415 inputfile=sx->GetString();
416 cout << " F2K input file : " << inputfile.Data() << endl;
418 cout << " Maximum number of events to be processed : " << fMaxevt << endl;
419 cout << " Print frequency : " << fPrintfreq << endl;
422 cout << " ROOT output file : " << fOutfile->GetName() << endl;
423 cout << " Output characteristics : splitlevel = " << fSplit << " buffersize = " << fBsize << endl;
429 for (Int_t ifile=0; ifile<ninfiles; ifile++)
431 TObjString* sx=(TObjString*)fInfiles->At(ifile);
434 inputfile=sx->GetString();
435 if (inputfile=="") continue;
437 // Open the input file in the default ascii format (autodetection) for reading
438 fInput=rdmc_mcopen(inputfile.Data(),"r",RDMC_DEFAULT_ASCII_F);
442 cout << " *IceF2k Exec* No input file found with name : " << inputfile.Data() << endl;
446 // Initialise the event structure
447 rdmc_init_mevt(&fEvent);
449 // Read the file header information
450 rdmc_rarr(fInput,&fHeader);
452 // Fill the database with geometry, calib. etc... parameters
453 // for all the devices
456 // Set the fit definitions according to the F2000 header info
459 // Set the trigger definitions according to the F2000 header info
462 while (!rdmc_revt(fInput,&fHeader,&fEvent))
464 if (fMaxevt>-1 && nevt>=fMaxevt) break;
466 // Reset the complete Event structure
469 evt->SetRunNumber(fEvent.nrun);
470 evt->SetEventNumber(fEvent.enr);
471 evt->SetMJD(fEvent.mjd,fEvent.secs,fEvent.nsecs);
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;
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);
492 // Invoke all available sub-tasks (if any)
498 if (!(nevt%fPrintfreq)) evt->HeaderData();
501 // Write the complete structure to the output Tree
502 if (otree) otree->Fill();
504 // Update event counter
507 if (fMaxevt>-1 && nevt>=fMaxevt) break;
510 // Flush possible memory resident data to the output file
511 if (fOutfile) fOutfile->Write();
513 // Remove the IceEvent object from the environment
514 // and delete it as well
521 ///////////////////////////////////////////////////////////////////////////
522 void IceF2k::FillOMdbase()
524 // Fill the database with geometry, calib. etc... parameters
525 // for all the devices.
537 Int_t adccal=fHeader.is_calib.adc;
538 Int_t tdccal=fHeader.is_calib.tdc;
539 Int_t totcal=fHeader.is_calib.tot;
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");
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");
559 TF1 ftotcal("ftotcal","x*[0]");
560 TF1 ftotdecal("ftotdecal","x/[0]");
561 ftotcal.SetParName(0,"BETA-TOT");
562 ftotdecal.SetParName(0,"BETA-TOT");
570 fOmdb=new AliObjMatrix();
571 fOmdb->SetNameTitle("OMDBASE","The OM geometry, calib. etc... database");
576 Double_t pos[3]={0,0,0};
577 for (Int_t i=0; i<fHeader.nch; i++)
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);
594 dev->SetPosition(pos,"car");
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);
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);
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);
618 dev->SetDecalFunction(&fadcdecal,1);
622 dev->SetCalFunction(&fadccal,1);
627 dev->SetDecalFunction(&ftdcdecal,2);
631 dev->SetCalFunction(&ftdccal,2);
636 dev->SetDecalFunction(&ftotdecal,3);
640 dev->SetCalFunction(&ftotcal,3);
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);
648 fOmdb->EnterObject(i+1,1,dev);
651 ///////////////////////////////////////////////////////////////////////////
652 void IceF2k::SetFitdefs()
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 :
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
688 // This memberfunction is based on the original idea/code by Adam Bouchta.
690 if (fHeader.n_fit<=0)
706 fFitdefs=new AliDevice();
709 fFitdefs->SetName("FitDefinitions");
710 fFitdefs->SetHitCopy (1);
715 for (Int_t i=0; i<fHeader.n_fit; i++)
717 s.SetUniqueID(fHeader.def_fit[i].id);
718 s.SetName(TString(fHeader.def_fit[i].tag));
720 for (Int_t j=0; j<fHeader.def_fit[i].nwords; j++)
722 s.SetSlotName(TString(fHeader.def_fit[i].words[j]),j+1);
730 ///////////////////////////////////////////////////////////////////////////
731 void IceF2k::SetTrigdefs()
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 :
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
769 if (fHeader.n_trigger<=0)
785 fTrigdefs=new AliDevice();
788 fTrigdefs->SetName("TrigDefinitions");
789 fTrigdefs->SetHitCopy (1);
794 for (Int_t i=0; i<fHeader.n_trigger; i++)
796 s.SetUniqueID(fHeader.def_trig[i].id);
797 s.SetName(TString(fHeader.def_trig[i].tag));
799 for (Int_t j=0; j<fHeader.def_trig[i].nwords; j++)
801 s.SetSlotName(TString(fHeader.def_trig[i].words[j]),j+1);
805 fTrigdefs->AddHit(s);
809 ///////////////////////////////////////////////////////////////////////////
810 void IceF2k::PutMcTracks()
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.
816 IceEvent* evt=(IceEvent*)GetMainObject();
817 if (!evt || fEvent.ntrack<=0) return;
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.
825 if (!fMctracks) return;
827 // Loop over all the tracks and add them to the current event
835 for (Int_t i=0; i<fEvent.ntrack; i++)
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");
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");
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");
860 // MC tracks are indicated by negative track id's
861 tid=fEvent.gen[i].tag;
864 idf2k=fEvent.gen[i].id;
868 idpdg=idf2k+10000000;
870 else if (idf2k <= 48)
872 idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
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;
884 // Check for the user selected MC track storage
885 if (fMctracks==1) // Store only muon and muon-neutrino tracks
887 if (abs(idpdg)!=13 && abs(idpdg)!=14) continue;
889 else if (fMctracks==2) // Store all lepton tracks
891 if (abs(idpdg)<11 || abs(idpdg)>16) continue;
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.);
903 // Create the pointers to each particle's parent particle.
906 for (Int_t itk=1; itk<=evt->GetNtracks (); itk++)
908 AliTrack* tx=evt->GetTrack(itk);
915 for (Int_t j=0; j<fEvent.ntrack; j++)
917 tid=fEvent.gen[j].tag;
918 if (-abs(tid) == txid) parid=fEvent.gen[j].parent;
921 if (parid<0) continue;
923 AliTrack* tpar=evt->GetIdTrack(-abs(parid));
927 tx->SetParentTrack(tpar);
930 ///////////////////////////////////////////////////////////////////////////
931 void IceF2k::PutRecoTracks()
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.
937 IceEvent* evt=(IceEvent*)GetMainObject();
938 if (!evt || fEvent.nfit<=0) return;
940 // Loop over all the tracks and add them to the current event
948 for (Int_t i=0; i<fEvent.nfit; i++)
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");
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");
967 if (fEvent.rec[i].e > 0)
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;
973 else // Give the track a nominal momentum of 1 GeV/c
975 vec[0]=fEvent.rec[i].px;
976 vec[1]=fEvent.rec[i].py;
977 vec[2]=fEvent.rec[i].pz;
979 p.SetVector (vec,"car");
982 // Use the fit number as track id
983 tid=fEvent.rec[i].tag;
986 idf2k=fEvent.rec[i].id;
990 idpdg=idf2k+10000000;
992 else if (idf2k <= 48)
994 idpdg=fPdg->ConvertGeant3ToPdg(idf2k);
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;
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.);
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++)
1015 fitdata->SetSignal(fEvent.fresult[i].val[jval],jval+1);
1018 // Store the various fit parameters for this track
1019 t.SetFitDetails(fitdata);
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.
1028 AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1029 Int_t nrec=evt->GetNtracks(1);
1032 AliTrack* tx=evt->GetTrack(evt->GetNtracks());
1033 if (tx) tx->AddTrackHypothesis(t);
1036 ///////////////////////////////////////////////////////////////////////////
1037 void IceF2k::PutHits()
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.
1042 IceEvent* evt=(IceEvent*)GetMainObject();
1045 // Loop over all the hits and add them to the current event
1048 s.SetSlotName("ADC",1);
1049 s.SetSlotName("LE",2);
1050 s.SetSlotName("TOT",3);
1053 if (fOmdb) maxchan=fHeader.nch;
1059 Float_t adcfirst=0; // Adc value of the first hit of an OM
1060 for (Int_t i=0; i<fEvent.nhits; i++)
1062 chan=fEvent.h[i].ch+1;
1063 if (chan>maxchan) continue; // Channels 9001, 9002 etc are trigger channels
1065 // Get corresponding device from the current event structure
1066 omx=(IceAOM*)evt->GetIdDevice(chan);
1071 omx=(IceAOM*)fOmdb->GetObject(chan,1);
1072 evt->AddDevice(omx);
1077 om.SetUniqueID(chan);
1080 omx=(IceAOM*)evt->GetIdDevice(chan);
1085 adc=fEvent.h[i].amp;
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.
1099 s.SetUniqueID(fEvent.h[i].id);
1101 s.SetSignal((fEvent.h[i].t-fToffset),2);
1102 s.SetSignal(fEvent.h[i].tot,3);
1106 sx=omx->GetHit(omx->GetNhits());
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);
1119 if (fcal) fcal->SetParameter(3,adc);
1120 if (fdecal) fdecal->SetParameter(3,adc);
1124 if (fcal) fcal->SetParameter(3,1.e20);
1125 if (fdecal) fdecal->SetParameter(3,1.e20);
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.
1133 tx=evt->GetIdTrack(tid); // Reco tracks
1134 if (!tx) tx=evt->GetIdTrack(-tid); // MC tracks
1135 if (tx) sx->AddTrack(*tx);
1139 if (tid == -2) sx->SetNameTitle("N","Noise");
1140 if (tid == -3) sx->SetNameTitle("A","Afterpulse");
1144 // Loop over all the waveforms and add the histo(s) to the corresponding OM's
1150 for (Int_t iwf=0; iwf<fEvent.nwf; iwf++)
1152 chan=fEvent.wf[iwf].om;
1153 if (chan<=0 || chan>maxchan) continue; // Skip trigger channels
1155 // Get corresponding device from the current event structure
1156 omx=(IceAOM*)evt->GetIdDevice(chan);
1161 omx=(IceAOM*)fOmdb->GetObject(chan,1);
1162 evt->AddDevice(omx);
1167 om.SetUniqueID(chan);
1170 omx=(IceAOM*)evt->GetIdDevice(chan);
1175 hname="BASELINE-WF";
1176 hname+=omx->GetNwaveforms()+1;
1177 omx->AddNamedSlot(hname);
1178 omx->SetSignal(fEvent.wf[iwf].baseline,hname);
1180 // Fill the waveform histogram
1184 hname+=omx->GetNwaveforms()+1;
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);
1193 for (Int_t jbin=1; jbin<=fEvent.wf[iwf].ndigi; jbin++)
1195 histo.SetBinContent(jbin,fEvent.wf[iwf].baseline-fEvent.wf[iwf].digi[jbin-1]);
1198 omx->SetWaveform(&histo,omx->GetNwaveforms()+1);
1201 // Set bi-directional links between hits and reco track hypotheses.
1202 // Note : Reco tracks are recognised by their positive id.
1204 TObjArray* rectracks=evt->GetTracks(1);
1205 for (Int_t jtk=0; jtk<rectracks->GetEntries(); jtk++)
1207 tx=(AliTrack*)rectracks->At(jtk);
1210 for (Int_t jhyp=1; jhyp<=tx->GetNhypotheses(); jhyp++)
1212 AliTrack* hypx=tx->GetTrackHypothesis(jhyp);
1213 if (!hypx) continue;
1215 // Loop over all combinations of F2K fits and used OM hits
1216 for (Int_t k=0; k<fEvent.nfit_uses; k++)
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);
1226 ///////////////////////////////////////////////////////////////////////////
1227 void IceF2k::PutTrigger()
1229 // Get the trigger info from the F2000 file into the IcePack structure.
1231 if (!fTrigdefs) return;
1233 IceEvent* evt=(IceEvent*)GetMainObject();
1234 if (!evt || fEvent.ntrig<=0) return;
1237 trig.SetNameTitle("Trigger","Amanda/IceCube event triggers");
1243 for (Int_t i=0; i<fEvent.ntrig; i++)
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();
1252 s.SetName(trigname);
1253 s.SetUniqueID(id+1);
1254 for (Int_t jval=0; jval<fEvent.ptrig[i].nval; jval++)
1256 slotname=tdef->GetSlotName(jval+1);
1257 s.SetSlotName(slotname,jval+1);
1258 s.SetSignal(fEvent.ptrig[i].val[jval],jval+1);
1263 // Store the trigger data into the IceEvent structure
1264 evt->AddDevice(trig);
1266 ///////////////////////////////////////////////////////////////////////////