#include <stdio.h>
#include <stdlib.h>
#include <string>
#include <map>
#include <utility>
#include <sstream>
#include "midas.h"
#include <TH2.h>
#include <TDirectory.h>
#include "TGlobalData.h"
#include "TSetupData.h"
#include "TPulseIsland.h"
Go to the source code of this file.
Functions | |
INT | MFastSlowCorrelator_init (void) |
INT | MFastSlowCorrelator (EVENT_HEADER *, void *) |
int | FillCoincidenceTimeHistogram (double tFast, double aFast, vector< TPulseIsland * > pulses, int start, TH1 *h, TH2 *hAmp, double t_low, double t_high, bool PP, TH1 *h_PP, TH2 *hAmp_PP, TH1 *hSC) |
bool | IsPileupProtected (double t, std::vector< TPulseIsland * > p, double window, int *lPP) |
Variables | |
HNDLE | hDB |
TGlobalData * | gData |
Object to hold data used and produced by modules throughout alcapana stage of analysis. | |
TSetupData * | gSetup |
Hardware information about digitizers and detectors to be used during alcapana stage of analysis. | |
ANA_MODULE | MFastSlowCorrelator_module |
static TH1 ** | hTime |
static TH2 ** | hAmplitudeCorrelation |
static TH1 ** | hSlowCounts |
static TH1 ** | hTime_PP |
static TH2 ** | hAmplitudeCorrelation_PP |
static int * | last |
static std::vector< std::string > | FastDetectors |
static std::vector< std::string > | SlowDetectors |
static int | nr_fast_detectors |
static int | nr_slow_detectors |
static int | nr_detectors |
int FillCoincidenceTimeHistogram | ( | double | tFast, | |
double | aFast, | |||
std::vector< TPulseIsland * > | pulses, | |||
int | start, | |||
TH1 * | h, | |||
TH2 * | hAmp, | |||
double | t_low, | |||
double | t_high, | |||
bool | PP, | |||
TH1 * | h_PP, | |||
TH2 * | hAmp_PP, | |||
TH1 * | hSC | |||
) |
Definition at line 216 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MMuSCTimeDifferences().
00219 { 00220 int i=0; 00221 for(int p = start; p < (int)pulses.size(); ++p){ 00222 00223 double tDiff = pulses[p]->GetPulseTime() - tFast; 00224 if(tDiff < t_low){ 00225 start = p; 00226 continue; 00227 } 00228 if(tDiff > t_high) break; 00229 h->Fill(tDiff); 00230 hAmp->Fill(aFast, pulses[p]->GetPulseHeight()); 00231 i++; 00232 if(PP){ 00233 h_PP->Fill(tDiff); 00234 hAmp_PP->Fill(aFast, pulses[p]->GetPulseHeight()); 00235 } 00236 } 00237 hSC->Fill(i); 00238 return start; 00239 }
bool IsPileupProtected | ( | double | t, | |
std::vector< TPulseIsland * > | p, | |||
double | window, | |||
int * | lPP | |||
) |
Definition at line 241 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MMuSCTimeDifferences().
00241 { 00242 if(t < window || t > 1.12E8-window) return false; // very basic bookending 00243 for(int i = *lPP; i < p.size(); ++i){ 00244 double t2 = p[i]->GetPulseTime(); 00245 if(t2 < t-window){ 00246 *lPP = i; 00247 } 00248 if(t2 > t+window) break; 00249 if(t != t2 && t2 > t-window && t2 < t+window) return false; 00250 } 00251 return true; 00252 }
INT MFastSlowCorrelator | ( | EVENT_HEADER * | pheader, | |
void * | pevent | |||
) |
This method processes one MIDAS block, producing a vector of TOctalFADCIsland objects from the raw Octal FADC data.
Definition at line 178 of file MFastSlowCorrelator.cpp.
References FastDetectors, FillCoincidenceTimeHistogram(), TGlobalData::fPulseIslandToChannelMap, TSetupData::GetBankName(), hAmplitudeCorrelation, hAmplitudeCorrelation_PP, hSlowCounts, hTime, hTime_PP, IsPileupProtected(), last, nr_detectors, and SlowDetectors.
00178 { 00179 double t_low = -10000; 00180 double t_high = 10000; 00181 00182 // Some typedefs 00183 typedef map<string, vector<TPulseIsland*> > TStringPulseIslandMap; 00184 00185 // Fetch a reference to the gData structure that stores a map 00186 // of (bank_name, vector<TPulseIsland*>) pairs 00187 TStringPulseIslandMap& pulse_islands_map = 00188 gData->fPulseIslandToChannelMap; 00189 00190 // Let's loop over fast channel times 00191 int lPP = 0; 00192 00193 for(int i=0; i<nr_detectors; ++i) last[i] = 0; 00194 for (int i=0; i<nr_detectors; ++i){ 00195 std::string FastBank = gSetup->GetBankName(FastDetectors[i]); 00196 std::string SlowBank = gSetup->GetBankName(SlowDetectors[i]); 00197 lPP = 0; 00198 for(std::vector<TPulseIsland*>::iterator it_Fast = pulse_islands_map[FastBank].begin(); 00199 it_Fast != pulse_islands_map[FastBank].end(); ++it_Fast){ 00200 00201 double tFast = (*it_Fast)->GetPulseTime(); 00202 double aFast = (*it_Fast)->GetPulseHeight(); 00203 00204 bool PP = IsPileupProtected(tFast, pulse_islands_map[FastBank], 10000, &lPP); 00205 00206 last[i] = FillCoincidenceTimeHistogram(tFast,aFast,pulse_islands_map[SlowBank], 00207 last[i], hTime[i], hAmplitudeCorrelation[i], t_low, t_high, 00208 PP, hTime_PP[i], hAmplitudeCorrelation_PP[i],hSlowCounts[i]); 00209 } 00210 } 00211 00212 00213 return SUCCESS; 00214 }
INT MFastSlowCorrelator_init | ( | ) |
This method initializes histograms.
Definition at line 78 of file MFastSlowCorrelator.cpp.
References FastDetectors, TSetupData::fBankToDetectorMap, hAmplitudeCorrelation, hAmplitudeCorrelation_PP, hSlowCounts, hTime, hTime_PP, last, nr_detectors, nr_fast_detectors, nr_slow_detectors, and SlowDetectors.
00079 { 00080 gDirectory->mkdir("FastSlowCorrelations")->cd(); 00081 // Let's pair fast and slow channels 00082 for(std::map<std::string, std::string>::iterator it = gSetup->fBankToDetectorMap.begin(); 00083 it != gSetup->fBankToDetectorMap.end(); ++it){ 00084 if(it->first != std::string("") && it->first != std::string("blank")){ 00085 if (it->second.compare("SiR1-sum-F")==0) continue; 00086 size_t StrPos = it->second.find("-F"); 00087 if(StrPos!=std::string::npos){ 00088 FastDetectors.push_back(it->second); 00089 // printf("%s\n",it->second.c_str()); 00090 std::string paired_slow = it->second.substr(0,StrPos); 00091 paired_slow+="-S"; 00092 for(std::map<std::string, std::string>::iterator it_slow = gSetup->fBankToDetectorMap.begin(); 00093 it_slow != gSetup->fBankToDetectorMap.end(); ++it_slow){ 00094 if(it_slow->first != std::string("") && it_slow->first != std::string("blank")){ 00095 if(it_slow->second.compare(paired_slow)==0){ 00096 SlowDetectors.push_back(it_slow->second); 00097 // printf("%s\n",it_slow->second.c_str()); 00098 } 00099 } 00100 } 00101 } 00102 } 00103 } 00104 00105 nr_fast_detectors = FastDetectors.size(); 00106 nr_slow_detectors = SlowDetectors.size(); 00107 if (nr_slow_detectors != nr_fast_detectors){ 00108 printf("Fast and Slow channel numbers don't match!\n"); 00109 return 1; 00110 } 00111 nr_detectors = nr_slow_detectors; 00112 00113 hTime = new TH1*[nr_detectors]; 00114 hTime_PP = new TH1*[nr_detectors]; 00115 hAmplitudeCorrelation = new TH2*[nr_detectors]; 00116 hAmplitudeCorrelation_PP = new TH2*[nr_detectors]; 00117 hSlowCounts = new TH1*[nr_detectors]; 00118 last = new int[nr_detectors]; 00119 00120 std::string hName; 00121 std::string hTitle; 00122 00123 for(int i = 0; i<nr_detectors; ++i){ 00124 last[i] = 0; 00125 00126 hName = std::string("h_"); 00127 hName += FastDetectors[i]; 00128 hName += "_"; 00129 hName += SlowDetectors[i]; 00130 hName += "_TDiff"; 00131 hTitle = std::string("Tdiff for "); 00132 hTitle += SlowDetectors[i]; 00133 hTitle += " - "; 00134 hTitle += FastDetectors[i]; 00135 hTime[i] = new TH1F(hName.c_str(), hTitle.c_str(), 500, -10000, 10000); 00136 hTime[i]->GetXaxis()->SetTitle("Time difference [ns]"); 00137 hName += "_PP"; 00138 hTitle += " with PP"; 00139 hTime_PP[i] = new TH1F(hName.c_str(), hTitle.c_str(), 500, -10000, 10000); 00140 hTime_PP[i]->GetXaxis()->SetTitle("Time difference [ns]"); 00141 00142 hName = std::string("h_"); 00143 hName += FastDetectors[i]; 00144 hName += "_"; 00145 hName += SlowDetectors[i]; 00146 hName += "Amplitude_Correlation"; 00147 hTitle = std::string("Tdiff for "); 00148 hTitle += SlowDetectors[i]; 00149 hTitle += " - "; 00150 hTitle += FastDetectors[i]; 00151 printf("%s %s\n",hName.c_str(),hTitle.c_str()); 00152 hAmplitudeCorrelation[i] = new TH2F(hName.c_str(), hTitle.c_str(), 250, 0, 1000, 00153 250, 0, 1000); 00154 hAmplitudeCorrelation[i]->GetXaxis()->SetTitle("Fast Amplitude"); 00155 hAmplitudeCorrelation[i]->GetYaxis()->SetTitle("Slow Amplitude"); 00156 hName += "_PP"; 00157 hTitle += " with PP"; 00158 hAmplitudeCorrelation_PP[i] = new TH2F(hName.c_str(), hTitle.c_str(), 250, 0, 1000, 00159 250, 0, 1000); 00160 hAmplitudeCorrelation_PP[i]->GetXaxis()->SetTitle("Fast Amplitude"); 00161 hAmplitudeCorrelation_PP[i]->GetYaxis()->SetTitle("Slow Amplitude"); 00162 00163 hName = std::string("h_"); 00164 hName += "SlowCounts_"; 00165 hName += FastDetectors[i]; 00166 hTitle = std::string("SlowCounts correspond to "); 00167 hTitle += FastDetectors[i]; 00168 hSlowCounts[i] = new TH1I(hName.c_str(),hTitle.c_str(),10,0,10); 00169 } 00170 gDirectory->cd("/"); 00171 printf("Initialization Complete!\n"); 00172 return SUCCESS; 00173 }
std::vector<std::string> FastDetectors [static] |
Definition at line 64 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
TH2** hAmplitudeCorrelation [static] |
Definition at line 59 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
TH2** hAmplitudeCorrelation_PP [static] |
Definition at line 62 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
HNDLE hDB |
TH1** hSlowCounts [static] |
Definition at line 60 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
TH1** hTime [static] |
Definition at line 58 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), MFastSlowCorrelator_init(), and PlotTAP_Time::ProcessEntry().
TH1** hTime_PP [static] |
Definition at line 61 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
int* last [static] |
Definition at line 63 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), MFastSlowCorrelator_init(), and SimpleHistograms::ProcessEntry().
ANA_MODULE MFastSlowCorrelator_module |
{ "MFastSlowCorrelator", "Ran Hong", MFastSlowCorrelator, NULL, NULL, MFastSlowCorrelator_init, NULL, NULL, 0, NULL, }
Definition at line 44 of file MFastSlowCorrelator.cpp.
int nr_detectors [static] |
Definition at line 68 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().
int nr_fast_detectors [static] |
Definition at line 66 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator_init().
int nr_slow_detectors [static] |
Definition at line 67 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator_init().
std::vector<std::string> SlowDetectors [static] |
Definition at line 65 of file MFastSlowCorrelator.cpp.
Referenced by MFastSlowCorrelator(), and MFastSlowCorrelator_init().