• Main Page
  • Related Pages
  • Namespaces
  • Data Structures
  • Files
  • File List
  • Globals

onlinemonitor/spadicmonitor/TSpadicProc.cxx (r4864/r3724)

Go to the documentation of this file.
00001 #include "TSpadicProc.h"
00002 
00003 #include "TGo4Analysis.h"
00004 #include "TROOT.h"
00005 #include "TSystem.h"
00006 #include "TPrincipal.h"
00007 
00008 #ifndef WITHOUT_SPADIC
00009   // still need this for subevent types:
00010   #include "roc/Board.h"
00011   #include "spadic/Message.h"
00012 #endif
00013 
00014 #define SPADIC_DEBUG 0
00015 
00016 #include "Riostream.h"
00017 using namespace std;
00018 
00019 // 25.10.2011 replace Susibo 03 with 10
00020   Int_t SpaPosMap[MAX_SPADIC] = {  0,  1, 11, 10, 15,  2, 12,  6, 18,  4, 16, 17,  5,  0,  0,  0,  0,  0,  0,  0 };  // map of 12 used Spadics           
00021 //  // 19.10.2011                   // MS                           // FF           // unused
00022 //  Int_t SpMap[MAX_SPADIC]   = {  0,  1,  2, 11, 12,  3,  6, 15, 18,  4, 16, 17,  5, 10, 14,  7,  8,  9, 13, 19 };  // map of Spadics in TRD order
00023 //  Int_t SuMap[MAX_SPADIC]   = {  0,  1,  5,  3,  9, 12,  7, 17, 18, 15, 13,  2,  6, 16, 14,  4, 10, 11,  8, 19 };  // histo mapping
00024   //                               0,  1,  2,  3,  4,  5,  6,  7,  8,  9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19
00025 //Int_t SpaPosMap[MAX_SPADIC] = {  0,  1, 11,  3, 15,  2, 12,  6, 18,  4, 16, 17,  5,  0,  0,  0,  0,  0,  0,  0 };  // map of 12 used Spadics           
00026 
00027 TSpadicProc::TSpadicProc(const char* name) :
00028   TCBMBeamtimeProc(name),
00029   fOutputEvent(0),
00030   fDoStopAnalysis(kFALSE)
00031 {
00032   TGo4Log::Info("TTRDSpadicProc: Create instance %s", GetName());
00033    
00034   fPar = (TSpadicParam*) MakeParameter("SpadicPar", "TSpadicParam", "set_SpadicPar.C");
00035 
00036   if (fPar) fPar->SetConfigSpadics();
00037 
00038   /* simple SPADIC online JAM:*/
00039   TString obname;
00040   TString obtitle;
00041   TString foldername;
00042 
00043   UInt_t sid, ch;
00044 
00045   for (sid = 0; sid < MAX_SPADIC; sid++){
00046  
00047         fSpadic_NbOfOverflow[sid] = 0;
00048         fSpadic_NbOfChOverThreshold[sid] = 0;
00049         fSpadic_NoiseDist[sid] = 0;
00050         fSpadic_NoiseDist2D[sid] = 0;
00051         fSpadic_PedelPos[sid] = 0;
00052         fSpadic_PedelPos2D[sid] = 0;
00053         fSpadic_PRF[sid] = 0;
00054         fSpadic_trace2D[sid] = 0;
00055     for (ch = 0; ch < NUM_SPADIC_CHA; ch++){
00056        fSpadic_trace[sid][ch] = 0;
00057        fSpadic_trace_clean[sid][ch] = 0;
00058        fSpadic_ADCdist[sid][ch] = 0;
00059     }
00060     fSpadic_spectrum[sid] = 0;
00061     fSpadic_ClusterWidth[sid] = 0;
00062     fSpadic_shape[sid] = 0;
00063     fSpadic_shape2Draw[sid] = 0;
00064     fSpadic_shape2D[sid] = 0;
00065     fSpadic_shapecnt[sid] = 0;
00066     fSpadic_PadMaxDist[sid] = 0;
00067     fSpadic_peak[sid] = 0;
00068     fSpadic_meanpos[sid] = 0;
00069 
00070          if (TSpadicEvent::ConfigSpadics[sid] == 0) continue;
00071 
00072          TGo4Log::Info("Create histograms for SPADIC %d", sid);
00073 
00074     obname.Form("Spadic_%d/Spadic%d_NbOfOverflow", sid, sid);
00075     obtitle.Form("Spadic %d Nb of Events in Overflow", sid);
00076     fSpadic_NbOfOverflow[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5,  "Channel Number", "Nb of Overflows");
00077     
00078     obname.Form("Spadic_%d/Spadic%d_NbOfChOverThreshold", sid, sid);
00079     obtitle.Form("Spadic %d Nb of bins over Threshold", sid);
00080     fSpadic_NbOfChOverThreshold[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5,  "Channel Number", "Nb of Events over Threshold");
00081     
00082     obname.Form("Spadic_%d/Spadic%d_NoiseDist", sid, sid);
00083     obtitle.Form("RMS of ADC distribution of spadic %d", sid);
00084     fSpadic_NoiseDist[sid] = MakeTH1('D', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5,  "Channel Number", "RMS of ADC distribution [ADC values]");
00085 
00086     obname.Form("Spadic_%d/Spadic%d_NoiseDist2D", sid, sid);
00087     obtitle.Form("RMS of ADC distribution of spadic %d 2D", sid);
00088     fSpadic_NoiseDist2D[sid] = MakeTH2('D', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5, 305, -50, 255,  "Channel Number", "RMS of ADC distribution [ADC values]");
00089     fSpadic_NoiseDist2D[sid]->SetContour(99); 
00090 
00091     obname.Form("Spadic_%d/Spadic%d_PedelPos", sid, sid);
00092     obtitle.Form("Pedestel (Mean of ADC value distribution) of spadic %d", sid);
00093     fSpadic_PedelPos[sid] = MakeTH1('D', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5,  "Channel Number", "Mean of ADC distribution [ADC values]");
00094 
00095     obname.Form("Spadic_%d/Spadic%d_PedelPos2D", sid, sid);
00096     obtitle.Form("Pedestel (Mean of ADC value distribution) of spadic %d 2D", sid);
00097     fSpadic_PedelPos2D[sid] = MakeTH2('D', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5, 305, -50, 255,  "Channel Number", "Mean of ADC distribution [ADC values]");
00098     fSpadic_PedelPos2D[sid]->SetContour(99); 
00099 
00100     obname.Form("Spadic_%d/Spadic%d_PRF", sid, sid);
00101     obtitle.Form("PRF of spadic %d", sid);
00102     fSpadic_PRF[sid] = MakeTH2('D', obname.Data(), obtitle.Data(), 400, -1.5, 1.5, 100, 0, 1);
00103     fSpadic_PRF[sid]->SetContour(99); 
00104     
00105     obname.Form("Spadic_%d/Spadic%d_trace2D", sid, sid);
00106     obtitle.Form("Trace2D Susibo_%02d", sid);
00107     fSpadic_trace2D[sid] = MakeTH2('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE, NUM_SPADIC_CHA, -0.5, NUM_SPADIC_CHA-0.5, "time (40ns)","Channel ID");
00108     fSpadic_trace2D[sid]->SetContour(99);
00109  
00110     for (ch = 0; ch < NUM_SPADIC_CHA; ch++){
00111       obname.Form("Spadic_%d/Channel%d/Spadic%d_ch%d_trace", sid, ch, sid, ch);
00112       obtitle.Form("Trace %2d  %2d", sid, ch);
00113       fSpadic_trace[sid][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE, "time (40ns)");
00114       
00115       obname.Form("Spadic_%d/Channel%d/Spadic%d_ch%d_trace_clean", sid, ch, sid, ch);
00116       obtitle.Form("Trace_clean %2d  %2d", sid, ch);
00117       fSpadic_trace_clean[sid][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE, "time (40ns)");
00118 
00119       obname.Form("Spadic_%d/Channel%d/Spadic%d_ch%d_ADCdist", sid, ch, sid, ch);
00120       obtitle.Form("ADC value distribution %2d  %2d", sid, ch);
00121       fSpadic_ADCdist[sid][ch] = MakeTH1('I', obname.Data(), obtitle.Data(), 256, 0, 256, "ADC value", "Counts");
00122     }
00123     
00124     obname.Form("Spadic_%d/Spadic%d_Spectrum", sid, sid);
00125     obtitle.Form("Spectrum of signal channels of spadic %d", sid);
00126     fSpadic_spectrum[sid] = MakeTH1('D', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE * NUM_SPADIC_CHA * 255 / 10, 0, SPADIC_TRACE_SIZE * NUM_SPADIC_CHA * 255, "Integrated ADC values of signal channels");
00127 
00128     obname.Form("Spadic_%d/Spadic%d_ClusterWidth", sid, sid);
00129     obtitle.Form("Cluster width of spadic %d", sid);
00130     fSpadic_ClusterWidth[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), 9, -0.5, 8.5, "Cluster width");
00131 
00132     obname.Form("Spadic_%d/Spadic%d_SHAPE", sid, sid);
00133     obtitle.Form("Average pulse shape of spadic %d", sid);
00134     fSpadic_shape[sid] = MakeTH1('D', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE);
00135 
00136     obname.Form("Spadic_%d/Spadic%d_SHAPE2D", sid, sid);
00137     obtitle.Form("Average pulse shape of spadic %d", sid);
00138     fSpadic_shape2D[sid] = MakeTH2('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE,305,-50,255);
00139     fSpadic_shape2D[sid]->SetContour(99);
00140 
00141     obname.Form("Spadic_%d/Spadic%d_SHAPE2Draw", sid, sid);
00142     obtitle.Form("Average pulse shape of spadic %d raw signal", sid);
00143     fSpadic_shape2Draw[sid] = MakeTH2('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE,305,-50,255);
00144     fSpadic_shape2Draw[sid]->SetContour(99);
00145 
00146     obname.Form("Spadic_%d/Spadic%d_SHAPEcnt", sid, sid);
00147     obtitle.Form("Counts for average pulse shape of spadic %d", sid);
00148     fSpadic_shapecnt[sid] = MakeTH1('D', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0, SPADIC_TRACE_SIZE);
00149       
00150        
00151     obname.Form("Spadic_%d/Spadic%d_PadMaxDist", sid, sid);
00152     obtitle.Form("Spadic %d peak height", sid);
00153     fSpadic_PadMaxDist[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), NUM_SPADIC_CHA, 0., NUM_SPADIC_CHA, "Pad with maximum signal");  
00154 
00155        
00156     //peak heightfrom bin
00157     obname.Form("Spadic_%d/Spadic%d_PeakHeight", sid, sid);
00158     obtitle.Form("Spadic %d peak height", sid);
00159     fSpadic_peak[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), 255, 0., 255., "PeakHeight");
00160 
00161     // mean position from window
00162     obname.Form("Spadic_%d/Spadic%d_PeakPos", sid, sid);
00163     obtitle.Form("Spadic %d peak position", sid);
00164     fSpadic_meanpos[sid] = MakeTH1('I', obname.Data(), obtitle.Data(), SPADIC_TRACE_SIZE, 0., SPADIC_TRACE_SIZE, "time (40ns)");
00165          
00166   }
00167   if (GetPicture("SpadicsSignalStatistic")==0) {
00168     TGo4Picture* pic = new TGo4Picture("SpadicsSignalStatistic", "signal statistic");  
00169     pic->SetDivision(2, 2);       
00170     for (sid = 1; sid < MAX_SPADIC; sid++) {
00171       pic->Pic(0, 0)->AddObject(fSpadic_NbOfChOverThreshold[sid]);
00172       pic->Pic(0, 1)->AddObject(fSpadic_NbOfOverflow[sid]);
00173       pic->Pic(1, 0)->AddObject(fSpadic_NoiseDist[sid]);
00174       pic->Pic(1, 1)->AddObject(fSpadic_PedelPos[sid]);
00175     }  
00176     AddPicture(pic,"Onlinemonitor");
00177   }
00178 
00179   // reduced to 12                                                                                                                                   
00180   if (GetPicture("SpadicsSpectra")==0) {
00181     TGo4Picture* pic = new TGo4Picture("SpadicsSpectra", "Integrated ADC values of signal channels");
00182     int numy = 4;
00183     int numx = (fPar->numUsedSpadics - 1) / numy;
00184     if (numx*numy < (fPar->numUsedSpadics - 1)) numx++;
00185     pic->SetDivision(numy, numx);
00186 
00187     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00188       int nx = (uid-1) / numy;
00189       int ny = (uid-1) % numy;
00190          
00191       pic->Pic(ny, nx)->AddObject(fSpadic_spectrum[SpaPosMap[uid]]);
00192     }
00193     AddPicture(pic,"Onlinemonitor");
00194   }
00195 
00196   if (GetPicture("SpadicsTraceD2")==0) {
00197     TGo4Picture* pic = new TGo4Picture("SpadicsTraceD2", "Input: Superposition of all shapes");
00198     int numy = 4;
00199     int numx = (fPar->numUsedSpadics-1) / numy;
00200     if (numx*numy < (fPar->numUsedSpadics - 1))
00201       numx++;
00202     pic->SetDivision(numx, numy);
00203 
00204     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00205       int nx = (uid-1) / numy;
00206       int ny = (uid-1) % numy;
00207          
00208       pic->Pic(nx, ny)->SetAutoScale(kFALSE);
00209       pic->Pic(nx, ny)->SetRangeZ(0, 255);
00210 
00211       pic->Pic(nx, ny)->AddObject(fSpadic_trace2D[SpaPosMap[uid]]);
00212     }
00213     AddPicture(pic,"Onlinemonitor");
00214   }
00215 
00216   if (GetPicture("SpadicsTrace2D")==0) {
00217     TGo4Picture* pic = new TGo4Picture("SpadicsTrace2D", "Input: Superposition of all shapes");
00218     int numy = 4;
00219     int numx = (fPar->numUsedSpadics-1) / numy;
00220     if (numx*numy < (fPar->numUsedSpadics - 1))
00221       numx++;
00222     pic->SetDivision(numy, numx);
00223 
00224     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00225       int nx = (uid-1) / numy;
00226       int ny = (uid-1) % numy;
00227          
00228       pic->Pic(ny, nx)->SetAutoScale(kFALSE);
00229       pic->Pic(ny, nx)->SetRangeZ(0, 255);
00230          
00231       pic->Pic(ny, nx)->AddObject(fSpadic_trace2D[SpaPosMap[uid]]);
00232     }
00233     AddPicture(pic,"Onlinemonitor");
00234   }
00235 
00236   if (GetPicture("SpadicsTrace")==0) {
00237     TGo4Picture* pic = new TGo4Picture("SpadicsTrace", "Input: Superposition of all shapes");
00238     int numy = 4;
00239     int numx = (fPar->numUsedSpadics-1) / numy;
00240     if (numx*numy < (fPar->numUsedSpadics - 1))
00241       numx++;
00242     pic->SetDivision(numy, numx);
00243 
00244     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00245       int nx = (uid-1) / numy;
00246       int ny = (uid-1) % numy;
00247          
00248       pic->Pic(ny, nx)->SetAutoScale(kFALSE);
00249       pic->Pic(ny, nx)->SetRangeY(-20, 255);
00250          
00251       for (ch = 0; ch < NUM_SPADIC_CHA; ch++) 
00252         pic->Pic(ny, nx)->AddObject(fSpadic_trace[SpaPosMap[uid]][ch]);
00253     }
00254     AddPicture(pic,"Onlinemonitor");
00255   }
00256       
00257   if (GetPicture("SpadicsTraceClean")==0) {
00258     TGo4Picture* pic = new TGo4Picture("SpadicsTraceClean", "Output: Superposition of all shapes");
00259     int numy = 4;
00260     int numx = (fPar->numUsedSpadics-1) / numy;
00261     if (numx*numy < (fPar->numUsedSpadics - 1))
00262       numx++;
00263     pic->SetDivision(numy, numx);
00264 
00265     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00266       int nx = (uid-1) / numy;
00267       int ny = (uid-1) % numy;
00268          
00269       pic->Pic(ny, nx)->SetAutoScale(kFALSE);
00270       pic->Pic(ny, nx)->SetRangeY(-20, 255);
00271          
00272       for (ch = 0; ch < NUM_SPADIC_CHA; ch++) 
00273         pic->Pic(ny, nx)->AddObject(fSpadic_trace_clean[SpaPosMap[uid]][ch]);
00274     }
00275     AddPicture(pic,"Onlinemonitor");
00276   }
00277     
00278   if (GetPicture("SpadicsPRF")==0) {
00279     TGo4Picture* pic = new TGo4Picture("SpadicsPRF", "PRF for noise correction quality test");
00280     int numy = 4;
00281     int numx = (fPar->numUsedSpadics-1) / numy;
00282     if (numx*numy < (fPar->numUsedSpadics - 1))
00283       numx++;
00284     pic->SetDivision(numy, numx);
00285 
00286     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00287       int nx = (uid-1) / numy;
00288       int ny = (uid-1) % numy;
00289 
00290       pic->Pic(ny, nx)->AddObject(fSpadic_PRF[SpaPosMap[uid]]);
00291     }
00292     AddPicture(pic,"Onlinemonitor");
00293   }
00294 
00295   if (GetPicture("SpadicsMaxAmplitude")==0) {
00296     TGo4Picture* pic = new TGo4Picture("SpadicsMaxAmplitude", "maximum amplitude");
00297     int numy = 4;
00298     int numx = (fPar->numUsedSpadics-1) / numy;
00299     if (numx*numy < (fPar->numUsedSpadics - 1))
00300       numx++;
00301     pic->SetDivision(numy, numx);
00302 
00303     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00304       int nx = (uid-1) / numy;
00305       int ny = (uid-1) % numy;
00306 
00307       pic->Pic(ny, nx)->AddObject(fSpadic_peak[SpaPosMap[uid]]);
00308     }
00309     AddPicture(pic,"Onlinemonitor");
00310   }
00311 
00312   if (GetPicture("SpadicsHitTime")==0) {
00313     TGo4Picture* pic = new TGo4Picture("SpadicsHitTime", "time of maximum amplitude");
00314     int numy = 4;
00315     int numx = (fPar->numUsedSpadics-1) / numy;
00316     if (numx*numy < (fPar->numUsedSpadics - 1))
00317       numx++;
00318     pic->SetDivision(numy, numx);
00319 
00320     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00321       int nx = (uid-1) / numy;
00322       int ny = (uid-1) % numy;
00323 
00324       pic->Pic(ny, nx)->AddObject(fSpadic_meanpos[SpaPosMap[uid]]);
00325     }
00326     AddPicture(pic,"Onlinemonitor");
00327   }
00328 
00329   if (GetPicture("SpadicsBeamProfile")==0) {
00330     TGo4Picture* pic = new TGo4Picture("SpadicsBeamProfile", "Pad max distribution");
00331     int numy = 4;
00332     int numx = (fPar->numUsedSpadics-1) / numy;
00333     if (numx*numy < (fPar->numUsedSpadics - 1))
00334       numx++;
00335     pic->SetDivision(numy, numx);
00336 
00337     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00338       int nx = (uid-1) / numy;
00339       int ny = (uid-1) % numy;
00340 
00341       pic->Pic(ny, nx)->AddObject(fSpadic_PadMaxDist[SpaPosMap[uid]]);
00342     }
00343     AddPicture(pic,"Onlinemonitor");
00344   }
00345 
00346   if (GetPicture("SpadicsClusterWidth")==0) {
00347     TGo4Picture* pic = new TGo4Picture("SpadicsClusterWidth", "Cluster width distribution");
00348     int numy = 4;
00349     int numx = (fPar->numUsedSpadics-1) / numy;
00350     if (numx*numy < (fPar->numUsedSpadics - 1))
00351       numx++;
00352     pic->SetDivision(numy, numx);
00353 
00354     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00355       int nx = (uid-1) / numy;
00356       int ny = (uid-1) % numy;
00357 
00358       pic->Pic(ny, nx)->AddObject(fSpadic_ClusterWidth[SpaPosMap[uid]]);
00359     }
00360     AddPicture(pic,"Onlinemonitor");
00361   }
00362 
00363   if (GetPicture("SpadicsNbOfOverflow")==0) {
00364     TGo4Picture* pic = new TGo4Picture("SpadicsNbOfOverflow", "Number of events in overflow per channel");
00365     int numy = 4;
00366     int numx = (fPar->numUsedSpadics-1) / numy;
00367     if (numx*numy < (fPar->numUsedSpadics - 1))
00368       numx++;
00369     pic->SetDivision(numy, numx);
00370 
00371     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00372       int nx = (uid-1) / numy;
00373       int ny = (uid-1) % numy;
00374 
00375       pic->Pic(ny, nx)->AddObject(fSpadic_NbOfOverflow[SpaPosMap[uid]]);
00376     }
00377     AddPicture(pic,"Onlinemonitor");
00378   }
00379 
00380   if (GetPicture("SpadicsNbOfChOverThreshold")==0) {
00381     TGo4Picture* pic = new TGo4Picture("SpadicsNbOfChOverThreshold", "Number of tbins over threshold per channel (Occupancy)");
00382     int numy = 4;
00383     int numx = (fPar->numUsedSpadics-1) / numy;
00384     if (numx*numy < (fPar->numUsedSpadics - 1))
00385       numx++;
00386     pic->SetDivision(numy, numx);
00387 
00388     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00389       int nx = (uid-1) / numy;
00390       int ny = (uid-1) % numy;
00391 
00392       pic->Pic(ny, nx)->AddObject(fSpadic_NbOfChOverThreshold[SpaPosMap[uid]]);
00393     }
00394     AddPicture(pic,"Onlinemonitor");
00395   }
00396 
00397   if (GetPicture("SpadicsADCdist")==0) {
00398     TGo4Picture* pic = new TGo4Picture("SpadicsADCdist", "Distribution of ADC values, for noise analysis");
00399     int numy = 4;
00400     int numx = (fPar->numUsedSpadics-1) / numy;
00401     if (numx*numy < (fPar->numUsedSpadics - 1))
00402       numx++;
00403     pic->SetDivision(numy, numx);
00404 
00405     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00406       int nx = (uid-1) / numy;
00407       int ny = (uid-1) % numy;
00408 
00409       for (ch = 0; ch < NUM_SPADIC_CHA; ch++) 
00410         pic->Pic(ny, nx)->AddObject(fSpadic_ADCdist[SpaPosMap[uid]][ch]);
00411     }
00412     AddPicture(pic,"Onlinemonitor");
00413   }
00414 
00415   if (GetPicture("SpadicsNoiseDist")==0) {
00416     TGo4Picture* pic = new TGo4Picture("SpadicsNoiseDist", "Noise-Distribution (RMS) of ADC distribution");
00417     int numy = 4;
00418     int numx = (fPar->numUsedSpadics-1) / numy;
00419     if (numx*numy < (fPar->numUsedSpadics - 1))
00420       numx++;
00421     pic->SetDivision(numy, numx);
00422 
00423     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00424       int nx = (uid-1) / numy;
00425       int ny = (uid-1) % numy;
00426 
00427       pic->Pic(ny, nx)->AddObject(fSpadic_NoiseDist[SpaPosMap[uid]]);
00428     }
00429     AddPicture(pic,"Onlinemonitor");
00430   }
00431 
00432   if (GetPicture("SpadicsNoiseDist2D")==0) {
00433     TGo4Picture* pic = new TGo4Picture("SpadicsNoiseDist2D", "Noise-Distribution (RMS) of ADC distribution 2D");
00434     int numy = 4;
00435     int numx = (fPar->numUsedSpadics-1) / numy;
00436     if (numx*numy < (fPar->numUsedSpadics - 1))
00437       numx++;
00438     pic->SetDivision(numy, numx);
00439 
00440     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00441       int nx = (uid-1) / numy;
00442       int ny = (uid-1) % numy;
00443 
00444       pic->Pic(ny, nx)->AddObject(fSpadic_NoiseDist2D[SpaPosMap[uid]]);
00445     }
00446     AddPicture(pic,"Onlinemonitor");
00447   }
00448 
00449   if (GetPicture("SpadicsPedelPos")==0) {
00450     TGo4Picture* pic = new TGo4Picture("SpadicsPedelPos", "Pedestel Position (Mean) of ADC distribution");
00451     int numy = 4;
00452     int numx = (fPar->numUsedSpadics-1) / numy;
00453     if (numx*numy < (fPar->numUsedSpadics - 1))
00454       numx++;
00455     pic->SetDivision(numy, numx);
00456 
00457     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00458       int nx = (uid-1) / numy;
00459       int ny = (uid-1) % numy;
00460 
00461       pic->Pic(ny, nx)->AddObject(fSpadic_PedelPos[SpaPosMap[uid]]);
00462     }
00463     AddPicture(pic,"Onlinemonitor");
00464   }
00465 
00466   if (GetPicture("SpadicsPedelPos2D")==0) {
00467     TGo4Picture* pic = new TGo4Picture("SpadicsPedelPos2D", "Pedestel Position (Mean) of ADC distribution 2D");
00468     int numy = 4;
00469     int numx = (fPar->numUsedSpadics-1) / numy;
00470     if (numx*numy < (fPar->numUsedSpadics - 1))
00471       numx++;
00472     pic->SetDivision(numy, numx);
00473 
00474     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00475       int nx = (uid-1) / numy;
00476       int ny = (uid-1) % numy;
00477 
00478       pic->Pic(ny, nx)->AddObject(fSpadic_PedelPos2D[SpaPosMap[uid]]);
00479     }
00480     AddPicture(pic,"Onlinemonitor");
00481   }
00482   
00483   if (GetPicture("SpadicsShape")==0) {
00484     TGo4Picture* pic = new TGo4Picture("SpadicsShape", "maximum signal shape");
00485     int numy = 4;
00486     int numx = (fPar->numUsedSpadics-1) / numy;
00487     if (numx*numy < (fPar->numUsedSpadics - 1))
00488       numx++;
00489     pic->SetDivision(numy, numx);
00490 
00491     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00492       int nx = (uid-1) / numy;
00493       int ny = (uid-1) % numy;
00494 
00495       pic->Pic(ny, nx)->SetAutoScale(kFALSE);
00496       pic->Pic(ny, nx)->SetRangeY(-20, 255);
00497       pic->Pic(ny, nx)->AddObject(fSpadic_shape[SpaPosMap[uid]]);
00498     }
00499     AddPicture(pic,"Onlinemonitor");
00500   }
00501   
00502   if (GetPicture("SpadicsShape2D")==0) {
00503     TGo4Picture* pic = new TGo4Picture("SpadicsShape2D", "maximum signal shape 2D");
00504     int numy = 4;
00505     int numx = (fPar->numUsedSpadics-1) / numy;
00506     if (numx*numy < (fPar->numUsedSpadics - 1))
00507       numx++;
00508     pic->SetDivision(numy, numx);
00509 
00510     for (Int_t uid = 1; uid < fPar->numUsedSpadics; uid++) {
00511       int nx = (uid-1) / numy;
00512       int ny = (uid-1) % numy;
00513 
00514       pic->Pic(ny, nx)->SetAutoScale(kFALSE);
00515       pic->Pic(ny, nx)->SetRangeY(-50, 255);
00516       pic->Pic(ny, nx)->AddObject(fSpadic_shape2D[SpaPosMap[uid]]);
00517     }
00518     AddPicture(pic,"Onlinemonitor");
00519   }
00520 
00521   TGo4Log::Info("TSpadicProc: Histograms created");
00522 }
00523 
00524 
00525   TSpadicProc::~TSpadicProc()
00526   {
00527   }
00528 
00529   void TSpadicProc::InitEvent(TGo4EventElement* outevnt){
00530     // test if we are in beamtime or standalone mode:
00531     // since output event object is never discarded within processor lifetime,
00532     // we just search for subevent by name once to speed up processing  
00533     if(fOutputEvent==0){
00534         TCBMBeamtimeEvent* btevent=dynamic_cast<TCBMBeamtimeEvent*>(outevnt);
00535         if(btevent){
00536             fOutputEvent=dynamic_cast<TSpadicEvent*>(btevent->GetSubEvent("SPADIC"));
00537           }     
00538         else{
00539 
00540             fOutputEvent= dynamic_cast<TSpadicEvent*>(outevnt);
00541           }
00542         if(fOutputEvent==0){
00543           GO4_STOP_ANALYSIS_MESSAGE("**** TSpadicProc: Fatal error: output event is not a TSpadicEvent!!! STOP GO4");
00544         }
00545 
00546       }
00547   }
00548 
00549 
00550 #ifdef WITHOUT_SPADIC
00551 
00552   void TSpadicProc::ProcessSubevent(TGo4MbsSubEvent* subevt){
00553     // do nothing
00554   }
00555 
00556 
00557   void TSpadicProc::ProcessSpadic(TGo4MbsSubEvent* psubevt){
00558     // do nothing
00559   }
00560 
00561 #else
00562 
00563 
00564   void TSpadicProc::ProcessSubevent(TGo4MbsSubEvent* subevt){
00565 #if __GO4BUILDVERSION__ > 40500
00566           Bool_t oldevent=IsKeepInputEvent();
00567           SetKeepInputEvent(kFALSE); // need this to reset our flag in case rocproc had requested it.
00568           if(oldevent) return; // ignore this subevent if we already had it
00569 #endif
00570 
00571 
00572           if (subevt->GetProcid() == roc::proc_TRD_Spadic){
00573 
00574                   Bool_t offSpill = kFALSE;
00575                   if(GetTriggerNumber()==8)
00576                   {
00577                           offSpill=kTRUE; // pulser events are flaged and can be ignored within 'ProcessSpadic(subevt,offSpill)'
00578                           //ignore pulser events in spill pause for the moment
00579                           //fOutputEvent->SetValid(kFALSE);
00580                           //return;
00581                   }
00582 
00583                   // TRD Spadic data
00584                   fDoStopAnalysis=kFALSE;
00585                   ProcessSpadic(subevt,offSpill);
00586                   if (fDoStopAnalysis){
00587 
00588                           for (UInt_t sid = 1; sid < MAX_SPADIC; sid++) {
00589                                   if (TSpadicEvent::ConfigSpadics[sid] == 0) continue;
00590 
00591                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_trace2D[sid]);
00592                                   for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++) {
00593                                           TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_trace[sid][ch]);
00594                                           TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_trace_clean[sid][ch]);
00595                                           TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_ADCdist[sid][ch]);
00596                                   }
00597                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_NoiseDist[sid]);
00598                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_PedelPos[sid]);
00599                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_NoiseDist2D[sid]);
00600                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_PedelPos2D[sid]);
00601                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_NbOfChOverThreshold[sid]);
00602                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_NbOfOverflow[sid]);
00603 
00604                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_PRF[sid]);
00605                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_shape[sid]);
00606                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_shape2D[sid]);
00607           TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_shape2Draw[sid]);
00608                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_meanpos[sid]);
00609                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_peak[sid]);
00610                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_ClusterWidth[sid]);
00611                                   TGo4Analysis::Instance()->SendObjectToGUI(fSpadic_PadMaxDist[sid]);
00612                           }
00613                           TGo4Analysis::Instance()->SetRunning(kFALSE);
00614                   }
00615                   return;
00616           }
00617   }
00618 
00619 
00620 void TSpadicProc::ProcessSpadic(TGo4MbsSubEvent* psubevt, Bool_t offSpill){
00621   
00622   UInt_t crateid=psubevt->GetSubcrate();  
00623   UInt_t sid =crateid;
00624   /*
00625     if(fPar->idMappingEnabled)
00626     {
00627     for(int i=0;i<fPar->numSpadics;++i)
00628     {
00629     if(crateid==fPar->idSpadic[i])
00630     {
00631     sid= i;
00632     break;
00633     }
00634     }
00635     }
00636   */
00637   if (sid >= MAX_SPADIC) {
00638     TGo4Log::Info("TSpadicProc: Skipping subevent with invalid Susibo id %u,"
00639                   " maximum allowed is %u", sid, MAX_SPADIC - 1);
00640     return;
00641   }
00642   TSpadicData* spdata=dynamic_cast<TSpadicData*>(fOutputEvent->getEventElement(sid));
00643   if(spdata==0){
00644      TGo4Log::Info("TSpadicProc: Skipping subevent with not configured Susibo id %u,"
00645                         ", please check spadic outputevent setup!", sid);
00646     return;
00647   }
00648   spdata->SetValid(kTRUE); // mark as existing in data stream
00649 
00650   spadic::Message mess((uint8_t*) psubevt->GetDataField());
00651   if (!mess.CheckMessage()) {
00652     TGo4Log::Info("TSpadicProc: Skipping subevent with invalid message, sid: %u", sid);
00653     return;
00654   }
00655   if (SPADIC_DEBUG) {
00656     cout << " - Message sid:" << sid << ", status:" << mess.GetStatusNumber();
00657     cout << ", evid:" << mess.GetEventIDNumber() << ", ts:" << mess.GetTimeStamp() << endl;
00658   }
00659 
00660   Double_t AverTrace[SPADIC_TRACE_SIZE] = {0.0};
00661   spdata->fSpadicOffSpillEvent=offSpill;
00662   for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00663     // generate sample trace histograms:
00664 
00665     //TH1* tracehis = fSpadic_trace[sid][ch];                                        // the original sample
00666     for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){
00667       uint8_t val = mess.Sample(ch, bin);
00668       AverTrace[bin] += val;
00669       //tracehis->SetBinContent(bin + 1, (int) val);
00670 
00671       if (fSpadic_trace[sid][ch] == 0) printf("Zero trace for spadic %d\n", sid);
00672       fSpadic_trace[sid][ch]->SetBinContent(bin + 1, (int) val);             // // // copy of original sampels, this spectra will be corrected
00673       fSpadic_trace2D[sid]->SetBinContent( (int)bin + 1,  (int)ch + 1, (int) val); 
00674       spdata->fSpadicPulse[ch][bin]=val;
00675       if (val >= 255) {
00676         spdata->fSpadicOverflows[ch][bin] = true;  
00677         if (!spdata->fSpadicOverFlowEvent)    
00678           spdata->fSpadicOverFlowEvent = true;  
00679       }                       // // // mark overflowting time bins
00680       if (val <= 0)   {
00681         spdata->fSpadicUndershoots[ch][bin] = true;
00682         if (!spdata->fSpadicUndershootsEvent)
00683           spdata->fSpadicUndershootsEvent = true;
00684       } 
00685     } // for bin
00686   }// for ch 
00687 
00688   if (fPar->Run_SimpleAnalysis){
00689     Float_t AveragePerBin[SPADIC_TRACE_SIZE] = {0.0};
00690     Float_t Pedestal[NUM_SPADIC_CHA] = {0.0};
00691     Float_t IntegratedSignal = 0.0;
00692     
00693     for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00694       for (UInt_t bin = 0; bin < 5; bin++){
00695         Pedestal[ch] += (uint8_t)mess.Sample(ch, bin);
00696       }
00697       Pedestal[ch] /= 5.0;
00698     }
00699     for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE; bin++){
00700       for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00701         AveragePerBin[bin] += (uint8_t)mess.Sample(ch, bin) - Pedestal[ch];
00702       }
00703       AveragePerBin[bin] /= Float_t(NUM_SPADIC_CHA);
00704     }
00705     Int_t     maxch = -1; 
00706     Float_t maxA = 0;  
00707     for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00708       for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE; bin++){
00709         spdata->fSpadicCompensated[ch][bin] = (uint8_t)mess.Sample(ch, bin) - Pedestal[ch] - AveragePerBin[bin];
00710         fSpadic_trace_clean[sid][ch]->SetBinContent(bin + 1, spdata->fSpadicCompensated[ch][bin]); 
00711         IntegratedSignal = spdata->fSpadicCompensated[ch][bin];
00712         if (spdata->fSpadicCompensated[ch][bin] > maxA){
00713           maxA = spdata->fSpadicCompensated[ch][bin];
00714           maxch = ch;
00715         }
00716       }
00717     }
00718     spdata->fSpadicHighestChannel = maxch;
00719     fSpadic_spectrum[sid]->Fill(IntegratedSignal);
00720     if (maxch > 0)
00721       spdata->fSpadicSignalCh[maxch-1] = true;
00722 
00723     spdata->fSpadicSignalCh[maxch] = true;
00724 
00725     if (maxch < NUM_SPADIC_CHA-1)
00726       spdata->fSpadicSignalCh[maxch+1] = true;
00727   }
00728 
00729   else {
00730     // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Find Signal Channels
00731     for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++) 
00732       AverTrace[bin] /= Float_t(NUM_SPADIC_CHA);
00733 
00734     Double_t maxsum = 0;                      // absolut integral maximum
00735     Int_t     maxch = -1;                     // channel with maximum signal integral
00736     Double_t sumCh[NUM_SPADIC_CHA] = {0};     // signal integral per channel
00737     Int_t     maxT[NUM_SPADIC_CHA] = {-1};    // maximum amplitude time bin
00738     Double_t  maxA[NUM_SPADIC_CHA] = {-100};  // maximum amplitude
00739     Double_t firstTb[NUM_SPADIC_CHA] = {0};   // baseline reference only for clustering
00740 
00741     for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00742       Double_t sum = 0;
00743       firstTb[ch] = (mess.Sample(ch, 0) - AverTrace[0]);
00744 
00745       for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){ 
00746         sum += (mess.Sample(ch, bin) - AverTrace[bin]);
00747    
00748         if ((mess.Sample(ch, bin) - AverTrace[bin]) > maxA[ch]){
00749           maxA[ch] = mess.Sample(ch, bin) - AverTrace[bin];
00750           maxT[ch] = bin;
00751         }
00752       }
00753       sumCh[ch] = sum;
00754       if (sum > maxsum) {
00755         maxsum = sum;
00756         maxch = ch;
00757       }
00758     }
00759     spdata->fSpadicHighestChannel = maxch;
00760     if (maxch == -1){
00761       printf(" <<<< no signal found for spadic %d\n",sid);
00762     }
00763     else {
00764       fSpadic_PadMaxDist[sid]->Fill(maxch);
00765       fSpadic_meanpos[sid]->Fill(maxT[maxch]);
00766       fSpadic_peak[sid]->Fill(maxA[maxch]);
00767       // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
00768       // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<   CLUSTERIZER
00769       Int_t nSpadicSignalCh = 0;
00770       if (fPar->Run_TimeClustering) {
00771         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00772           if (maxT[ch] >= maxT[maxch] - 5 && maxT[ch] <= maxT[maxch] + 5){ // if signal is correlated in time domain
00773             spdata->fSpadicSignalCh[ch] = true;
00774             nSpadicSignalCh++;
00775           }
00776         }
00777       }
00778       if (fPar->Run_AmplitudeClustering) {
00779         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00780           if (maxA[ch] > 15 + firstTb[ch]){                               //   if signal is above a minimum threshold with respect to the first time bin of this channel
00781             spdata->fSpadicSignalCh[ch] = true;
00782             nSpadicSignalCh++;    
00783           }
00784         }
00785       }
00786       fSpadic_ClusterWidth[sid]->Fill(nSpadicSignalCh);   
00787       // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
00789       if (fPar->Run_RawSignalStatistics){
00790         Float_t pedestal[NUM_SPADIC_CHA] = {0.0};
00791 
00792         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00793           if (fSpadic_trace[sid][ch]->GetMaximum() >= 255) 
00794             fSpadic_NbOfOverflow[sid]->Fill(ch);
00795 
00796           if (fSpadic_trace[sid][ch]->GetMaximum() > fPar->HitThreshold) 
00797             fSpadic_NbOfChOverThreshold[sid]->Fill(ch);
00798       
00799           //if (fSpadic_trace[sid][ch]->GetMaximum() < fPar->HitThreshold){                          //// fill RMS ditribution only with not hit channels
00800           if (!spdata->fSpadicSignalCh[ch]){
00801             for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){
00802               pedestal[ch] += fSpadic_trace[sid][ch]->GetBinContent(bin);
00803               fSpadic_ADCdist[sid][ch]->Fill(fSpadic_trace[sid][ch]->GetBinContent(bin));
00804             }
00805             pedestal[ch] /= SPADIC_TRACE_SIZE;
00806             fSpadic_PedelPos2D[sid]->Fill(ch,pedestal[ch]);
00807           } 
00808 
00809           fSpadic_NoiseDist[sid]->SetBinContent(ch+1, fSpadic_ADCdist[sid][ch]->GetRMS());          
00810           fSpadic_NoiseDist2D[sid]->Fill(ch, fSpadic_ADCdist[sid][ch]->GetRMS());
00811           fSpadic_PedelPos[sid]->SetBinContent(ch+1, fSpadic_ADCdist[sid][ch]->GetMean());          
00812           //fSpadic_PedelPos2D[sid]->Fill(ch, fSpadic_ADCdist[sid][ch]->GetMean());
00813           //}
00814         }
00815       }
00817       // <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<   Noise Compensation
00818       if (fPar->Run_SimpleNoiseReduction) {
00819         Int_t MaxCharge = 0;
00820         for (UInt_t bin = 5; bin < 20; bin++){
00821           for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00822             if(mess.Sample(ch, bin) > MaxCharge) {
00823               MaxCharge = mess.Sample(ch, bin);
00824               maxch = ch;
00825             }
00826           }
00827         }
00828         for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){ 
00829           if (maxch != 0 && maxch != NUM_SPADIC_CHA-1) 
00830             AverTrace[bin] = (mess.Sample(0, bin) + mess.Sample(NUM_SPADIC_CHA-1, bin)) * 0.5;
00831           if (maxch == NUM_SPADIC_CHA-1)
00832             AverTrace[bin] = mess.Sample(0, bin);
00833           if (maxch == 0)
00834             AverTrace[bin] = mess.Sample(NUM_SPADIC_CHA-1, bin);
00835         }
00836       }
00837   
00838       else {
00839         Int_t nOBins[NUM_SPADIC_CHA] = {0};
00840         Double_t AverTraceO[NUM_SPADIC_CHA] = {0.0};
00841         for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){ 
00842           AverTrace[bin] = 0.0;  
00843           for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){
00844             if (!spdata->fSpadicSignalCh[ch]){
00845               AverTrace[bin] += mess.Sample(ch, bin);
00846               if(spdata->fSpadicOverflows[ch][bin]){
00847                 nOBins[ch]++;
00848                 AverTraceO[ch] += mess.Sample(ch, bin);
00849               }
00850             }
00851           }
00852           if(nSpadicSignalCh<NUM_SPADIC_CHA) // Debug check validity of array? JAM
00853             {
00854               AverTrace[bin] /= (NUM_SPADIC_CHA - nSpadicSignalCh);
00855             }
00856           else
00857             {
00858               //cout <<"WWWWWWWWWW arning: nSpadicSignalCh="<<nSpadicSignalCh<<" exceeds NUM_SPADIC_CHA "<<NUM_SPADIC_CHA << endl;
00859             }
00860 
00861         }
00862 
00863         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++)
00864           AverTraceO[ch] /= nOBins[ch];
00865       }
00866       Double_t pedestal[NUM_SPADIC_CHA] = {0.0};
00867       //Int_t pedestalTB = 3;
00868       Double_t intensCh[NUM_SPADIC_CHA] = {0.0};
00869       if (fPar->Run_PedestleCorrection) {
00870         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){ 
00871           for (Int_t bin = 0; bin < fPar->pedestalTB; bin++){ 
00872             pedestal[ch] += mess.Sample(ch, bin) - AverTrace[bin];
00873           }
00874           if (fPar->pedestalTB > 0 && fPar->pedestalTB < SPADIC_TRACE_SIZE)
00875             pedestal[ch] /= fPar->pedestalTB;
00876           else
00877             pedestal[ch] = 0;                                          
00878         }
00879       }
00880 
00881       for (UInt_t bin = 0; bin < SPADIC_TRACE_SIZE-1; bin++){ 
00882         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){ 
00883           spdata->fSpadicCompensated[ch][bin] = mess.Sample(ch, bin);
00884 
00885           if(!spdata->fSpadicOverflows[ch][bin]){
00886             spdata->fSpadicCompensated[ch][bin] -= AverTrace[bin] - pedestal[ch];
00887           }
00888           else {
00889             //spdata->fSpadicCompensated[ch][bin] -= AverTraceO[ch];
00890             spdata->fSpadicCompensated[ch][bin] -= AverTrace[bin] - pedestal[ch];
00891           }
00892 
00893           fSpadic_trace_clean[sid][ch]->SetBinContent(bin + 1, spdata->fSpadicCompensated[ch][bin]); 
00894           intensCh[ch] += spdata->fSpadicCompensated[ch][bin];
00895           if (maxch>=0 && (ch == (unsigned) maxch)){
00896             double val = spdata->fSpadicCompensated[ch][bin];
00897             double weight = fSpadic_shapecnt[sid]->GetBinContent(bin+1);
00898             double value0 = fSpadic_shape[sid]->GetBinContent(bin+1);
00899             //     printf("ch:%d Val:%e  weight:%e, value0:%e \n",ch,val,weight,value0);
00900             //     cout <<"SpadicCompensated:"<<spdata->fSpadicCompensated[ch][bin] << endl;
00901             //     cout <<"AverTrace:"<<AverTrace[bin] << endl;
00902             //     cout <<"pedestal:"<<pedestal[ch] << endl;
00903             //    printf("Setting fSpadic_shape sid:%d, ptr:0x%p name:%s bin:%d val=%e\n",sid,fSpadic_shape[sid],fSpadic_shape[sid]->GetName(),bin,(value0 * weight + val) / (weight + 1.));
00904             //
00905 
00906       fSpadic_shape2Draw[sid]->Fill(bin,fSpadic_trace[sid][ch]->GetBinContent(bin));  // AATEST
00907       
00908             fSpadic_shape2D[sid]->Fill(bin, spdata->fSpadicCompensated[ch][bin]);
00909             fSpadic_shape[sid]->SetBinContent(bin+1, (value0 * weight + val) / (weight + 1.));
00910             fSpadic_shapecnt[sid]->SetBinContent(bin+1, weight+1.);
00911           }
00912         }
00913       }
00914       Double_t intens = 0;
00915       if(fPar->UseClusteredSignal)
00916         for (UInt_t ch = 0; ch < NUM_SPADIC_CHA; ch++){ 
00917           if(spdata->fSpadicSignalCh[ch])
00918             intens += intensCh[ch];
00919         }
00920       else
00921         intens = intensCh[maxch-1] + intensCh[maxch] + intensCh[maxch+1];
00922 
00923       fSpadic_spectrum[sid]->Fill(intens);
00924       // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
00925       // >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> PRF  
00926       if (fPar->Run_PadResponseFunction) {
00927         if (maxch >= 1 && maxch <= NUM_SPADIC_CHA - 1){
00928           Double_t dxPos = 0.5 * log((intensCh[maxch+1] / intensCh[maxch-1])) / log((pow(intensCh[maxch],2) / (intensCh[maxch+1] * intensCh[maxch-1])));
00929           Double_t chargePercent = intensCh[maxch-1] / (intensCh[maxch-1] + intensCh[maxch] + intensCh[maxch+1]);
00930           fSpadic_PRF[sid]->Fill(-dxPos-1,chargePercent);
00931           //PRF1D->Fill(-dxPos-pWidth,chargePercent);
00932           chargePercent = intensCh[maxch] / (intensCh[maxch-1] + intensCh[maxch] + intensCh[maxch+1]);
00933           fSpadic_PRF[sid]->Fill(dxPos,chargePercent);
00934           //PRF1D->Fill(dxPos,chargePercent);
00935           chargePercent = intensCh[maxch+1] / (intensCh[maxch-1] + intensCh[maxch] + intensCh[maxch+1]);
00936           fSpadic_PRF[sid]->Fill(-dxPos+1,chargePercent);
00937           //PRF1D->Fill(-dxPos+pWidth,chargePercent);
00938         } 
00939       }
00940     }
00941   }
00942 }
00943 #endif
00944 
00945 

Generated on Tue Dec 10 2013 04:52:23 for ROCsoft by  doxygen 1.7.1