ARA ROOT v3.10 Software

calibration/ATRI/voltageCalib/adcSampleBlock.cxx

00001 
00002 
00003 
00004 
00005 
00006 //AraRoot Includes
00007 #include "UsefulAtriStationEvent.h"
00008 #include "RawAtriStationEvent.h"
00009 #include "araSoft.h"
00010 
00011 
00012 //Root Includes
00013 #include "TTree.h"
00014 #include "TFile.h"
00015 #include "TH1.h"
00016 #include "TTree.h"
00017 #include "TMath.h"
00018 #include "TCanvas.h"
00019 
00020 
00021 //Standard Includes
00022 #include <iostream>
00023 #include <fstream>
00024 
00025 
00026 #define MAX_ADC 4096
00027 
00028 //Prototype Functions
00029 Int_t adcSample(char*, Int_t, Int_t, Double_t);
00030 
00031 
00032 int main(int argc, char **argv)
00033 {
00034   Int_t runNum=0, pedNum=0, dda=0, chan=0;
00035   char baseName[FILENAME_MAX];
00036 
00037   if(argc<5) {
00038     std::cerr << "Usage: " << argv[0] << " <baseDir> <runNum> <pedNum> <voltage> \n";
00039     return 1;
00040   }
00041   sprintf(baseName, argv[1]);
00042   runNum=atoi(argv[2]);
00043   pedNum=atoi(argv[3]);
00044   Double_t voltage = atof(argv[4]);
00045   return adcSample(baseName, runNum, pedNum, voltage);
00046   //return test(baseName, runNum, pedNum, voltage);
00047 
00048 }
00049 
00050 Int_t adcSample(char* baseDirName, Int_t runNum, Int_t pedNum, Double_t voltage){
00051   printf("%s\n", baseDirName);
00052   printf("%i %i %f\n", runNum, pedNum, voltage);
00053 
00054 
00055   char runFileName[FILENAME_MAX];
00056   char pedFileName[FILENAME_MAX];
00057   sprintf(runFileName, "%s/root/run%i/event%i.root", baseDirName, runNum, runNum);
00058   sprintf(pedFileName, "%s/raw_data/run_%06i/pedestalValues.run%06i.dat", baseDirName, pedNum, pedNum);
00059 
00060   fprintf(stderr, "%s\n", runFileName);
00061   fprintf(stderr, "%s\n", pedFileName);
00062   
00063   TFile *fp = new TFile(runFileName);
00064   if(!fp) {
00065     std::cerr << "Can't open file\n";
00066     return -1;
00067   }
00068   TTree *eventTree = (TTree*) fp->Get("eventTree");
00069   if(!eventTree) {
00070     std::cerr << "Can't find eventTree\n";
00071     return -1;
00072   }
00073   RawAtriStationEvent *evPtr=0;
00074   eventTree->SetBranchAddress("event",&evPtr);
00075   Long64_t numEntries=eventTree->GetEntries();
00076 
00077   std::cout << "Number of entries in file is " << numEntries << std::endl;
00078 
00079   Long64_t starEvery=numEntries/80;
00080   if(starEvery==0) starEvery++;
00081   
00082   Int_t stationId=0;
00083   eventTree->GetEntry(0);
00084   stationId= evPtr->stationId;
00085   AraEventCalibrator *calib = AraEventCalibrator::Instance();
00086   calib->setAtriPedFile(pedFileName, stationId);
00087   
00088   
00089   
00090   //General output stuff
00091   char outFileName[FILENAME_MAX];
00092   sprintf(outFileName, "%s/root/run%i/adcSampleBlock.root", baseDirName, runNum);
00093   TFile *outFile = new TFile(outFileName, "RECREATE");
00094 
00095   //  std::cout << "Output file " << outFileName << std::endl;
00096 
00097   //Variables what we need
00098   Int_t chanIndex;
00099   Int_t capArray;
00100   Int_t block;
00101   Int_t logicalBlock;
00102   Int_t sample;
00103   Int_t thisCapArray;
00104   Int_t dda=0, chan=0;
00105   Double_t ****maxPos = new Double_t***[DDA_PER_ATRI];
00106   Double_t ****maxNeg = new Double_t***[DDA_PER_ATRI];
00107 
00108   for(dda=0;dda<DDA_PER_ATRI;dda++){
00109     maxPos[dda] = new Double_t**[RFCHAN_PER_DDA];
00110     maxNeg[dda] = new Double_t**[RFCHAN_PER_DDA];
00111     for(chan=0;chan<RFCHAN_PER_DDA;chan++){
00112       maxPos[dda][chan] = new Double_t*[SAMPLES_PER_BLOCK];
00113       maxNeg[dda][chan] = new Double_t*[SAMPLES_PER_BLOCK];
00114       for(sample=0;sample<SAMPLES_PER_BLOCK;sample++){
00115         maxPos[dda][chan][sample] = new Double_t[BLOCKS_PER_DDA];
00116         maxNeg[dda][chan][sample] = new Double_t[BLOCKS_PER_DDA];
00117         for(block=0;block<BLOCKS_PER_DDA;block++){
00118           maxPos[dda][chan][sample][block]=-4096;
00119           maxNeg[dda][chan][sample][block]=4096;
00120         }
00121       }
00122     }
00123   }
00124   
00125   // //  numEntries=1;//FIXME
00126   for(int entry=10;entry<numEntries;entry++){
00127     if(entry%starEvery==0) std::cerr <<"*";
00128     eventTree->GetEntry(entry);
00129     UsefulAtriStationEvent realEvent(evPtr, AraCalType::kJustPed);
00130     capArray = evPtr->blockVec[0].getCapArray(); //capArray of first block
00131     
00132     for(dda=0;dda<DDA_PER_ATRI;dda++){
00133       for(chan=0;chan<RFCHAN_PER_DDA;chan++){
00134         chanIndex=chan+RFCHAN_PER_DDA*dda;
00135         TGraph *gr = realEvent.getGraphFromElecChan(chanIndex);
00136         Int_t numSamples = gr->GetN();
00137         Int_t numBlocks = numSamples/SAMPLES_PER_BLOCK;
00138         Double_t *yVals=gr->GetY();
00139 
00140         for(block=0; block<numBlocks; block++){ 
00141           if(block%2) thisCapArray=1-capArray;
00142           else thisCapArray=capArray;
00143           logicalBlock=evPtr->blockVec[block*DDA_PER_ATRI+dda].getBlock();
00144           for(sample=0;sample<SAMPLES_PER_BLOCK;sample++){
00145             if(yVals[sample+block*SAMPLES_PER_BLOCK]>maxPos[dda][chan][sample][logicalBlock]) maxPos[dda][chan][sample][logicalBlock]=yVals[sample+block*SAMPLES_PER_BLOCK];
00146             if(yVals[sample+block*SAMPLES_PER_BLOCK]<maxNeg[dda][chan][sample][logicalBlock]) maxNeg[dda][chan][sample][logicalBlock]=yVals[sample+block*SAMPLES_PER_BLOCK];
00147           }
00148         }//block
00149         
00150         delete gr;
00151       }//chan
00152     }//dda
00153 
00154   }//entry
00155   std::cerr << "\n";
00156 
00157   TTree *maxPosNegTree = new TTree("maxPosNegTree", "Tree of maximum positive and negative values");
00158   maxPosNegTree->Branch("dda", &dda, "dda/I");
00159   maxPosNegTree->Branch("chan", &chan, "chan/I");
00160   maxPosNegTree->Branch("sample", &sample, "sample/I");
00161   maxPosNegTree->Branch("block", &block, "block/I");
00162   Double_t maxNegValue=0, maxPosValue=0;
00163   maxPosNegTree->Branch("maxPosValue", &maxPosValue, "maxPosValue/D");
00164   maxPosNegTree->Branch("maxNegValue", &maxNegValue, "maxNegValue/D");
00165   maxPosNegTree->Branch("voltage", &voltage, "voltage/D");
00166   
00167   for(dda=0;dda<DDA_PER_ATRI;dda++){
00168     for(chan=0;chan<RFCHAN_PER_DDA;chan++){
00169       for(sample=0;sample<SAMPLES_PER_BLOCK;sample++){
00170         for(block=0;block<BLOCKS_PER_DDA;block++){
00171           maxPosValue=maxPos[dda][chan][sample][block];
00172           maxNegValue=maxNeg[dda][chan][sample][block];
00173           maxPosNegTree->Fill();
00174         }
00175       }
00176     }
00177   }
00178 
00179   maxPosNegTree->Write();
00180 
00181   outFile->Write();
00182 
00183   for(dda=0;dda<DDA_PER_ATRI;dda++){
00184     for(chan=0;chan<RFCHAN_PER_DDA;chan++){
00185       for(sample=0;sample<SAMPLES_PER_BLOCK;sample++){
00186         delete [] maxPos[dda][chan][sample];
00187         delete [] maxNeg[dda][chan][sample];
00188       }
00189       delete [] maxPos[dda][chan];
00190       delete [] maxNeg[dda][chan];
00191     }
00192     delete [] maxPos[dda];
00193     delete [] maxNeg[dda];
00194  }
00195  
00196 
00197   return 0;
00198 }
00199 
00200 

Generated on Tue Jul 16 16:58:02 2013 for ARA ROOT v3.10 Software by doxygen 1.4.7