////////////////////////////////////////////////////////// // // Tracer.C - "makeclass" for root "tree1" trees generated // by the uconn.cc utility. // // author: richard.t.jones at uconn.edu // version: june 1, 2014 // // example usage: // $ root -l .root // root> .L Tracer.C // root> TTree *t=tree1; // root> Tracer tr(t); // root> tr.Show(6,0,0,2); // to see event 2, channel 6 // // This class has been automatically generated on // Tue Feb 11 14:01:19 2014 by ROOT version 5.34/09 // from TTree tree1/Sipm test // found on file: data/simple_91.root ////////////////////////////////////////////////////////// #define Tracer_cxx #include "Tracer.h" #include #include #include #include // See the comments in pulse_analyzer.C for details. // You must define exactly one of the following: // PHEIGHT_TRIG_SCANNING // PHEIGHT_TRIG_SUMMER // PHEIGHT_TRIG_PULSER #define PHEIGHT_TRIG_PULSER TH1D* Tracer::Loop(int channel, int pulse, int fadc) { // In a ROOT session, you can do: // Root > .L Tracer.C // Root > Tracer t // Root > t.GetEntry(12); // Fill t data members with entry number 12 // Root > t.Show(); // Show values of entry 12 // Root > t.Show(16); // Read and show values of entry 16 // Root > t.Loop(); // Loop on all entries // // This is the loop skeleton where: // jentry is the global entry number in the chain // ientry is the entry number in the current Tree // Note that the argument to GetEntry must be: // jentry for TChain::GetEntry // ientry for TTree::GetEntry and TBranch::GetEntry // // To read only selected branches, Insert statements like: // METHOD1: // fChain->SetBranchStatus("*",0); // disable all branches // fChain->SetBranchStatus("branchname",1); // activate branchname // METHOD2: replace line // fChain->GetEntry(jentry); //read all branches // b_branchname->GetEntry(ientry); //read only this branch if (fChain == 0) return 0; if (pulse != 0) return 0; char title[80]; sprintf(title,"pulse integral spectrum, channel %d", channel); TH1D *h1d = new TH1D("spec",title,500,-2.,498.); Long64_t nentries = fChain->GetEntriesFast(); pulse_analyzer pana(fadc_raw[fadc][0], 16, fWindowWidth); Long64_t nbytes = 0, nb = 0; for (Long64_t jentry=2; jentry < nentries; jentry++) { Long64_t ientry = LoadTree(jentry); if (ientry < 0) break; nb = fChain->GetEntry(jentry); nbytes += nb; if (jentry/100000*100000 == jentry) { std::cout << jentry << std::endl; } overflow_fill(fadc_raw[fadc][channel],fWindowWidth); #if PHEIGHT_TRIG_SCANNING h1d->Fill(pana.pheight_trig_scanning(channel)); #elif PHEIGHT_TRIG_SUMMER h1d->Fill(pana.pheight_trig_summer(channel)); #else h1d->Fill(pana.pheight_trig_pulser(channel)); #endif } TCanvas *c1 = (TCanvas*)gROOT->FindObject("c1"); if (c1 == 0) { c1 = new TCanvas("c1","c1",0,20,600,500); } h1d->Draw(); return h1d; } TH1D* Tracer::Spec(int channel, int pulse, int fadc) { // Generate a histogram containing the pulse height spectrum // for the given channel of the given fadc, using the data // collected in raw mode (pulse=0), and return the result // as a pointer to a new TH1D object, which the user code is // responsible to delete when it is no longer needed. if (fChain == 0) return 0; if (pulse != 0) return 0; TH1D *h1d; char title[80]; sprintf(title,"triggered pulse spectrum, channel %d", channel); if (channel < 6) { h1d = new TH1D("spec",title,500,-2.,498.); } else { h1d = new TH1D("spec",title,500,-2.,498.); } Long64_t nentries = fChain->GetEntriesFast(); pulse_analyzer pana(fadc_raw[fadc][0], 16, fWindowWidth); Long64_t nbytes = 0, nb = 0; for (Long64_t jentry=2; jentry < nentries; jentry++) { Long64_t ientry = LoadTree(jentry); if (ientry < 0) break; nb = fChain->GetEntry(jentry); nbytes += nb; if (jentry/100000*100000 == jentry) { std::cout << jentry << std::endl; } overflow_fill(fadc_raw[fadc][channel],fWindowWidth); #if PHEIGHT_TRIG_SCANNING h1d->Fill(pana.pheight_trig_scanning(channel)); #elif PHEIGHT_TRIG_SUMMER h1d->Fill(pana.pheight_trig_summer(channel)); #else h1d->Fill(pana.pheight_trig_pulser(channel)); #endif } TCanvas *c1 = (TCanvas*)gROOT->FindObject("c1"); if (c1 == 0) { c1 = new TCanvas("c1","c1",0,20,600,500); } h1d->Draw(); return h1d; } TProfile* Tracer::Prof(int channel, int pulse, int fadc) { // Generate a profile histogram containing the pulse height series // data vs event number, for the given channel of the given fadc, // using data collected in raw mode (pulse=0), and return the result // as a pointer to a new TProfile object, which the user code is // responsible to delete when it is no longer needed. if (fChain == 0) return 0; if (pulse != 0) return 0; TProfile *prof; char title[80]; sprintf(title,"pulse height vs event no., channel %d", channel); Long64_t nentries = fChain->GetEntriesFast(); prof = new TProfile("prof",title,nentries,0,nentries,0.,4000.); pulse_analyzer pana(fadc_raw[fadc][0], 16, fWindowWidth); Long64_t nbytes = 0, nb = 0; for (Long64_t jentry=2; jentry < nentries; jentry++) { Long64_t ientry = LoadTree(jentry); if (ientry < 0) break; nb = fChain->GetEntry(jentry); nbytes += nb; if (jentry/100000*100000 == jentry) { std::cout << jentry << std::endl; } overflow_fill(fadc_raw[fadc][channel],fWindowWidth); #if PHEIGHT_TRIG_SCANNING prof->Fill(jentry, pana.pheight_trig_scanning(channel)); #elif PHEIGHT_TRIG_SUMMER prof->Fill(jentry, pana.pheight_trig_summer(channel)); #else prof->Fill(jentry, pana.pheight_trig_pulser(channel)); #endif } TCanvas *c1 = (TCanvas*)gROOT->FindObject("c1"); if (c1 == 0) { c1 = new TCanvas("c1","c1",0,20,600,500); } prof->Draw(); return prof; } int Tracer::overflow_fill(Int_t *fadc_array, int fadc_size) { // It seems that when raw fadc sequence data are being reported, one // sees pulse heights that overflow the maximum input range get written // out as zero. Assume any zero values in the raw fadc data array to be // overflows, and replace them with the saturation vaue 4096. int nfilled = 0; for (int i=0; i < fadc_size; ++i) { if (fadc_array[i] == 0) { fadc_array[i] = 4096; ++nfilled; } } return nfilled; } TTree *Tracer::Tree(char *tname, int pulse, int fadc) { // From the raw fadc sequence data (pulse=0) contained in the input tree, // compute pulse heights for each channel of board fadc, and save them // in a tree with just one number per channel per event. It is the user's // responsibility to delete this tree when it is no longer neeed. if (fChain == 0) return 0; if (pulse != 0) return 0; TTree *tree = new TTree(tname, tname); double pheight[16]; tree->Branch("pheight",&pheight,"pheight[16]/D"); Long64_t nentries = fChain->GetEntriesFast(); pulse_analyzer pana(fadc_raw[fadc][0], 16, fWindowWidth); Long64_t nbytes = 0, nb = 0; for (Long64_t jentry=2; jentry < nentries; jentry++) { Long64_t ientry = LoadTree(jentry); if (ientry < 0) break; nb = fChain->GetEntry(jentry); nbytes += nb; if (jentry/100000*100000 == jentry) { std::cout << jentry << std::endl; } for (int channel = 1; channel <= 16; ++channel) { overflow_fill(fadc_raw[fadc][channel],fWindowWidth); #if PHEIGHT_TRIG_SCANNING pheight[channel] = pana.pheight_trig_scanning(channel); #elif PHEIGHT_TRIG_SUMMER pheight[channel-1] = pana.pheight_trig_summer(channel); #else pheight[channel-1] = pana.pheight_trig_pulser(channel); #endif } tree->Fill(); } TCanvas *c1 = (TCanvas*)gROOT->FindObject("c1"); if (c1 == 0) { c1 = new TCanvas("c1","c1",0,20,600,500); } tree->Draw("pheight[1]:pheight[6]","pheight[1] > 10"); return tree; } #include