midas2root.cxx

Go to the documentation of this file.
00001 // Example Program for converting MIDAS format to ROOT format.
00002 //
00003 // T. Lindner (Jan 2016) 
00004 //
00005 // Example is for the CAEN V792 ADC module
00006 
00007 #include <stdio.h>
00008 #include <iostream>
00009 #include <time.h>
00010 #include <vector>
00011 
00012 #include "TRootanaEventLoop.hxx"
00013 #include "TFile.h"
00014 #include "TTree.h"
00015 
00016 #include "TAnaManager.hxx"
00017 
00018 #ifdef USE_V792
00019 #include "TV792Data.hxx"
00020 #endif
00021 
00022 class Analyzer: public TRootanaEventLoop {
00023 
00024 public:
00025 
00026   // An analysis manager.  Define and fill histograms in 
00027   // analysis manager.
00028   TAnaManager *anaManager;
00029 
00030   // The tree to fill.
00031   TTree *fTree;
00032 
00033 #ifdef USE_V792
00034   // CAEN V792 tree variables
00035   int nchannels;
00036   int adc_value[32];
00037 #endif
00038 
00039 
00040   Analyzer() {
00041 
00042   };
00043 
00044   virtual ~Analyzer() {};
00045 
00046   void Initialize(){
00047 
00048 
00049   }
00050   
00051   
00052   void BeginRun(int transition,int run,int time){
00053     
00054     // Create a TTree
00055     fTree = new TTree("midas_data","MIDAS data");
00056 
00057 #ifdef USE_V792    
00058     fTree->Branch("nchannels",&nchannels,"nchannels/I");
00059     fTree->Branch("adc_value",adc_value,"adc_value[nchannels]/I");
00060 #endif
00061 
00062   }   
00063 
00064 
00065   void EndRun(int transition,int run,int time){
00066         printf("\n");
00067   }
00068 
00069   
00070   
00071   // Main work here; create ttree events for every sequenced event in 
00072   // Lecroy data packets.
00073   bool ProcessMidasEvent(TDataContainer& dataContainer){
00074 
00075     int id = dataContainer.GetMidasEvent().GetSerialNumber();
00076     if(id%10 == 0) printf(".");
00077 
00078 
00079 #ifdef USE_V792    
00080     TV792Data *data = dataContainer.GetEventData<TV792Data>("ADC0");
00081     if(data){
00082       nchannels = 32;
00083       for(int i = 0; i < nchannels;i++) adc_value[i] = 0;
00084       
00085       /// Get the Vector of ADC Measurements.
00086       std::vector<VADCMeasurement> measurements = data->GetMeasurements();
00087       for(unsigned int i = 0; i < measurements.size(); i++){ // loop over measurements
00088         
00089         int chan = measurements[i].GetChannel();
00090         uint32_t adc = measurements[i].GetMeasurement();
00091         
00092         if(chan >= 0 && chan < nchannels)
00093           adc_value[chan] = adc;
00094       }
00095     }
00096 #endif
00097 
00098     fTree->Fill();
00099 
00100     return true;
00101 
00102   };
00103   
00104   // Complicated method to set correct filename when dealing with subruns.
00105   std::string SetFullOutputFileName(int run, std::string midasFilename)
00106   {
00107     char buff[128]; 
00108     Int_t in_num = 0, part = 0;
00109     Int_t num[2] = { 0, 0 }; // run and subrun values
00110     // get run/subrun numbers from file name
00111     for (int i=0; ; ++i) {
00112       char ch = midasFilename[i];
00113         if (!ch) break;
00114         if (ch == '/') {
00115           // skip numbers in the directory name
00116           num[0] = num[1] = in_num = part = 0;
00117         } else if (ch >= '0' && ch <= '9' && part < 2) {
00118           num[part] = num[part] * 10 + (ch - '0');
00119           in_num = 1;
00120         } else if (in_num) {
00121           in_num = 0;
00122           ++part;
00123         }
00124     }
00125     if (part == 2) {
00126       if (run != num[0]) {
00127         std::cerr << "File name run number (" << num[0]
00128                   << ") disagrees with MIDAS run (" << run << ")" << std::endl;
00129         exit(1);
00130       }
00131       sprintf(buff,"output_%.6d_%.4d.root", run, num[1]);
00132       printf("Using filename %s\n",buff);
00133     } else {
00134       sprintf(buff,"output_%.6d.root", run);
00135     }
00136     return std::string(buff);
00137   };
00138 
00139 
00140 
00141 
00142 
00143 }; 
00144 
00145 
00146 int main(int argc, char *argv[])
00147 {
00148 
00149   Analyzer::CreateSingleton<Analyzer>();
00150   return Analyzer::Get().ExecuteLoop(argc, argv);
00151 
00152 }
00153 

Generated on 12 Feb 2016 for ROOT Analyzer by  doxygen 1.6.1