NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
RnoConvert.cxx
Go to the documentation of this file.
1
17
19
44
45#include "RnoConvert.h"
46#include "Riostream.h"
47
48ClassImp(RnoConvert); // Class implementation to enable ROOT I/O
49
51RnoConvert::RnoConvert(const char* name,const char* title) : NcJob(name,title)
52{
59
60 fSplit=99;
61 fBsize=32000;
62 fMaxevt=-1;
63 fPrintfreq=0;
65 fOutfile=0;
66 fHdr=0;
67 fDs=0;
68 fWf=0;
69 fComb=0;
70 fPed=0;
71 fData=0;
74}
75
77{
83
84 if (fData)
85 {
86 delete fData;
87 fData=0;
88 }
89
90 if (fOutfile)
91 {
92 delete fOutfile;
93 fOutfile=0;
94 }
95}
96
98{
106
107 fMaxevt=n;
108}
109
110void RnoConvert::SetPrintFreq(Int_t m,Int_t level)
111{
127
128 if (level<0 || level>2)
129 {
130 cout << " *" << ClassName() << "::SetPrintFreq* Unsupported input : level=" << level << " --> Print level set to 0." << endl;
131 cout << endl;
132 level=0;
133 }
134
135 fPrintfreq=m;
136 fPrintlevel=level;
137}
138
140{
147
148 if (split>=0) fSplit=split;
149}
150
152{
159
160 if (bsize>=0) fBsize=bsize;
161}
162
163void RnoConvert::AddInputFile(TString file,TString type)
164{
184
185 // Expand the path name of the provided input file
186 file=gSystem->ExpandPathName(file.Data());
187
188 if (type!="hdr" && type!="ds" && type!="wf" && type!="combined" && type!="pedestal")
189 {
190 printf(" *%-s::AddInputFile* Unsupported type : %-s. File %-s not added. \n",ClassName(),type.Data(),file.Data());
191 return;
192 }
193
194 TChain* chain=0;
195 if (type=="hdr") chain=fHdr;
196 if (type=="ds") chain=fDs;
197 if (type=="wf") chain=fWf;
198 if (type=="combined") chain=fComb;
199 if (type=="pedestal") chain=fPed;
200
201 if (!chain)
202 {
203 if (type!="pedestal")
204 {
205 chain=new TChain(type.Data());
206 }
207 else
208 {
209 chain=new TChain();
210 }
211 if (type=="hdr") fHdr=chain;
212 if (type=="ds") fDs=chain;
213 if (type=="wf") fWf=chain;
214 if (type=="combined") fComb=chain;
215 if (type=="pedestal") fPed=chain;
216 }
217
218 chain->Add(file.Data());
219
220 printf(" *%-s::AddInputFile* Added RNO-G %-s data input file : %-s \n",ClassName(),type.Data(),file.Data());
221}
222
224{
230
231 if (fOutfile) delete fOutfile;
232 fOutfile=ofile;
233}
234
236{
245
246 // Expand the path name of the specified output file
247 name=gSystem->ExpandPathName(name.Data());
248
249 if (fOutfile) delete fOutfile;
250 fOutfile=new TFile(name.Data(),"RECREATE","RNO-G data in RnoEvent structure");
251}
252
253void RnoConvert::SetSelectLevels(Int_t min,Int_t max)
254{
265
266 fMinSelectLevel=min;
267 fMaxSelectLevel=max;
268}
269
271{
277
278 return fMinSelectLevel;
279}
280
282{
288
289 return fMaxSelectLevel;
290}
291
292void RnoConvert::ListInput(TString type,Option_t* opt)
293{
312
313 if (!fHdr && !fDs && !fWf && !fComb && !fPed)
314 {
315 printf("\n *%-s::ListInput* No input file of type %-s has been attached. \n",ClassName(),type.Data());
316 return;
317 }
318
319 TString s=opt;
320 if (s=="") s="Default";
321
322 // Printout for all stored data
323 if (type=="*")
324 {
325 printf("\n *%-s::ListInput* Overview of all the stored input data with option : %-s \n",ClassName(),s.Data());
326 if (fHdr) fHdr->Print(opt);
327 if (fDs) fDs->Print(opt);
328 if (fWf) fWf->Print(opt);
329 if (fComb) fComb->Print(opt);
330 if (fPed) fPed->Print(opt);
331 return;
332 }
333
334 // Selected printout
335 TChain* chain=0;
336 if (type=="hdr") chain=fHdr;
337 if (type=="ds") chain=fDs;
338 if (type=="wf") chain=fWf;
339 if (type=="combined") chain=fComb;
340 if (type=="pedestal") chain=fPed;
341
342 if (chain)
343 {
344 printf("\n *%-s::ListInput* Overview of the %-s input data with option : %-s \n",ClassName(),type.Data(),s.Data());
345 chain->Print(opt);
346 }
347 else
348 {
349 printf("\n *%-s::ListInput* No input file of type %-s has been attached. \n",ClassName(),type.Data());
350 return;
351 }
352}
353
355{
361
362 return fOutfile;
363}
364
366{
372
373 if (fData)
374 {
375 delete fData;
376 fData=0;
377 }
378
379 // Check whether loaded input data exists
380 if (!fHdr && !fDs && !fWf && !fComb && !fPed) return;
381
382 fData=new TChain("Data");
383
384 // Build the total input chain
385 Bool_t first=kTRUE;
386
387 if (fHdr)
388 {
389 if (first)
390 {
391 fData=new TChain("Data");
392 fData->Add(fHdr);
393 first=kFALSE;
394 }
395 else
396 {
397 if (fData) fData->AddFriend(fHdr);
398 }
399 }
400
401 if (fDs)
402 {
403 if (first)
404 {
405 fData=new TChain("Data");
406 fData->Add(fDs);
407 first=kFALSE;
408 }
409 else
410 {
411 if (fData) fData->AddFriend(fDs);
412 }
413 }
414
415 if (fWf)
416 {
417 if (first)
418 {
419 fData=new TChain("Data");
420 fData->Add(fWf);
421 first=kFALSE;
422 }
423 else
424 {
425 if (fData) fData->AddFriend(fWf);
426 }
427 }
428
429 if (fComb)
430 {
431 if (first)
432 {
433 fData=new TChain("Data");
434 fData->Add(fComb);
435 first=kFALSE;
436 }
437 else
438 {
439 if (fData) fData->AddFriend(fComb);
440 }
441 }
442
443 if (fPed)
444 {
445 if (first)
446 {
447 fData=new TChain("Data");
448 fData->Add(fPed);
449 first=kFALSE;
450 }
451 else
452 {
453 if (fData) fData->AddFriend(fPed);
454 }
455 }
456}
457
458void RnoConvert::Exec(Option_t* opt)
459{
480
481 // Create the main input data chain for processing
483
484 if (!fData)
485 {
486 cout << " *" << ClassName() << "::Exec* No data input file(s) specified." << endl;
487 return;
488 }
489
490 // Create output tree if necessary
491 TTree* otree=0;
492 if (fOutfile)
493 {
494 otree=new TTree("T","RNO-G data converted to RnoEvent structures");
495 otree->SetDirectory(fOutfile);
496 }
497
498 // An initial RNO-G detector structure
499 RnoDetector det;
500
501 // Create the RnoEvent structure
502 RnoEvent* evt=new RnoEvent();
503 evt->SetOwner(1);
504
505 // Branch in the tree for the event structure
506 if (otree) otree->Branch("Events","RnoEvent",&evt,fBsize,fSplit);
507
508 // Initialise the job working environment
509 SetMainObject(evt);
510
511 // Some output for the user's convenience
512 cout << endl;
513 cout << " *" << ClassName() << "::Exec* Overview of scheduled processing." << endl;
514 cout << " ***" << endl;
515 cout << " *** Start processing of job " << GetName() << " ***" << endl;
516 cout << " ***" << endl;
517 cout << " Maximum number of events to be processed : " << fMaxevt << " (<0 means all events)" << endl;
518 cout << " Print frequency : " << fPrintfreq << " (<=0 means no printout)" << endl;
519 cout << " Print level : " << fPrintlevel << endl;
520 if (fOutfile)
521 {
522 cout << " RnoEvent output file : " << fOutfile->GetName() << endl;
523 cout << " Output characteristics : splitlevel = " << fSplit << " buffersize = " << fBsize << endl;
524 cout << " Required event selection level interval for output : [" << fMinSelectLevel << "," << fMaxSelectLevel << "]" << endl;
525 }
526 cout << endl;
527
529
530 // The number of entries in the input chain
531 Int_t nen=fData->GetEntries();
532
533 // Check for the maximum number of events to be processed
534 if (fMaxevt>-1 && nen>fMaxevt) nen=fMaxevt;
535
537 // Variables in the RNO-G data tree //
539
540 // Header info
541 Int_t run=0;
542 Int_t event=0;
543 Int_t station=0;
544 Double_t trigtime=0;
545
546 // DAQ info
547 NcDevice daq("DAQ","DAQ status");
548 daq.AddNamedSlot("RADIANT-update");
549 daq.AddNamedSlot("FLOWER-update");
550 daq.AddNamedSlot("Sampling-rate");
551
552 // Trigger info
553 NcTagger trigger("Trigger","Trigger tags");
554
555 // The waveforms of all channels are provided as Short_t radiant_data[24][2048]
556 // but will be readout linearly and stored as an NcSample for each channel, as indicated below.
557 NcSample signal("Signals","Radiant signals");
558 signal.SetStoreMode(1);
559
560 // The pedestals of all channels are provided as UShort_t pedestals[24][4096]
561 // but will be readout linearly and stored as an NcSample for each channel, as indicated below.
562 NcSample pedestal("Pedestals","Pedestal values");
563 pedestal.SetStoreMode(1);
564
566 // Loop over the entries in the input data chain //
568
569 TLeaf* lx=0;
570 TLeaf* lradiant=0; // Pointer to the Radiant waveform data
571 TLeaf* lpedestal=0; // Pointer to the Pedestal waveform data
572 Int_t idx=0;
573 Int_t nsamples=0; // The sampling buffer length
574 NcDevice* devx=0;
575 Int_t idsample=0; // Index for the NcSample in the NcDevice storage
576 TString name="";
577 TString namex1="";
578 TString namex2="";
579 TString namex3="";
580 Float_t value=0;
581 Int_t ivalue=0;
582 Bool_t flag=kFALSE;
583 Int_t iradiant=0; // RADIANT firmware update sequence number
584 Int_t iflower=0; // FLOWER firmware update sequence number
585 Float_t fsample; // The DAQ sampling rate in Hz
586 Int_t nwritten=0; // The number of events written to output
587 Int_t perc=0; // Percentage of the number of processed events
588
589 for (Int_t ient=0; ient<nen; ient++)
590 {
591 // Reset the detector structure
592 det.Reset();
593 idsample=0;
594
595 fData->GetEntry(ient);
596
597 // Loop over all the leaves and extract the relevant data for this entry.
598 // This approach makes the functionality independent of the Tree/Branch structure.
599 run=0;
600 event=0;
601 station=0;
602 trigtime=0;
603 nsamples=0;
604 trigger.Reset();
605 lradiant=0;
606 lpedestal=0;
607 TIterator* iter=fData->GetIteratorOnAllLeaves();
608 while((lx=(TLeaf*)iter->Next()))
609 {
610 name=lx->GetName();
611
612 // Header data
613 if (!run && name=="run_number") run=lx->GetValue();
614 if (!event && name=="event_number") event=lx->GetValue();
615 if (!station && name=="station_number") station=lx->GetValue();
616 if (!trigtime && name=="trigger_time") trigtime=lx->GetValue();
617 if (!trigtime && name=="when") trigtime=lx->GetValue(); // In case of a pedestal data file
618 if (!nsamples && name=="buffer_length") nsamples=lx->GetValue();
619
620 // Trigger data
621 if (name.Contains("trigger_info"))
622 {
623 name.ReplaceAll("trigger_info.","");
624 value=lx->GetValue();
625 ivalue=TMath::Nint(value);
626 namex1="";
627 namex2="";
628 namex3="";
629
630 // The low threshold time window and number of coincidences
631 if (name.Contains("lt_info.window")) namex1="lt-window";
632 if (name.Contains("lt_info.num_coinc")) namex1="lt-ncoinc";
633
634 if (namex1!="")
635 {
636 trigger.AddNamedSlot(namex1);
637 trigger.SetSignal(value,namex1);
638 continue;
639 }
640
641 // The various radiant (=surface) time windows and number of coincidences
642 if (name.Contains("radiant_info.RF_window"))
643 {
644 namex1="LPDA-up-window";
645 namex2="LPDA-down-window";
646 value=lx->GetValue(0);
647 trigger.AddNamedSlot(namex1);
648 trigger.SetSignal(value,namex1);
649 value=lx->GetValue(1);
650 trigger.AddNamedSlot(namex2);
651 trigger.SetSignal(value,namex2);
652 continue;
653 }
654 if (name.Contains("radiant_info.RF_ncoinc"))
655 {
656 namex1="LPDA-up-ncoinc";
657 namex2="LPDA-down-ncoinc";
658 value=lx->GetValue(0);
659 trigger.AddNamedSlot(namex1);
660 trigger.SetSignal(value,namex1);
661 value=lx->GetValue(1);
662 trigger.AddNamedSlot(namex2);
663 trigger.SetSignal(value,namex2);
664 continue;
665 }
666
667 // The actual trigger tags
668 if (name.Contains("_trigger"))
669 {
670 // Settings of the various radiant (=surface) triggers
671 if (name.Contains("which_radiant"))
672 {
673 namex1="LPDA-up_trigger";
674 namex2="LPDA-down_trigger";
675 namex3="radiant-unknown_trigger";
676 if (ivalue<-100) // None of the surface triggers
677 {
678 trigger.SetPass(namex1,kFALSE);
679 trigger.SetPass(namex2,kFALSE);
680 trigger.SetPass(namex3,kFALSE);
681 }
682 if (ivalue==-1) // Unknown radiant trigger
683 {
684 trigger.SetPass(namex1,kFALSE);
685 trigger.SetPass(namex2,kFALSE);
686 trigger.SetPass(namex3,kTRUE);
687 }
688 if (ivalue==0) // Only upward surface trigger
689 {
690 trigger.SetPass(namex1,kTRUE);
691 trigger.SetPass(namex2,kFALSE);
692 trigger.SetPass(namex3,kFALSE);
693 }
694 if (ivalue==1) // Only downward surface trigger
695 {
696 trigger.SetPass(namex1,kFALSE);
697 trigger.SetPass(namex2,kTRUE);
698 trigger.SetPass(namex3,kFALSE);
699 }
700 }
701 else
702 {
703 flag=kFALSE;
704 if (ivalue) flag=kTRUE;
705 trigger.SetPass(name,flag);
706 }
707 }
708 } // End of trigger data
709
710 // Pointers to the waveform data
711 if (name=="radiant_data") lradiant=lx;
712 if (name=="pedestals") lpedestal=lx;
713 } // End of loop over the leaves
714
715 // Create this station in the detector structure
716 RnoStation* stax=det.GetStation(station,kTRUE);
717
718 if (!stax) break;
719
720 // DAQ info
721 iradiant=2; // The RADIANT update version
722 iflower=0; // The FLOWER update version
723 fsample=3.2e9; // The sampling rate in Hz
724 if (station==11)
725 {
726 if (run<571) iradiant=1;
727 if (run<474) iradiant=0;
728 if (run>=571) iradiant=2;
729 }
730 if (station==21)
731 {
732 if (run<753) iradiant=1;
733 if (run<646) iradiant=0;
734 }
735 if (station==22)
736 {
737 if (run<656) iradiant=1;
738 if (run<574) iradiant=0;
739 }
740
741 daq.SetSignal(iradiant,"RADIANT-update");
742 daq.SetSignal(iflower,"FLOWER-update");
743 daq.SetSignal(fsample,"Sampling-rate");
744
745 stax->AddDevice(daq);
746
747 stax->AddDevice(trigger);
748
749 // Readout the signal waveforms of all Radiant channels
750 if (!nsamples) nsamples=2048;
751 lx=lradiant;
752 if (lx)
753 {
754 idsample++;
755 idx=0;
756 for (Int_t i=0; i<24; i++)
757 {
758 // Access the corresponding channel of this station
759 name="Ch";
760 name+=i;
761 devx=stax->GetDevice(name,kTRUE);
762
763 if (!devx) continue;
764
765 signal.Reset();
766 signal.SetNames("ADC");
767 for (Int_t j=0; j<nsamples; j++)
768 {
769 // Retrieve the value of radiant[i][j]
770 value=lx->GetValue(idx);
771 signal.Enter(value);
772 idx++;
773 }
774 devx->SetSample(&signal,idsample);
775 }
776 }
777
778 // Readout the pedestals of all Radiant channels
779 nsamples=4096;
780 lx=lpedestal;
781 if (lx)
782 {
783 idsample++;
784 idx=0;
785 for (Int_t i=0; i<24; i++)
786 {
787 // Access the corresponding channel of this station
788 name="Ch";
789 name+=i;
790 devx=stax->GetDevice(name,kTRUE);
791
792 if (!devx) continue;
793
794 pedestal.Reset();
795 pedestal.SetNames("ADC");
796 for (Int_t j=0; j<nsamples; j++)
797 {
798 // Retrieve the value of pedestals[i][j]
799 value=lx->GetValue(idx);
800 pedestal.Enter(value);
801 idx++;
802 }
803 devx->SetSample(&pedestal,idsample);
804 }
805 }
806
807 // Transfer the RNO-G data into the RnoEvent structure
808 evt->Reset();
809 evt->SetUnixTime(trigtime,"A");
810 evt->SetRunNumber(run);
811 evt->SetEventNumber(event);
812 evt->SetDetector(det);
813 UInt_t id=station+100*run;
814 evt->SetUniqueID(id);
815
816 // Invoke all available sub-tasks (if any)
817 CleanTasks();
818 ExecuteTasks(GetName());
819
820 // Provide a printout every "Printfreq" events
821 if (fPrintfreq>0)
822 {
823 if ((!(ient%fPrintfreq) || (ient+1)==nen))
824 {
825 perc=100*(ient+1)/nen;
826 if (fPrintlevel==0 || fPrintlevel==2)
827 {
828 cout << " *** Processed input entry : " << ient << " run : " << run << " event : " << event << " (" << perc << "%)" << endl;
829 }
830 if (fPrintlevel>0)
831 {
832 evt->HeaderData();
833 cout << endl;
834 }
835 }
836 }
837
838 // Write the event to the output file (if the event select level is o.k.)
839 Int_t select=evt->GetSelectLevel();
840 if (select<fMinSelectLevel) continue;
841 if (fMaxSelectLevel>=fMinSelectLevel && select>fMaxSelectLevel) continue;
842
843 if (otree) otree->Fill();
844 nwritten++;
845 } // End of loop over the entries
846
847 // Flush possible memory resident data to the output file
848 if (fOutfile)
849 {
850 fOutfile->Write();
851 cout << endl;
852 cout << " *" << ClassName() << "::Exec* Number of (selected) events written to output : " << nwritten << endl;
853 }
854
855 // Remove the RnoEvent object from the environment
856 // and delete it as well
857 if (evt)
858 {
859 RemoveObject(evt);
860 delete evt;
861 }
862
863 // Delete the input data chain
864 if (fData)
865 {
866 delete fData;
867 fData=0;
868 }
869}
870
ClassImp(RnoConvert)
void AddNamedSlot(TString s)
NcDevice * GetDevice(Int_t i) const
void AddDevice(NcDevice &d)
virtual void Reset(Int_t mode=0)
Signal (Hit) handling of a generic device.
Definition NcDevice.h:14
virtual void Reset(Int_t mode=0)
Definition NcDevice.cxx:222
void SetEventNumber(Int_t evt)
Definition NcEvent.cxx:573
virtual void HeaderData()
Definition NcEvent.cxx:996
Int_t GetSelectLevel() const
Definition NcEvent.cxx:659
void SetDetector(NcDetector d)
Definition NcEvent.cxx:1036
virtual void SetOwner(Bool_t own=kTRUE)
Definition NcEvent.cxx:490
void SetRunNumber(Int_t run)
Definition NcEvent.cxx:562
void SetMainObject(TObject *obj)
Definition NcJob.cxx:305
void RemoveObject(TObject *obj)
Definition NcJob.cxx:384
NcJob(const char *name="NcJob", const char *title="")
Definition NcJob.cxx:139
void ListEnvironment()
Definition NcJob.cxx:205
Sampling and statistics tools for various multi-dimensional data samples.
Definition NcSample.h:28
void Reset()
Definition NcSample.cxx:255
void Enter(Double_t x)
Definition NcSample.cxx:357
void SetNames(TString name1="X", TString name2="Y", TString name3="Z", TString name4="T")
Definition NcSample.cxx:327
void SetStoreMode(Int_t mode=1, Int_t nmax=0, Int_t i=0)
virtual void SetSignal(Double_t sig, Int_t j=1)
Definition NcSignal.cxx:516
void SetSample(NcSample *sample, Int_t j=1)
Handling of generic event classification tags.
Definition NcTagger.h:12
void SetPass(TString name, Bool_t flag)
Definition NcTagger.cxx:168
Int_t SetUnixTime(Double_t sec, TString utc="A", Int_t leap=0, Double_t dut=0)
Job for conversion of RNO-G Root data into RnoEvent data structures.
Definition RnoConvert.h:23
TChain * fPed
Definition RnoConvert.h:52
TFile * GetOutputFile()
Int_t GetMaxSelectLevel() const
virtual ~RnoConvert()
Int_t fMinSelectLevel
Definition RnoConvert.h:54
Int_t GetMinSelectLevel() const
void AddInputFile(TString file, TString type)
Int_t fMaxSelectLevel
Definition RnoConvert.h:55
void SetPrintFreq(Int_t m, Int_t level=0)
Int_t fSplit
Definition RnoConvert.h:42
void SetSelectLevels(Int_t min, Int_t max)
TChain * fComb
Definition RnoConvert.h:51
TFile * fOutfile
Definition RnoConvert.h:47
Int_t fPrintfreq
Definition RnoConvert.h:45
Int_t fPrintlevel
Definition RnoConvert.h:46
TChain * fWf
Definition RnoConvert.h:50
RnoConvert(const char *name="RnoConvert", const char *title="")
void ListInput(TString type="*", Option_t *opt="")
TChain * fData
Definition RnoConvert.h:53
virtual void Exec(Option_t *opt)
void SetSplitLevel(Int_t split)
void SetOutputFile(TFile *ofile)
TChain * fDs
Definition RnoConvert.h:49
void SetMaxEvents(Int_t n)
void SetBufferSize(Int_t bsize)
Int_t fMaxevt
Definition RnoConvert.h:44
void CreateMainChain()
TChain * fHdr
Definition RnoConvert.h:48
Int_t fBsize
Definition RnoConvert.h:43
Handling of RNO-G detector data.
Definition RnoDetector.h:14
RnoStation * GetStation(Int_t id, Bool_t create)
Handling of RNO-G event data.
Definition RnoEvent.h:14
virtual void Reset()
Definition RnoEvent.cxx:75
Handling of RNO-G event data.
Definition RnoStation.h:19