Subversion Repositories f9daq

Rev

Rev 301 | Blame | Compare with Previous | Last modification | View Log | RSS feed

  1. #include "TFile.h"
  2. #include "TH1.h"
  3. #include "TH2.h"
  4. #include "TFile.h"
  5. #include "TGraph.h"
  6. #include "TMultiGraph.h"
  7. #include "TCanvas.h"
  8. #include "TStyle.h"
  9. #include "TLegend.h"
  10. #include "TFile.h"
  11. //#include "TROOT.h"
  12.  
  13. #include <stdio.h>
  14. #include <time.h>
  15.  
  16. #define DATA_PATH  "./data/"
  17. #define CFRAC 0.4
  18. #define DWIDTH 1024
  19. #define NCH 4
  20. #include "drs.h"
  21.  
  22. ClassImp (drs)
  23.  
  24.  
  25.  TCanvas *c_wfm;
  26.         TFile *froot;
  27.  
  28.         TH1F *h1_chmax[NCH];
  29.         TH1F *h1_chmaxt[NCH];
  30.         TH1F *h1_chmin[NCH];
  31.         TH1F *h1_chmint[NCH];
  32.         TH1F *h1_letime[NCH];
  33.         TH1F *h1_cftime[NCH];
  34.         TH1F *h1_cftdif[NCH];
  35.         TH1F *h1_qdc[NCH];
  36.         TH2F *h2_dpo[NCH];
  37.         TH2F *h2_tdcadc[NCH];
  38.         TH2F *h2_ctdcadc[NCH];
  39.         TH2F *h2_ctalk;
  40.                 TH2F *h2_pos;
  41.                 TH1F *h1_qsumx;
  42.                 TH1F *h1_qsumy;
  43.         TGraph *gr_wfm[NCH];
  44.                 TLegend *leg;
  45.                 TMultiGraph *wfm;
  46.  
  47. int drs::histoinit( int nch, DrsChannel* d,  int trgch, double ymin, double ymax){
  48.         char histname[128];
  49.             c_wfm= new TCanvas("c_wfm","WFM",50,50,700,700);
  50.        
  51.         c_wfm->Divide(1,1);
  52.         c_wfm->cd(1);
  53.         const double maxt=200;
  54.         leg = new TLegend(0.8, 0.8, 1, 1);
  55.                 wfm = new TMultiGraph();
  56.         for (int i=0; i<nch; i++) {
  57.                 sprintf(histname,"chmax_%02d",i);
  58.                 h1_chmax[i]=new TH1F(histname,Form("%d"),1500,d[i].ymin,d[i].ymax);
  59.                 sprintf(histname,"chmaxt_%02d",i);
  60.                 h1_chmaxt[i]=new TH1F(histname,histname,DWIDTH,-0.1,maxt);
  61.                 sprintf(histname,"chmin_%02d",i);
  62.                 h1_chmin[i]=new TH1F(histname,histname,1500,d[i].ymin,d[i].ymax);
  63.                 sprintf(histname,"chmint_%02d",i);
  64.                 h1_chmint[i]=new TH1F(histname,histname,DWIDTH,-0.1,maxt);
  65.                 sprintf(histname,"letime_%02d",i);
  66.                 h1_letime[i]=new TH1F(histname,histname,500,d[i].twin[0],d[i].twin[1]);
  67.                 sprintf(histname,"cftime_%02d",i);
  68.                 h1_cftime[i]=new TH1F(histname,histname,500,d[i].twin[0],d[i].twin[1]);
  69.                 sprintf(histname,"cftdif_%02d",i);
  70.                 h1_cftdif[i]=new TH1F(histname,histname,2000,d[i].twin[0]-d[trgch].twin[0]-40,d[i].twin[1]-d[trgch].twin[0]+20);
  71.                 sprintf(histname,"qdc_%02d",i);
  72.                                 float dx = 7*(d[i].ymax - d[i].ymin);
  73.                
  74.                                 double ymin, ymax;
  75.                                 if(d[i].edge){
  76.                                         ymin=d[i].ymin*8;
  77.                                         ymax=d[i].ymax*8;
  78.                                        
  79.                                 }
  80.                                 else{
  81.                                         ymin=-d[i].ymax*8;
  82.                                         ymax=-d[i].ymin*8;                                     
  83.                                        
  84.                                 }
  85.                                
  86.                                 h1_qdc[i]=new TH1F(histname,histname,800, ymin , ymax);
  87.                 sprintf(histname,"dpo_%02d",i);
  88.                 h2_dpo[i]=new TH2F(histname,histname,DWIDTH,-0.1,maxt,1500,d[i].ymin,d[i].ymax);
  89.                 sprintf(histname,"tdcadc_%02d",i);
  90.                 h2_tdcadc[i]=new TH2F(histname,histname,250,d[i].ymin,d[i].ymax,400,d[i].twin[0]-d[trgch].twin[0],d[i].twin[1]-d[trgch].twin[0]);
  91.                 sprintf(histname,"ctdcadc_%02d",i);
  92.                 h2_ctdcadc[i]=new TH2F(histname,histname,250,d[i].ymin,d[i].ymax,400,d[i].twin[0]-d[trgch].twin[0],d[i].twin[1]-d[trgch].twin[0]);
  93.                 gr_wfm[i]=new TGraph(DWIDTH);
  94.                 //gr_wfm[i]->SetMinimum(d[i].ymin);
  95.                 //gr_wfm[i]->SetMaximum(d[i].ymax);
  96.                 gr_wfm[i]->SetLineColor(i+1);
  97.                 gr_wfm[i]->SetMarkerColor(i+1);
  98.                 sprintf(histname,"ch %d",i);
  99.                 leg->AddEntry(gr_wfm[i], histname, "l");
  100.                                 wfm->Add(gr_wfm[i]);
  101.                 //if (i) gr_wfm[i]->Draw("LP"); else gr_wfm[i]->Draw("ALP");
  102.                 //gr_wfm[i]->GetHistogram()->SetTitle("WaveForms");
  103.                                 //gr_wfm[i]->GetHistogram()->SetDirectory(0);
  104.         }
  105.                 wfm->SetMinimum(ymin);
  106.         wfm->SetMaximum(ymax);
  107.         wfm->Draw("ALP");
  108.                
  109.         //wfm->SetTitle("WaveForms");
  110.             //wfm->SetDirectory(0);
  111.         leg->Draw();
  112.         h2_ctalk=new TH2F("ctalk","ctalk",250,0.,0.5,250,0.,0.5);
  113.         h2_pos=new TH2F("cpos","cpos",40,-1,1,40,-1,1);
  114.                 h1_qsumx=new TH1F("qsumx","qsumx",100,-1,6);
  115.                 h1_qsumy=new TH1F("qsumy","qsumy",100,-1,6);
  116. return 0;              
  117. }
  118.  
  119. drs::drs( const char *FileName, int nch, DrsChannel* d, int trgch,int first, int neve, int updfrq,double ymin, double ymax)
  120.  
  121. {
  122.         int debug=0;
  123.         int i,ich,iev;
  124.         int stat;
  125.         int nev =  first + neve;        
  126. //  char FileName[300]="./data/n2.dat";
  127.         FILE *fin;
  128.  
  129.         char recid[5]="SAMO";
  130.         unsigned short *srecid = (unsigned short *) &recid[0];
  131.  
  132.         struct evrec {
  133.                 unsigned long evn;
  134.                 unsigned short tYear,tMonth,tDay,tHour,tMin,tSec,tmSec,rc;
  135.         } evrec;
  136.  
  137.         float dtdata[NCH][DWIDTH];
  138.         float tcdata[NCH][DWIDTH];
  139.         float tccalib[NCH];
  140.         float scaler[NCH];
  141.         float *t;
  142.         float *y;
  143.         float *ay;
  144.  
  145.         short lwfm[NCH];
  146.         int nwfm[NCH];
  147.         int board;
  148.         int version;
  149.  
  150.  
  151.         unsigned short ichdat[NCH][DWIDTH];
  152.         float chdat[NCH][DWIDTH];
  153.         float achdat[NCH][DWIDTH];
  154.         float chmax[NCH][2];
  155.         float chmin[NCH][2];
  156.         float letime[NCH],cftime[NCH],qdc[NCH];
  157.  
  158.  
  159. //  gROOT->Reset();
  160.  
  161.         gStyle->SetOptStat(1111111);
  162. //  gStyle->SetPalette(52);
  163.        
  164.                 histoinit(nch, d, trgch, ymin,ymax);
  165.                 char rootname[0xff];
  166.         sprintf(rootname, "%s_x%d_y%d.root", FileName, 0, 0);
  167.                
  168.                 froot = new TFile(rootname,"RECREATE");
  169.        
  170.  
  171.         fin=fopen(FileName,"rb");
  172.         if (fin==NULL) {
  173.                 printf("Error opening file %s\n",FileName);
  174.                 return;
  175.         }
  176.         printf("%s, %lu\n",recid,sizeof(evrec));
  177.  
  178.  
  179.         for (ich=0; ich<nch; ich++) {
  180.                 nwfm[ich]=0;
  181.                 for (i=0; i<DWIDTH; i++) {
  182.                         achdat[ich][i]=0;
  183.                 }
  184.         }
  185.         unsigned int trgcell = 0;
  186.         int event = -1;
  187.  
  188.         iev=0;
  189.         time_t t0,told;
  190.         time(&t0);
  191.         told=t0-1;
  192.                 int posrec[9];
  193.         while(!feof(fin)) {
  194.                 if (iev==nev) break;
  195.                 stat=fread(recid,1,4,fin);
  196.                 if (debug) {
  197.                         if (recid[1]=='#') printf("%c%c%u\n",recid[0],recid[1],srecid[1]);
  198.                         else printf("%s\n",recid );
  199.                 }
  200.                 switch (recid[0]) {
  201.                                        
  202.                                         case 'P':{ // position record
  203.                                                                               int len;
  204.                                                                           stat=fread(&len,1,4,fin);
  205.                                                                                   if (len>0 && len<37) {
  206.                                                      stat=fread(posrec,1,len,fin);
  207.                                                                                                          if (recid[3]=='R'){
  208.                                                                                                                 if (iev>=first){
  209.                                                                                                                   froot->Write();
  210.                                                                                                                   froot->Close();
  211.                                                                                                                   sprintf(rootname, "%s_x%d_y%d.root", FileName, posrec[3], posrec[4]);
  212.                                                           froot = new  TFile(rootname,"RECREATE");
  213.                                                                                                                   histoinit(nch,d, trgch, ymin,ymax);
  214.                                                                                                             }  
  215.                                                                                                         printf("Position  len=%d x=%d y=%d ix=%d iy=%d\n", len, posrec[0], posrec[1], posrec[3], posrec[4]);
  216.                                                                                                      } else
  217.                                                         printf("Position record buffer length %d\n", len);                                                                                                               
  218.                                                                                   } else {
  219.                                                                                                         printf("Wrong buffer length %d\n", len);
  220.                                                                                   }
  221.                         break;
  222.                                                                 }
  223.                 case 'D': // DRS
  224.                         version = atoi(&recid[3]);
  225.                         break;
  226.                 case 'T': // TIME
  227.                         switch (recid[1]) {
  228.                         case 'I':
  229.                                 event=0;
  230.                                 if (version ==0) version++;
  231.                                 break;
  232.                         case '#': trgcell = srecid[1]; break;
  233.                         }
  234.                         break;
  235.                 case 'B':
  236.                         board = srecid[1];
  237.                         break;
  238.                 case 'C':
  239.                         sscanf(recid,"C%d",&ich);
  240.                         ich=ich-1;
  241.                         if (event) {
  242.                                 if (version>1) stat=fread(&scaler[ich],1,sizeof(float),fin);
  243.                                 stat=fread(ichdat[ich],1,sizeof(ichdat[ich]),fin);
  244.                         } else {
  245.  
  246.                                 stat=fread(dtdata[ich],1,sizeof(float)*DWIDTH,fin);
  247.  
  248.                         }
  249.                         break;
  250.                 case 'E':
  251.                         stat=fread(&evrec.evn,1,sizeof(evrec),fin);
  252.                         iev++;
  253.                         if (!version)  {
  254.                                 stat=fread(dtdata[0],1,sizeof(float)*DWIDTH,fin);
  255.                         }
  256.                         event=1;
  257.                         break;
  258.                 default:
  259.                         printf("Unknown header! %s\n", recid);
  260.                         break;
  261.                 }
  262.                 if (iev<first) continue;
  263.                 if (recid[0] == 'E') {
  264.                         for ( i=0; i<nch; i++) {
  265.                                 if (!lwfm[i]) continue;
  266.                                                                 double pulseheight =  (d[i].edge)? chmax[i][0]: chmin[i][0];
  267.                                                                        
  268.                                                                 if ((pulseheight >d[i].vcut[0])&&
  269.                                                                     (pulseheight <d[i].vcut[1]) ) {
  270.                                                                    
  271.                                         h1_cftdif[i]->Fill(cftime[i]-cftime[trgch]);
  272.                                 }
  273.                                 if (!d[i].edge) pulseheight = - pulseheight;
  274.                                 h2_tdcadc[i] ->Fill(pulseheight,letime[i]-cftime[trgch]);
  275.                                 h2_ctdcadc[i]->Fill(pulseheight,cftime[i]-cftime[trgch]);
  276.                                                                
  277.                                                                 //if (i==1) printf("ph=%f %f\n", pulseheight, letime[i]-cftime[trgch]);
  278.                         }
  279.                                                 double qsumx=qdc[0]+qdc[1];
  280.                                                 double qsumy=qdc[2]+qdc[3];
  281.                                                 double x0=(qsumx)?(qdc[0]-qdc[1])/qsumx:-2;
  282.                                                 double x1=(qsumy)?(qdc[2]-qdc[3])/qsumy:-2;
  283.                                                 h2_pos->Fill(x0,x1);
  284.                                                 h1_qsumx->Fill(qsumx);
  285.                                                 h1_qsumy->Fill(qsumy);
  286.                                                 double ph0 =  (d[0].edge)? chmax[0][0]: -chmin[0][0];
  287.                                                 double ph1 =  (d[1].edge)? chmax[1][0]: -chmin[1][0];
  288.                         h2_ctalk->Fill(ph0,ph1);
  289.                 }
  290.                                
  291.                 if (event<1) continue;
  292.                 if (recid[0]!='C') continue;
  293.                 if (ich>=nch) continue;
  294.  
  295.                 lwfm[ich]=0;
  296.                 qdc[ich]=d[ich].offset;
  297.                 chmax[ich][0]=-1; chmax[ich][1]=-1;
  298.                 chmin[ich][0]=1; chmin[ich][1]=-1;
  299.                 letime[ich]=-1.; letime[ich]=-1.;
  300.                 cftime[ich]=-1.; cftime[ich]=-1.;
  301.  
  302.  
  303.                 lwfm[ich]=1;
  304.                 nwfm[ich]+=1;
  305.  
  306.  
  307.  
  308.                 t  = &tcdata[ich][0];
  309.                 y  = &chdat[ich][0];
  310.                 ay = &achdat[ich][0];
  311.  
  312.                 float sum=0;
  313.                 for (i=0; i<DWIDTH; i++) {
  314.                         sum+= dtdata[ich][(i+trgcell)%DWIDTH];
  315.                         if (version) t[i]=sum;
  316.                         else t[i]=dtdata[0][i]; // old format
  317.                 }
  318.                 tccalib[ich]=t[(DWIDTH-trgcell)%DWIDTH];
  319.  
  320.                 for (i=0; i<DWIDTH; i++) {
  321.                         if (version) t[i]-=(tccalib[ich]-tccalib[0]); // assume ich 0 is connected and always first in the data
  322.                         y[i]=-0.5+(float)ichdat[ich][i]/0xFFFF;
  323.                         if (version) y[i]+=evrec.rc/1000.;
  324.                         h2_dpo[ich]->Fill(t[i],y[i]);
  325.                         gr_wfm[ich]->SetPoint(i,t[i],y[i]+0.05*(ich-1));
  326.                 }
  327.  
  328.  
  329.                 for (i=0; i<DWIDTH; i++) {
  330.                         // accumulated waveforms
  331.                         ay[i]+=y[i];
  332.  
  333.                         // signal minumum
  334.                         if (y[i]>chmax[ich][0]) {
  335.                                 chmax[ich][0]=y[i];
  336.                                 chmax[ich][1]=t[i];
  337.                         }
  338.  
  339.                         // signal maximum
  340.                         if (y[i]<chmin[ich][0]) {
  341.                                 chmin[ich][0]=y[i];
  342.                                 chmin[ich][1]=t[i];
  343.                         }
  344.                         float t0  = t[i];
  345.                         float t0p = t[i+1];
  346.                         float t0n = t[i-1];
  347.                         // charge in the range
  348.                         if ((t0>d[ich].adcgate[0])&&(t0<d[ich].adcgate[1])) {
  349.                                 float sign =  (d[ich].edge)? 1: -1;
  350.                                                                 qdc[ich]+=sign*y[i]*(t0p-t0n)/2.;
  351.                         }
  352.                 }
  353.  
  354.                 for (i=1; i<DWIDTH-1; i++) {
  355.                         float t0  = t[i];
  356.                         float t0p = t[i+1];
  357.                         float t0n = t[i-1];
  358.                         // leading edge in the range
  359.                         if ((t0>d[ich].twin[0])&&(t0<d[ich].twin[1])) {
  360.                                 if ( ((!d[ich].edge)&&(y[i]<d[ich].threshold) ||
  361.                                       ( d[ich].edge)&&(y[i]>d[ich].threshold ) )
  362.                                      &&(letime[ich]<0.)) {
  363.                                         letime[ich]=t0n+(t0-t0n)*(d[ich].threshold-y[i-1])/(y[i]-y[i-1]);
  364.                                 }
  365.                                 // constant fraction time in the range
  366.                                                                 if (cftime[ich]<0.){
  367.                                                                   // negative signals  
  368.                                   if ((!d[ich].edge) && (chmin[ich][0]<d[ich].threshold)&&(y[i]<chmin[ich][0]*d[ich].cfrac )){
  369.                                                                           cftime[ich]=t0n+(t0-t0n)*(chmin[ich][0]*d[ich].cfrac-y[i-1])/(y[i]-y[i-1]);
  370.                                                                   }
  371.                                                                   // positive signals                                                                          
  372.                                   if (( d[ich].edge) && (chmax[ich][0]>d[ich].threshold)&&(y[i]>chmax[ich][0]*d[ich].cfrac)){
  373.                                                                           cftime[ich]=t0n+(t0-t0n)*(chmax[ich][0]*d[ich].cfrac-y[i-1])/(y[i]-y[i-1]);    
  374.                                                                   }        
  375.                                                                 }  
  376.                         }
  377.                 }
  378.  
  379.                 if ((chmax[ich][1]>d[ich].twin[0])&&(chmax[ich][1]<d[ich].twin[1])) {
  380.                         h1_chmax[ich]->Fill(chmax[ich][0]);
  381.                         if (chmax[ich][0]>d[ich].threshold) h1_chmaxt[ich]->Fill(chmax[ich][1]);
  382.                 }
  383.                 if ((chmin[ich][1]>d[ich].twin[0])&&(chmin[ich][1]<d[ich].twin[1])) {
  384.                         h1_chmin[ich]->Fill(chmin[ich][0]);
  385.                         if (chmin[ich][0]<d[ich].threshold) h1_chmint[ich]->Fill(chmin[ich][1]);
  386.                 }
  387.                 h1_qdc[ich]->Fill(qdc[ich]);
  388.                 h1_letime[ich]->Fill(letime[ich]);
  389.                 h1_cftime[ich]->Fill(cftime[ich]);
  390.                 time(&t0);
  391.                
  392.                 if ((t0!=told) || ((updfrq>0)&&(!(iev%updfrq) ) ) ) {
  393.                         c_wfm->Update();
  394.                         for (i=0; i<DWIDTH; i++) {
  395.                                 ay[i]/=nwfm[ich];
  396.                                 gr_wfm[ich]->SetPoint(i,t[i],ay[i]+0.00*(ich-1));
  397.                         }
  398.                         printf("processing event %d\n",iev);
  399.                         told=t0;
  400.                 }
  401.         }
  402.  
  403.         printf("number of events: %d\n",iev);
  404.         fclose(fin);
  405.                 TString pdfname = Form("%s.pdf",FileName);
  406.         gStyle->SetOptStat(0);
  407.         TCanvas* c_dpo= new TCanvas("c_dpo","DPO",750,50,700,700);
  408.         c_dpo->Divide(2,2);
  409.  
  410.         for (i=0; i<nch; i++) {c_dpo->cd(i+1)->SetLogz(1); h2_dpo[i]->DrawCopy("colz");}
  411.  
  412.         gStyle->SetOptStat(1111111);
  413.         TCanvas* c_cftdif= new TCanvas("c_cftdif","Constant fraction time difference",750,50,700,700);
  414.         c_cftdif->Divide(2,2);
  415.         for (i=0; i<nch; i++) {c_cftdif->cd(i+1); h1_cftdif[i]->DrawCopy();}
  416.                
  417.                 /*
  418.                 TCanvas* c_cftcorr= new TCanvas("c_cftcorr","pulseheight vs timing",750,50,700,700);
  419.         c_cftcorr->Divide(2,2);
  420.         for (i=0; i<nch; i++) {c_cftcorr->cd(i+1); h2_tdcadc[i]->DrawCopy("colz");}
  421.                 */
  422.                 TCanvas* c_adc= new TCanvas("c_adc","adc",750,50,700,700);
  423.         c_adc->Divide(2,2);
  424.         for (i=0; i<nch; i++) {c_adc->cd(i+1); h1_qdc[i]->DrawCopy();}
  425.                
  426.                 /*
  427.                 TCanvas* c_qsum= new TCanvas("c_qsum","qsum",750,50,700,700);
  428.         c_qsum->Divide(1,2);
  429.         c_qsum->cd(1); h1_qsumx->DrawCopy();
  430.                 c_qsum->cd(2); h1_qsumy->DrawCopy();
  431.                
  432.                 gStyle->SetOptStat(0);
  433.                 TCanvas* c_pos= new TCanvas("c_pos","position",750,50,700,700);
  434.                
  435.                 c_pos->cd();
  436.                 h2_pos->DrawCopy("colz");
  437.                 */
  438.                
  439.                 //c_pos->Print(pdfname+"(","pdf");
  440.                 //c_wfm->Print(pdfname,"pdf");
  441.                 c_adc->Print(pdfname+"(","pdf");
  442.                 //c_qsum->Print(pdfname,"pdf");
  443.                 c_dpo->Print(pdfname,"pdf");
  444.                 c_cftdif->Print(pdfname + ")","pdf");
  445.                
  446.                 c_dpo->Write();
  447.                 c_adc->Write();
  448.                 c_cftdif->Write();
  449.                 froot->Write();
  450.                 froot->Close();
  451.  
  452.         for (i=0; i<nch; i++) {c_cftdif->cd(i+1); h1_cftdif[i]->DrawCopy();}
  453.                 froot->Write();
  454.                 froot->Close();
  455.  
  456.         return;
  457. }
  458.  
  459.  
  460. void drsana(int nev=1000, int updfrq=20, char *FileName="./data/x1.dat", int trgch=0){
  461.  
  462. //  float threshold[nch][2]={{0.,-0.02},{0.,0.},{0.,0.},{0.,-0.25}};
  463. //  float twin[nch][2]={{120.,130.},{100.,150.},{100.,150.},{70.,90.}};
  464. //  float adcgate[nch][2]={{12nch.,136.},{120.,140.},{120.,140.},{70.,90.}};
  465. //  float threshold[4][2]={{0.,-0.05},{0.,-0.05},{0.,-0.25},{0.,-0.25}};
  466. //  float twin[4][2]={{80.,100.},{80.,100.},{90.,110.},{50.,70.}};
  467. //  float adcgate[4][2]={{80.,110.},{80.,110.},{90.,110.},{50.,70.}};
  468.         double threshold[4][2]={{0.,-0.05},{0.,-0.05},{0.,-0.25},{0.,-0.25}};
  469.         float twin[4][2]={{100.,120.},{100.,120.},{110.,130.},{60.,80.}};
  470.         float adcgate[4][2]={{105.,125.},{105.,125.},{110.,130.},{60.,80.}};
  471.         double vcut[4][2]={{-0.49,-0.07},{-1,1},{-1,1},{-1,1}};
  472.         int edge[4]={1,1,1,1};// 0 -- rising, 1 -- falling
  473.         float offset[4]={0.25,0.25,0.25,0.25};//
  474.  
  475.         DrsChannel c[2];
  476.  
  477.         c[0].cfrac = CFRAC;
  478.         c[0].threshold = threshold[0][0];
  479.         c[0].twin[0] = twin[0][0];
  480.         c[0].twin[1] = twin[0][1];
  481.  
  482.         c[0].adcgate[0] = adcgate[0][0];
  483.         c[0].adcgate[1] = adcgate[0][1];
  484.  
  485.         c[0].vcut[0] = vcut[0][0];
  486.         c[0].vcut[1] = vcut[0][1];
  487.  
  488.         c[0].edge = edge[0];
  489.         c[0].edge = offset[0];
  490.  
  491.         c[1].cfrac = CFRAC;
  492.         c[1].threshold = threshold[1][0];
  493.         c[1].twin[0] = twin[1][0];
  494.         c[1].twin[1] = twin[1][1];
  495.  
  496.         c[1].adcgate[0] = adcgate[1][0];
  497.         c[1].adcgate[1] = adcgate[1][1];
  498.  
  499.         c[1].vcut[0] = vcut[1][0];
  500.         c[1].vcut[1] = vcut[1][1];
  501.  
  502.         c[1].edge = edge[1];
  503.                 c[1].offset = offset[1];
  504.         drs * d = new drs(FileName, 2, c, trgch, 0, nev, updfrq,  -0.55,0.55 );
  505. }
  506.  
  507.  
  508.  
  509. void sensl(int nev=1000, int updfrq=20, char *FileName="./data/x1.dat", int trgch=3){
  510.         // nev = -1 all events
  511.         double threshold[4]={-0.005,-0.05,-0.005,-0.05}; // threshold for determination of timing
  512.         float twin[4][2]={{20,40.},{980.,1150.},{100.,160.},{100.,160.}}; // window for calculation of the timing
  513.         float adcgate[4][2]={{20.,30.}, {980.,1150.},{112.,123.},{150.,163.}}; //integration window
  514.         double vcut[4][2]={{-1,1},{-1,1},{-1,1},{-1,1}};  // cut on the pulseheight for filling the timing histograms
  515.         int edge[4]={1,0,1,0};// 1 -- rising edge, 0 -- falling edge
  516.         double offset[4]={0.25,0.25, 0.025,0.25};// qdc offset
  517.                 double ymin[4]={-0.1,-0.5,-0.01,-0.6};// ymin offset
  518.         double ymax[4]={0.5,0.5,0.07,0.1};// ymin offset
  519.         DrsChannel c[NCH];
  520.         for (int k=0;k<NCH;k++){
  521.           c[k].cfrac = CFRAC;    
  522.           c[k].threshold = threshold[k];
  523.                   for (int j=0;j<2;j++){
  524.             c[k].twin[j] = twin[k][j];
  525.             c[k].adcgate[j] = adcgate[k][j];
  526.             c[k].vcut[j] = vcut[k][j];        
  527.             c[k].edge = edge[k];
  528.                         c[k].offset = offset[k];
  529.                         c[k].ymin = ymin[k];
  530.                         c[k].ymax = ymax[k];
  531.                   }
  532.         }
  533.         drs * d = new drs(FileName, NCH, c, trgch, 0, nev, updfrq , -1,1);
  534. }
  535.  
  536.  
  537.  
  538. void ana(TString FileName="./data/x1.dat"){
  539.             int nev=-1;
  540.                 int updfrq=1000;
  541.                 int trgch=0;
  542.         // nev = -1 all events
  543.         double threshold[4]={-0.005,-0.05,-0.005,-0.05}; // threshold for determination of timing
  544.         float twin[4][2]={{20,40.},{980.,1150.},{100.,160.},{100.,160.}}; // window for calculation of the timing
  545.         float adcgate[4][2]={{20.,30.}, {980.,1150.},{112.,123.},{150.,163.}}; //integration window
  546.         double vcut[4][2]={{-1,1},{-1,1},{-1,1},{-1,1}};  // cut on the pulseheight for filling the timing histograms
  547.         int edge[4]={1,0,1,0};// 1 -- rising edge, 0 -- falling edge
  548.         double offset[4]={0.25,0.25, 0.025,0.25};// qdc offset
  549.                 double ymin[4]={-0.1,-0.5,-0.01,-0.6};// ymin offset
  550.         double ymax[4]={0.5,0.5,0.07,0.1};// ymin offset
  551.         DrsChannel c[NCH];
  552.         for (int k=0;k<NCH;k++){
  553.           c[k].cfrac = CFRAC;    
  554.           c[k].threshold = threshold[k];
  555.                   for (int j=0;j<2;j++){
  556.             c[k].twin[j] = twin[k][j];
  557.             c[k].adcgate[j] = adcgate[k][j];
  558.             c[k].vcut[j] = vcut[k][j];        
  559.             c[k].edge = edge[k];
  560.                         c[k].offset = offset[k];
  561.                         c[k].ymin = ymin[k];
  562.                         c[k].ymax = ymax[k];
  563.                   }
  564.         }
  565.                 cout << "DRS4 " << FileName << endl;
  566.         drs * d = new drs(FileName.Data(), NCH, c, trgch, 0, nev, updfrq , -0.5,0.5);
  567. }
  568.