GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
Deadtime.cxx
Go to the documentation of this file.
1// To compile:
2// Note: GRSISort looks for .cxx extensions when compiling (for example it looks in the myAnalysis directory)
3// Alternatively you may use the following to compile:
4// g++ myanalysis.cxx -o myanalysis -std=c++0x -I$GRSISYS/GRSISort/include/ `grsi-config --cflags --all-libs --root`
5
6#include <fstream>
7#include <cstdlib>
8#include <cstring>
9#include <algorithm>
10#include <iostream>
11#include <functional>
12#include <iomanip>
13#include <string>
14#include <cmath>
15#include <cstdio>
16#include <ctime>
17
18#include "TF1.h"
19#include "TMath.h"
20#include "TH1.h"
21#include "TH1F.h"
22#include "THStack.h"
23#include "TString.h"
24#include "TCanvas.h"
25#include "TPad.h"
26#include "TFile.h"
27//#include "TFileIter.h"
28#include "TKey.h"
29#include "TTree.h"
30#include "TH2F.h"
31#include "TROOT.h"
32#include "TPPG.h"
33#include "TScaler.h"
34#include "TApplication.h"
35#include "TLeaf.h"
36
37#ifndef __CINT__
38#include "TGriffin.h"
39#include "TSceptar.h"
40#endif
41
42void Printaddress(int* channel);
43void MakeSpectra(const char*& filename, int& prog, const char*& fname, int& nsclr, int& ncycle, double* rate,
44 int* channel, int& index, const int* trun, double& thresh);
45void CheckFile(const char*& fname);
46void DoAnalysis(const char*& fname, int& nfile, double* rate, int& nsclr, int& patlen, int& ncycle, int* trun,
47 double& eor, const char*& hname, const char*& iname, const char*& jname, const char*& kname,
48 const char*& lname, const char*& mname, const char*& nname, int& nscaler);
49
50int main(int argc, char* argv[])
51{
52 TApplication theApp("Analysis", &argc, argv);
53 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~USER
54 //INPUTS
55 std::ifstream filelist("/home/mbowry/Documents/ROOT/Backup/newfilelist_hs.txt"); // input file(s) DATA
56 std::ifstream ODBlist("/home/mbowry/Documents/ROOT/Backup/ODBlist_hs.txt"); // input file(s) ODB (runinfo)
57 const char* fname = "viewscaler.root"; // output file (scaler spectra)
58 const char* hname = "random_counts.txt"; // output file (frequency histograms, source)
59 const char* iname = "combined_counts.txt"; // output file (frequency histograms, source+pulser)
60 const char* jname = "random_seed.txt"; // output file (check random number generation)
61 const char* kname = "RESULTS.txt"; // output file (deadtime results)
62 const char* lname = "asymmetric_error.txt"; // output file (error combination histogram)
63 const char* mname = "random_deadtime.txt"; // output file (random deadtime histogram)
64 const char* nname = "accepted_rand.txt"; // output file (accepted random rate histogram)
65 int nsclr = 9; // total number of pulser inputs
66 std::array<int, 9> addr = {0x0000, 0x000d, 0x0101, 0x0107, 0x020b,
67 0x020c, 0x0302, 0x030c, 0x030d}; // addresses with pulser inputs
68 std::array<double, 4> freq = {2.e3, 5.e3, 1.e4, 2.e4}; // precision pulser rates (2,5,10,20 kHz)
69 int patlen = 10; // pattern length in seconds
70 int ncycle = 1; // read out period of scalers (seconds);
71 int scaleri = 0; // scaler# START
72 int scalerf = 3; // scaler# END (0->3 = all scalers)
73 double thresh = 0.3; // threshold for rejection of spurious scaler events
74 double eor = 0.3; // threshold for end-of-run cut off
75 int specoff = 1; // temporary flag: if =1, only analysis performed
76 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~DEFINITIONS
77 std::array<int, 4> tdiff = {0, 0, 0, 0}; // run time extracted from ODB
78 int* td = tdiff.data();
79 int* p = addr.data();
80 double* q = freq.data();
81 int counter = 0;
82 int nds = 0;
83 int nscaler = 0;
84 int len = 70;
85 auto* line = new char[len];
86 auto* odb = new char[len];
87 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*
88 if(!(filelist.is_open())) {
89 std::cerr << "Failed to open filelist. Check path. " << std::endl;
90 exit(EXIT_FAILURE);
91 } else if(!(ODBlist.is_open())) {
92 std::cerr << "Failed to open ODBlist. Check path. " << std::endl;
93 exit(EXIT_FAILURE);
94 }
95 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*
96 // in the ODB file there are the lines:
97 // Start time binary = DWORD : 1445453916
98 // Stop time binary = DWORD : 1445454042
99 // the difference (stop-start) equals the run time in seconds. Useful for histogram binning purposes.
100 size_t posa = 0;
101 size_t posb = 0;
102 int sub = 28;
103 const char* starttime = "Start time binary";
104 const char* stoptime = "Stop time binary";
105 int tstart = 0;
106 int tend = 0;
107 int runtime = 0;
108 int nppg = 0;
109 std::string odbline;
110 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*
111 int line_count = std::count( // read in number of lines(files)
112 std::istreambuf_iterator<char>(filelist), std::istreambuf_iterator<char>(), '\n');
113 printf(DMAGENTA "-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/" RESET_COLOR
114 "\n");
115 printf(DGREEN "Started Deadtime v1 / 04-Mar-2016 / mbowry@triumf.ca" RESET_COLOR "\n");
116 printf(DBLUE "Calculates deadtime on a channel-by-channel basis using precision scaler data" RESET_COLOR "\n");
117 printf(DBLUE "Sub-runs must be merged (e.g. using gadd) before running this program" RESET_COLOR "\n");
118 printf(DBLUE "A list of fragments (run files) and their corresponding ODB files must be provided" RESET_COLOR "\n");
119 printf(DGREEN "Number of files loaded: %i" RESET_COLOR "\n", line_count);
120 printf(DBLUE "Spectra are saved in %s" RESET_COLOR "\n", fname);
121 printf(DBLUE "Deadtime results are recorded in %s" RESET_COLOR "\n", kname);
122 printf(DMAGENTA "-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/" RESET_COLOR
123 "\n");
124 printf(DBLUE "Working, be patient .. " RESET_COLOR "\n");
125 printf("\n");
126 filelist.clear();
127 filelist.seekg(0, std::ios::beg);
128 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*
129 Printaddress(p); // Credit to to E.Kwan for this part
130 p = addr.data();
131 while(counter < line_count) {
132 filelist.getline(line, len, '\n');
133 const char* filename = line;
134 //-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/ODB STUFF
135 ODBlist.getline(odb, len, '\n'); //<retreive ODB file name ('odb' assigned here)
136 std::ifstream InFile(odb); //<set the I/O stream to the current file
137 while(InFile.good()) {
138 getline(InFile, odbline);
139 posa = odbline.find(starttime);
140 posb = odbline.find(stoptime);
141 if(posa != std::string::npos && odbline.size() > static_cast<size_t>(sub)) {
142 odbline = odbline.substr(sub);
143 tstart = std::stoi(odbline, nullptr, 10);
144 } else if(posb != std::string::npos && odbline.size() > static_cast<size_t>(sub - 1)) {
145 odbline = odbline.substr((sub - 1));
146 tend = std::stoi(odbline, nullptr, 10);
147 }
148 }
149 runtime = tend - tstart;
150 nppg = floor(runtime / patlen);
151 tdiff[counter] = runtime;
152 std::cout << "\n";
153 if(runtime <= 600) {
154 printf(DBLUE "Read ODB %s : run time = %i seconds / %i transitions" RESET_COLOR "\n", odb, runtime, nppg);
155 } else if(runtime > 600) {
156 printf(DBLUE "Read ODB %s : run time = %i minutes / %i transitions" RESET_COLOR "\n", odb, (runtime / 60),
157 nppg);
158 }
159 //-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-/-
160 for(int index = scaleri; index < (scalerf + 1); index++) {
161 if(specoff == 1) {
162 printf(DMAGENTA "Creation of histograms suppressed .. performing analysis step only." RESET_COLOR "\n");
163 } else {
164 MakeSpectra(filename, counter, fname, nsclr, ncycle, q, p, index, td, thresh);
165 }
166 p = addr.data();
167 nds += 1;
168 }
169 counter++;
170 q++; // these two lines had a de-reference which makes no sense (and makes the compiler complain)
171 td++; // unless this was meant to increment the value the pointers are pointing to (in which case it should read
172 // (*q)++); VB
173 tstart = tend = 0;
174 }
175 delete[] line;
176 delete[] odb;
177 q = freq.data();
178 td = tdiff.data();
179 if(counter > 0) {
180 nscaler = nds / counter;
181 }
182 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*
183 CheckFile(fname);
184 DoAnalysis(fname, line_count, q, nsclr, patlen, ncycle, td, eor, hname, iname, jname, kname, lname, mname, nname,
185 nscaler);
186 printf(DBLUE "Spectra are saved in %s" RESET_COLOR "\n", fname);
187 printf(DBLUE "Deadtime results are %s" RESET_COLOR "\n", kname);
188 printf(DBLUE "Done. Control-c to exit." RESET_COLOR "\n");
189 theApp.Run(kTRUE);
190
191 return EXIT_SUCCESS;
192}
193
194void Printaddress(int* channel)
195{
196 printf(DBLUE "Addresses with both source and pulser inputs:" RESET_COLOR "\n");
197 for(int i = 0; i < 10; ++i) {
198 printf(DBLUE "%04x " RESET_COLOR, *channel);
199 channel++;
200 }
201 printf("\n");
202}
203
204void MakeSpectra(const char*& filename, int& prog, const char*& fname, int& nsclr, int& ncycle, double*, int* channel,
205 int& index, const int* const trun, double&)
206{
207 int nsc = nsclr;
208
209 // define spectra
210 auto** grif = new TH1D*[nsc];
211 // define file pointer
212 TFile* vs = nullptr;
213
214 // make spectra
215 auto* rf = new TFile(filename, "read");
216 auto* maple = static_cast<TTree*>(rf->Get("ScalerTree")); // Scaler data
217
218 int nofBins = *trun / ncycle;
219 double xaxis = 0;
220 double yaxis = 0;
221 double prev = 0;
222 double xpast = 0;
223 int j = 0;
224 int k = 0;
225 double clk = 20e-9; // 2 clock tics (20ns)
226
227 while(j < nsc) {
228 grif[j] =
229 new TH1D(Form("grif%d_0x%04x_%d", prog, *channel, index),
230 Form("Address 0x%04x, scaler %i vs time in cycle; time [s]; counts/%d s", *channel, index, ncycle),
231 nofBins, 0., *trun);
232 j++;
233 channel++; // used to have a de-reference, see comments above; VB
234 }
235
236 TScalerData* scaler = nullptr;
237 TScaler(maple).Clear();
238 maple->SetBranchAddress("TScalerData", &scaler);
239 Long64_t nentries = maple->GetEntries();
240
241 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*
242 // Build histograms directly from trees
243 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*
244 if(prog == 0 && index == 0) {
245 vs = new TFile(fname, "recreate");
246 } else {
247 vs = new TFile(fname, "update");
248 }
249 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*
250 for(int i = 0; i < nsc; i++) {
251 channel--; // used to have a de-reference, see comments above; VB
252 }
253 while(k < nsc) {
254 for(Long64_t e = 0; e < nentries; e++) {
255 maple->GetEntry(e);
256 if(scaler->GetAddress() == static_cast<UInt_t>(*channel)) {
257 xaxis = static_cast<double>(scaler->GetTimeStamp()) / 1e9;
258 // we check both the value of the scaler and the timestamp (ts difference should be = readout time)
259 if(prev != 0 && prev < scaler->GetScaler(index) && (xaxis - xpast) <= static_cast<double>(ncycle) + clk) {
260 yaxis = (scaler->GetScaler(index) - prev);
261 grif[k]->Fill(xaxis, yaxis);
262 }
263 prev = scaler->GetScaler(index);
264 xpast = xaxis;
265 }
266 }
267 k++;
268 channel++;
269 prev = xpast = 0;
270 }
271
272 printf(DGREEN "Created histograms for scaler %i : file = %s" RESET_COLOR "\n", index, filename);
273 // write hists
274 for(int i = 0; i < nsc; i++) {
275 grif[i]->Write();
276 }
277 vs->Close();
278}
279
280void CheckFile(const char*& fname)
281{
282 TFile f(fname);
283 if(f.IsZombie()) {
284 printf("\n");
285 printf(DRED "Error opening ROOT file %s" RESET_COLOR "\n", fname);
286 exit(-1);
287 } else {
288 printf("\n");
289 printf(DBLUE "ROOT file %s opens with no problems.." RESET_COLOR "\n", fname);
290 }
291}
292
293void DoAnalysis(const char*& fname, int& nfile, double* rate, int& nsclr, int& patlen, int&, int* trun, double& eor,
294 const char*& hname, const char*& iname, const char*& jname, const char*& kname, const char*& lname,
295 const char*& mname, const char*& nname, int& nscaler)
296{
297
298 auto* vs = new TFile(fname, "read");
299 std::ofstream ofile;
300 ofile.open("diagnostic.txt");
301 FILE* random = fopen(hname, "w");
302 FILE* combine = fopen(iname, "w");
303 FILE* rng = fopen(jname, "w");
304 FILE* deadtime = fopen(kname, "w");
305 FILE* asymerror = fopen(lname, "w");
306 FILE* randdt = fopen(mname, "w");
307 FILE* randw = fopen(nname, "w");
308 ofile.precision(4);
309
310 int nsc = nsclr;
311 int cnt = 0;
312 std::array<int, 2> ppgstat = {0, 0};
313 // generate random seed from system time for use in error analysis
314 time_t timer = time(nullptr);
315 srand(timer);
316
317 TFile f(fname);
318 TIter next(f.GetListOfKeys());
319 TKey* key = nullptr;
320 auto** spec = new TH1D*[nfile * (nsc * nscaler)];
321
322 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~USER EDITABLE
323 // limits for deadtime matrices [us]
324 // NOTE: The order of the limits (rp1,rp2 etc.) MUST be identical to the file order
325 std::array<double, 8> rp1 = {-15, 15, -15, 15, 0, 20, 85, 115}; // 2kHz, scaler0-3
326 std::array<double, 8> rp2 = {-15, 15, -15, 15, 0, 20, 85, 115}; // 5kHz, scaler0-3
327 std::array<double, 8> rp3 = {-15, 15, -15, 15, 0, 20, 170, 230}; // 10kHz, scaler0-3
328 std::array<double, 8> rp4 = {-15, 15, -15, 15, 0, 20, 200, 300}; // 20kHz, scaler0-3
329 double* lowrtau = rp1.data(); // internal pointers: points to correction coefficients in each array
330 double* upprtau = &rp1[1]; //
331 int cflag = 0; // counter
332 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~
333
334 std::cout << "\n";
335 // output headers
336 fprintf(deadtime, "#np=pulser,nr=source,nrp=source+pulser,tau=deadtime,erb=FWHM est. err,erf=full range "
337 "err,erp=standard err propagation");
338 fprintf(deadtime, "\n");
339 fprintf(deadtime, "%s\t\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s", "#spec", "chan", "np", "nr", "+/-",
340 "nrp", "+", "-", "tau", "+erb", "-erb", "+erf", "-erf", "+/-erp");
341 fprintf(deadtime, "\n");
342 fprintf(deadtime, "%s\t\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s", "#", "-", "kHz", "kHz", "kHz", "kHz",
343 "kHz", "kHz", "us", "us", "us", "us", "us", "us");
344 fprintf(deadtime, "\n");
345 fprintf(randdt, "%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s", "#rc", "rr", "rc-rr", "rtau", "lim1", "lim2", "flag",
346 "err(rc-rr)");
347 fprintf(randdt, "\n");
348
349 while((key = dynamic_cast<TKey*>(next())) != nullptr) {
350 int numpat = floor(*trun / patlen); // minimum number of expected transitions based on run time
351 const char* sname = key->GetName();
352 spec[cnt] = dynamic_cast<TH1D*>(vs->Get(sname));
353
354 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~*~*~*variables
355 double wcomb = 0;
356 double ncomb = 0;
357 double rcomb = 0;
358 double sdcomb = 0;
359 double wrand = 0;
360 double nrand = 0;
361 double rrand = 0;
362 double sdrand = 0;
363 double sum = 0;
364 double sumc = 0;
365 double maxl = 0;
366 double minl = 1e6;
367 double x = 0;
368 double w = 0;
369 double diff = 0;
370 double rcmin = 0;
371 int minb = 0;
372 double lbd = 0.;
373 double ubd = 0.;
374 double rmax = 0.;
375 int flag = 0;
376 double z = 0;
377 double y = 0;
378 double v = 0;
379 double dz = 0;
380 double dv = 0;
381 double d2z = 0;
382 int fbin = 10;
383 double max = 0;
384 double dlim = 0;
385 int nt = 0;
386 int ord = 0;
387 int sref = 0;
388 int pref = 0;
389 int chop = 2; //'chop'= ignore first/last two bins of each pattern
390 auto** ppg = new int*[numpat];
391 for(int i = 0; i < numpat; ++i) {
392 ppg[i] = new int[2];
393 };
394 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*find PPG transition boundaries
395 // run through and find boundaries (pulser on/off etc.)
396 //'ord' defines increase (1) or decrease (0) in rate
397 int xbins = *trun;
398 auto** trans = new double*[xbins];
399 for(int i = 0; i < xbins; ++i) {
400 trans[i] = new double[3];
401 };
402 auto** freq = new double*[fbin];
403 for(int i = 0; i < fbin; ++i) {
404 freq[i] = new double[2];
405 };
406 auto** bnd = new int*[numpat];
407 for(int i = 0; i < numpat; ++i) {
408 bnd[i] = new int[2];
409 };
410 // initialise the bnd matrix (prevents the program from hanging up later!)
411 for(int i = 0; i < numpat; i++) {
412 bnd[i][0] = 0;
413 bnd[i][1] = 0;
414 }
415 // find change in counts between bins
416 // we always look at the relative change transitioning 'up' (OFF->ON), even if we are transitioning 'down'
417 // we don't want to misidentify gaps in spectrum as transitions! (& vice versa)
418 for(int i = 0; i <= xbins; i++) {
419 x = spec[cnt]->GetBinContent(i);
420 if(w > 0 && x > 0 && x > w) {
421 diff = ((x - w) / w);
422 ord = 1;
423 } else if(w > 0 && x > 0 && x < w) {
424 diff = ((w - x) / x);
425 ord = 0;
426 } else if(x > 0 && w == 0) {
427 // this statement deals with gaps in spectrum
428 // we look at earlier bins and calculate difference compared to the current bin (up to t=0)
429 for(int j = (i - 1); j >= 0; j--) {
430 z = spec[cnt]->GetBinContent(j);
431 if(z > 0 && x > z) {
432 diff = ((x - z) / z);
433 ord = 1;
434 break;
435 }
436 if(z > 0 && x < z) {
437 diff = ((z - x) / x);
438 ord = 0;
439 break;
440 // what if there are no earlier bins >0? (rare, unless at very start of run)
441 }
442 }
443 } else if(w > 0 && x == 0) {
444 // we look at later bins and calculate difference compared to the current bin (up to t=trun)
445 // what if there are no later bins with counts >0? (rare, unless at very end of run)
446 for(int j = (i + 1); j <= *trun; j++) {
447 z = spec[cnt]->GetBinContent(j);
448 if(z > 0) {
449 diff = 0;
450 break;
451 }
452 }
453 if(z == 0) {
454 for(int j = *trun; j >= (*trun - static_cast<int>(0.7 * patlen)); j--) {
455 y = spec[cnt]->GetBinContent(j);
456 if(y > 0 && v > 0) {
457 dz = (abs(y - v) / v);
458 if(dz > 0 && dv > 0) {
459 d2z = abs(dz - dv);
460 }
461 }
462 if(d2z > eor) {
463 flag += 1;
464 diff = 0;
465 break;
466 }
467 v = spec[cnt]->GetBinContent(j);
468 dv = dz;
469 }
470 if(flag == 0) {
471 diff = (0.9 * rmax);
472 }
473 flag = 0;
474 }
475 } else if(x == w) {
476 diff = 0;
477 } else if(x == 0 && w == 0) {
478 diff = 0;
479 }
480 trans[i][0] = i;
481 trans[i][1] = diff;
482 trans[i][2] = ord;
483 w = spec[cnt]->GetBinContent(i);
484 if(abs(diff) > rmax) {
485 rmax = diff;
486 }
487 }
488 rmax = rmax + (0.1 * rmax);
489 // set up the frequency histogram
490 for(int i = 0; i <= xbins; i++) {
491 for(int j = 0; j < fbin; j++) {
492 if(i == 1) {
493 freq[j][0] = (0. + j * (rmax / fbin));
494 freq[j][1] = 0;
495 }
496 if(abs(trans[i][1]) > (0. + j * (rmax / fbin)) &&
497 abs(trans[i][1]) <= (rmax / fbin) + (j * (rmax / fbin))) {
498 freq[j][1] = freq[j][1] + 1;
499 }
500 }
501 }
502 // set transition limit using histogram
503 for(int j = 0; j < fbin; j++) {
504 if(freq[j][1] > max) {
505 max = freq[j][1];
506 dlim = (2 * ((rmax / fbin) + (j * (rmax / fbin))));
507 }
508 }
509 // std::cout<<"PPG transitions assumed above "<<(dlim*100.)<<" % change in rate.."<<std::endl;
510 for(int i = 0; i <= xbins; i++) {
511 if(abs(trans[i][1]) > dlim) {
512 nt += 1;
513 if(nt > numpat) {
514 printf(DRED "Warning: Found more PPG transition(s) than expected in spectrum %s (%i/%i)" RESET_COLOR
515 "\n",
516 sname, nt, numpat);
517 // std::cout<<"(see bin number "<<i<<" in spectrum "<<sname<<")"<<std::endl;
518 ppgstat[1] += 1;
519 break;
520 }
521 bnd[nt - 1][0] = i;
522 bnd[nt - 1][1] = static_cast<int>(trans[i][2]);
523 }
524 }
525 if(nt == numpat) {
526 printf(DGREEN "Found correct number of PPG transitions in spectrum %s (%i/%i)" RESET_COLOR "\n", sname, nt,
527 numpat);
528 ppgstat[0] += 1;
529 } else if(nt < numpat) {
530 printf(DRED "Warning: Found too few PPG transitions in spectrum %s (%i/%i)" RESET_COLOR "\n", sname, nt,
531 numpat);
532 ppgstat[1] += 1;
533 }
534 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
535 // DECIDE THE PPG ORDER (PULSER OFF/PULSER ON)
536 ppg[0][0] = 0;
537 for(int i = 0; i < numpat; i++) {
538 ppg[i + 1][0] = bnd[i][0];
539 ppg[i][1] = (bnd[i][0]) - 1;
540 }
541 if(bnd[0][1] == 0) { // pulser must be first (first transition is pulser OFF)
542 sref = 1;
543 pref = 0;
544 } else if(bnd[0][1] == 1) { // source must be first (first transition is pulser ON)
545 sref = 0;
546 pref = 1;
547 }
548 // diagnostic only
549 if(cnt % nsc == 0) {
550 for(int i = 0; i < numpat; i++) {
551 ofile << ppg[i][0] << "\t" << ppg[i][1] << std::endl;
552 }
553 ofile << "/" << std::endl;
554 }
555 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
556 // SOURCE ONLY
557 for(int i = sref; i < numpat; i += 2) {
558 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
559 x = spec[cnt]->GetBinContent(j);
560 if(x != 0) {
561 wrand = wrand + x;
562 nrand += 1;
563 }
564 }
565 }
566 rrand = (wrand / nrand); // mean
567 // diagnostic spectrum (dspec) parameters
568 int lim1 = static_cast<int>(rrand - (0.5 * dlim * rrand));
569 int lim2 = static_cast<int>(rrand + (0.5 * dlim * rrand));
570 int dsbin = (lim2 - lim1) / 20;
571 auto** dspec = new int*[dsbin];
572 for(int i = 0; i < dsbin; ++i) {
573 dspec[i] = new int[2];
574 };
575 for(int i = 0; i < dsbin; i++) {
576 dspec[i][0] = lim1 + (i * ((lim2 - lim1) / dsbin));
577 dspec[i][1] = 0;
578 }
579 // standard deviation / increment dspec
580 for(int i = sref; i < numpat; i += 2) {
581 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
582 x = spec[cnt]->GetBinContent(j);
583 if(cnt % nsc == 0) { // dspec for channel 0 only
584 for(int k = 0; k < dsbin; k++) {
585 if(x >= dspec[k][0] && x < dspec[k + 1][0]) {
586 dspec[k][1] += 1;
587 }
588 }
589 }
590 if(x != 0) {
591 sum = sum + pow((x - rrand), 2);
592 }
593 }
594 }
595 sdrand = sqrt(sum / (nrand - 1)); // standard deviation
596 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~FREQUENCY HISTOGRAM
597 if(cnt % nsc == 0) {
598 for(int i = 0; i < dsbin; i++) { // dspec
599 fprintf(random, "%i \t %i", dspec[i][0], dspec[i][1]);
600 fprintf(random, "\n");
601 }
602 fprintf(random, "/");
603 fprintf(random, "\n");
604 }
605 // std::cout<<"Mean = "<<rrand<<", standard deviation = "<<sdrand<<std::endl;
606 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
607 // SOURCE+PULSER
608 for(int i = pref; i < numpat; i += 2) {
609 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
610 x = spec[cnt]->GetBinContent(j);
611 if(x != 0) {
612 wcomb = wcomb + x;
613 ncomb += 1;
614 } else {
615 continue;
616 }
617 }
618 }
619 rcomb = (wcomb / ncomb); // mean
620 // now calculate max. likelihood
621 for(int i = pref; i < numpat; i += 2) {
622 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
623 x = spec[cnt]->GetBinContent(j);
624 if(x != 0) {
625 maxl = (0.5 * (pow((x - rcomb), 2))) / x;
626 if(maxl < minl) {
627 minl = maxl;
628 minb = j;
629 }
630 } else {
631 continue;
632 }
633 }
634 }
635 rcmin = spec[cnt]->GetBinContent(minb);
636 lbd = rcmin;
637 ubd = rcmin;
638 // find min/max values of parabola with respect to minimum
639 for(int i = pref; i < numpat; i += 2) {
640 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
641 x = spec[cnt]->GetBinContent(j);
642 if(x != 0) {
643 maxl = (0.5 * (pow((x - rcmin), 2))) / x;
644 // ofile<<j<<"\t"<<x<<"\t"<<maxl<<std::endl;
645 if(maxl <= (minl + 0.5) && x < rcmin && x < lbd) { // fixed
646 lbd = x;
647 }
648 if(maxl <= (minl + 0.5) && x > rcmin && x > ubd) { // fixed
649 ubd = x;
650 }
651 } else {
652 continue;
653 }
654 }
655 }
656 // diagnostic spectrum (dspec) parameters (re-define for source+pulser)
657 lim1 = static_cast<int>(lbd - (2.0 * abs(rcmin - lbd)));
658 lim2 = static_cast<int>(ubd + (2.0 * abs(rcmin - ubd)));
659 // std::cout<<lim1<<"\t"<<lim2<<"\t"<<rcmin<<std::endl;
660 dsbin = (lim2 - lim1) / 20; // dspec[dsbin][2]; What was this statement supposed to do? It has no effect; VB
661 for(int i = 0; i < dsbin; i++) {
662 dspec[i][0] = lim1 + (i * ((lim2 - lim1) / dsbin));
663 dspec[i][1] = 0;
664 }
665 // standard deviation / increment dspec
666 for(int i = pref; i < numpat; i += 2) {
667 for(int j = (ppg[i][0] + chop); j <= (ppg[i][1] - chop); j++) {
668 x = spec[cnt]->GetBinContent(j);
669 if(cnt % nsc == 0) { // dspec for channel 0 only
670 for(int k = 0; k < dsbin; k++) {
671 if(x >= dspec[k][0] && x < dspec[k + 1][0]) {
672 dspec[k][1] += 1;
673 }
674 }
675 }
676 if(x != 0) {
677 sumc = sumc + pow((x - rcomb), 2);
678 }
679 }
680 }
681 sdcomb = sqrt(sumc / (ncomb - 1)); // standard deviation (not strictly applicable to source+pulser data)
682 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~FREQUENCY HISTOGRAM
683 if(cnt % nsc == 0) {
684 for(int i = 0; i < dsbin; i++) { // dspec
685 fprintf(combine, "%i \t %i", dspec[i][0], dspec[i][1]);
686 fprintf(combine, "\n");
687 }
688 fprintf(combine, "/");
689 fprintf(combine, "\n");
690 }
691 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
692 // To find the error in (rcmin-rrand), we generate a likelihood curve using an iterative procedure
693 // Ref. "Asymmetric Statistical Errors", Roger Barlow, http://arxiv.org/pdf/physics/0406120.pdf
694 // WARNING: u must be given a reasonable range. This is assumed to be ~1.5*max(a1,a2)
695
696 double a1 = abs(rcmin - ubd);
697 double a2 = abs(rcmin - lbd);
698 double maxa = std::max(a1, a2);
699 double uhi = (1.5 * maxa);
700 double ulo = -uhi;
701 double du = (uhi - ulo) / 1e3;
702 double u = ulo;
703 int itr = 0;
704 int umin = 0;
705 double w1 = 0.;
706 double x1 = 0.;
707 double x2 = 0.;
708 int nrow = static_cast<int>((uhi - ulo) / du);
709 double w2 = pow((pow(sdrand, 2)), 2) / ((2. * pow(sdrand, 2))); // const.
710 minl = -1e6;
711 auto** array = new double*[nrow];
712 for(int i = 0; i < nrow; ++i) {
713 array[i] = new double[2];
714 }
715 double sigm = 0.;
716 double sigp = 0.;
717 double lnl = 0.;
718
719 while(itr < nrow) {
720 if(itr == 0) {
721 x1 = 0.;
722 x2 = 0.;
723 } else {
724 x1 = (u * w1) / (w1 + w2);
725 x2 = (u * w2) / (w1 + w2);
726 }
727 w1 = pow(((a1 * a2) + ((a1 - a2) * x1)), 2) / ((2. * (a1 * a2)) + ((a1 - a2) * x1));
728 lnl = -0.5 * ((pow(x1, 2) / ((a1 * a2) + ((a1 - a2) * x1))) + (pow(x2, 2) / (pow(sdrand, 2))));
729 if(lnl <= 0.) { // reject positive likelihoods (these tend occur at large +/- values of u)
730 if(itr > 0 && lnl > minl) { // mean value
731 minl = lnl;
732 umin = itr;
733 }
734 array[itr][0] = u;
735 array[itr][1] = lnl;
736 u += du;
737 itr++;
738 } else {
739 u += du;
740 }
741 }
742 // determine upper/lower limits
743 sigm = array[umin][0];
744 sigp = array[umin][0];
745 for(int i = 1; i < nrow; i++) {
746 if(array[i][1] >= (minl - 0.5) && array[i][0] < array[umin][0] && array[i][0] < sigm) {
747 sigm = array[i][0];
748 }
749 if(array[i][1] >= (minl - 0.5) && array[i][0] > array[umin][0] && array[i][0] > sigp) {
750 sigp = array[i][0];
751 }
752 }
753 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*PARABOLA HISTOGRAM
754 if(cnt % nsc == 0) {
755 fprintf(asymerror, "%s\t %s\t %s", "#chan", "Err[Hz]", "Ln(L)");
756 fprintf(asymerror, "\n");
757 for(int i = 1; i < nrow; i++) {
758 fprintf(asymerror, "%i \t %.4E \t %.4E", cnt % nsc, array[i][0], array[i][1]);
759 fprintf(asymerror, "\n");
760 }
761 fprintf(asymerror, "/");
762 fprintf(asymerror, "\n");
763 }
764 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
765 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
766 double tau = ((1.0 / rrand) * (1.0 - sqrt((rcmin - rrand) / (*rate)))) * 1.0e6; // deadtime (us)
767 double rtau = 0;
768 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
769 // FINAL ERROR ANALYSIS*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
770 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
771 // Idea: now that the error boundaries are known for rrand, rcmin and (rcmin-rrand), generate a random number
772 // which chooses a value within these boundaries - the resulting tau is plotted and the range gives the error in
773 // tau
774
775 int iter = 10000;
776 int bin = 10;
777 double tempa = 0.;
778 double tempb = 0.;
779 double var1 = 0.;
780 double var2 = 0.;
781 double var3 = 0.;
782 double l1 = 0.;
783 double l2 = 0.;
784 int wbin = 40;
785 int wrow = 0;
786 int wsize = static_cast<int>(pow(wbin, 2));
787 auto** randcheck = new double*[bin];
788 for(int i = 0; i < bin; ++i) {
789 randcheck[i] = new double[2];
790 }
791 auto** randtau = new double*[iter];
792 for(int i = 0; i < iter; ++i) {
793 randtau[i] = new double[2];
794 }
795 auto** wspec = new double*[wsize];
796 for(int i = 0; i < wsize; ++i) {
797 wspec[i] = new double[3];
798 }
799 int binsize = 3;
800
801 // Check the "randomness" of the random number generator **RCHECK**
802 for(int i = 0; i < bin; i++) {
803 randcheck[i][0] = 0 + (static_cast<double>(i) * (1 / static_cast<double>(bin)));
804 randcheck[i][1] = 0;
805 }
806 // initialise matrix to get final uncertainty
807 for(int i = 0; i < wbin; i++) {
808 for(int j = 0; j < wbin; j++) {
809 wspec[wrow][0] = ((rcmin - rrand) + sigm) + (static_cast<double>(i) * ((sigp - sigm) / static_cast<double>(wbin))); // x,
810 // (rcmin-rrand)
811 wspec[wrow][1] = (*lowrtau) + (static_cast<double>(j) * ((*upprtau - *lowrtau) / static_cast<double>(wbin))); // y, (rtau);
812 wspec[wrow][2] = 0; // z, frequency
813 wrow += 1;
814 }
815 }
816 // generate all possible guesses for rcmin and rrand to satisfy test value of (rcmin-rrand) in the limits l1,l2
817 int i = 0;
818 while(i < (nrow - binsize)) {
819 if(array[i][0] >= sigm && array[i][0] <= sigp) {
820 l1 = ((rcmin - rrand) + array[i][0]);
821 l2 = ((rcmin - rrand) + array[i + binsize][0]);
822
823 for(int j = 0; j < iter; j++) {
824 tempa = rand() / static_cast<double>(RAND_MAX);
825 var1 = ((tempa * (a1 + a2)) + (rcmin - a2));
826 tempb = rand() / static_cast<double>(RAND_MAX);
827 var2 = ((tempb * (2.e0 * sdrand)) + (rrand - sdrand));
828 var3 = (var1 - var2);
829 //~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~'good'
830 //events
831 if(var3 >= l1 && var3 <= l2) {
832 // calculate deadtime using var3(rcmin-rrand), var2(rrand)
833 rtau = ((1.0 / var2) * (1.0 - sqrt(var3 / (*rate)))) * 1.0e6;
834 // build wspec matrix
835 for(int k = 0; k < wsize; k++) {
836 if(var3 >= wspec[k][0] && var3 < wspec[k + wbin][0] && rtau >= wspec[k][1] &&
837 rtau < wspec[k + 1][1]) {
838 wspec[k][2] += 1; // this seems to work now
839 }
840 }
841 //~*~*~*~*~*~*~*~*~*~*~*~~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~
842 } else {
843 continue;
844 }
845 }
846 i += binsize;
847 } else {
848 i += 1;
849 continue;
850 }
851 }
852 // diagnostic
853 // printf(DMAGENTA "(comb-rand)= %f, sigp= %f, sigm= %f, rand= %f, sdrand= %f" RESET_COLOR
854 // "\n",(rcmin-rrand),sigp,sigm,rrand,sdrand);
855 //~*~*~*~*~*~*~*~*~*~*~*~~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~
856 // Determine final error bounds using wspec matrix
857 double gmax = 0.;
858 double hmax = 0.;
859 double erbp = 0.;
860 double erbm = 0.;
861 double erfp = 0.;
862 double erfm = 0.;
863 double frmin = 1e6;
864 double frmax = 0.;
865 double srmin = 1e6;
866 double srmax = 0.;
867
868 for(i = 0; i < wsize; i++) {
869 if(wspec[i][2] > gmax) {
870 gmax = wspec[i][2];
871 hmax = (gmax / 2);
872 }
873 }
874 // calculate max/min range (select above half the height of 'peak')
875 for(i = 0; i < wsize; i++) {
876 if(wspec[i][1] > srmax && wspec[i][2] > hmax) {
877 srmax = wspec[i][1];
878 }
879 if(wspec[i][1] < srmin && wspec[i][2] > hmax) {
880 srmin = wspec[i][1];
881 }
882 }
883 // calculate max/min range (only select above zero)
884 for(i = 0; i < wsize; i++) {
885 if(wspec[i][1] > frmax && wspec[i][2] > 0) {
886 frmax = wspec[i][1];
887 }
888 if(wspec[i][1] < frmin && wspec[i][2] > 0) {
889 frmin = wspec[i][1];
890 }
891 }
892 // final errors: b=best estimate, f=full range
893 erbm = abs(tau - srmin);
894 erbp = abs(tau - srmax);
895 erfm = abs(tau - frmin);
896 erfp = abs(tau - frmax);
897 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*STANDARD ERROR PROPAGATION
898 double p1 = sqrt(pow(sdcomb, 2) + pow(sdrand, 2));
899 double p1r = p1 / (rcmin - rrand);
900 double p2 = (rcmin - rrand) / (*rate);
901 // Here we assume that sqrt(rel.err.) ~ 0.5 * rel.err. This is approximately true where the +/- uncertainties are
902 // similar.
903 double p3 = (0.5 * p1r * p2);
904 double eprop = tau * (sqrt(pow((p3 / (1 - sqrt(p2))), 2) + pow((sdrand / rrand), 2)));
905 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*RANDOM TAU HISTOGRAM
906 // if(cnt%nsc==0){
907 fprintf(randdt, "%i \t %.2f \t %.2f \t %.2f \t %.2f", cnt, srmin, srmax, frmin, frmax);
908 // fprintf(randdt, "#//end channel 0");
909 fprintf(randdt, "\n");
910 //}
911 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~WIDTH HISTOGRAM
912 if(cnt % nsc == 0) {
913 for(i = 0; i < wsize; i++) {
914 fprintf(randw, "%.2f \t %.2f \t %.2f", wspec[i][0], wspec[i][1], wspec[i][2]);
915 fprintf(randw, "\n");
916 }
917 fprintf(randw, "//end channel 0");
918 fprintf(randw, "\n");
919 }
920 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~RCHECK HISTOGRAM
921 // print random number generator results to file for a single spectrum
922 /*if(cnt%nsc==0){
923 for (i=0; i<bin; i++){
924 fprintf(rng, "%f \t %f", randcheck[i][0], randcheck[i][1]);
925 fprintf(rng, "\n");
926 }
927 fprintf(rng, "//end channel 0");
928 fprintf(rng, "\n");
929 }*/
930 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~PRINT FINAL RESULTS TO FILE
931 fprintf(deadtime, "%s\t%i\t%.1f\t%.1f\t%.2f\t%.1f\t%.2f\t%.2f\t %.1f\t%.1f\t%.1f\t%.1f\t%.1f\t%.1f\n", sname,
932 (cnt % nsc), ((*rate) / 1e3), (rrand / 1e3), (sdrand / 1e3), (rcmin / 1e3), (abs(rcmin - ubd) / 1e3),
933 (abs(rcmin - lbd) / 1e3), tau, erbp, erbm, erfp, erfm, eprop);
934 fprintf(deadtime, "\n");
935 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
936 cnt++;
937 // change pulser rate and runtime accordingly for each file
938 if(cnt % (nsc * nscaler) == 0) {
939 rate++; // used to have a de-reference, see comments above; VB
940 trun++; // used to have a de-reference, see comments above; VB
941 }
942 // change limits accordingly for each scaler in each file
943 if(cnt % (nsc) == 0) {
944 for(i = 0; i < 2; i++) {
945 lowrtau++; // used to have a de-reference, see comments above; VB
946 upprtau++; // used to have a de-reference, see comments above; VB
947 }
948 }
949 if(cnt % (nsc) == 0 && cnt % (nsc * nscaler) == 0) {
950 cflag++;
951 if(cflag == 1) {
952 lowrtau = &rp2[0]; // NOLINT(readability-container-data-pointer)
953 upprtau = &rp2[1];
954 } else if(cflag == 2) {
955 lowrtau = &rp3[0]; // NOLINT(readability-container-data-pointer)
956 upprtau = &rp3[1];
957 } else if(cflag == 3) {
958 lowrtau = &rp4[0]; // NOLINT(readability-container-data-pointer)
959 upprtau = &rp4[1];
960 }
961 }
962 } // END SPECTRUM ANALYSIS LOOP
963 //*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~*~**~*~*~*~*~*~**~*~*~*~*~*~*
964 delete[] spec;
965 ofile.close();
966 fclose(random);
967 fclose(combine);
968 fclose(rng);
969 fclose(deadtime);
970 fclose(asymerror);
971 fclose(randdt);
972 fclose(randw);
973 int good = ppgstat[0];
974 printf("\n");
975 // info/warnings
976 if(good < cnt) {
977 printf(DRED "Correct # of PPG transitions obtained for %i out of %i spectra" RESET_COLOR "\n", good, cnt);
978 } else {
979 printf(DBLUE "Correct # of PPG transitions obtained for %i out of %i spectra" RESET_COLOR "\n", good, cnt);
980 }
981}
int main(int argc, char *argv[])
Definition Deadtime.cxx:50
void CheckFile(const char *&fname)
Definition Deadtime.cxx:280
void DoAnalysis(const char *&fname, int &nfile, double *rate, int &nsclr, int &patlen, int &ncycle, int *trun, double &eor, const char *&hname, const char *&iname, const char *&jname, const char *&kname, const char *&lname, const char *&mname, const char *&nname, int &nscaler)
Definition Deadtime.cxx:293
void MakeSpectra(const char *&filename, int &prog, const char *&fname, int &nsclr, int &ncycle, double *rate, int *channel, int &index, const int *trun, double &thresh)
Definition Deadtime.cxx:204
void Printaddress(int *channel)
Definition Deadtime.cxx:194
#define DRED
Definition Globals.h:18
#define DGREEN
Definition Globals.h:17
#define DBLUE
Definition Globals.h:15
#define RESET_COLOR
Definition Globals.h:5
#define DMAGENTA
Definition Globals.h:20
ULong64_t GetTimeStamp() const
Definition TScaler.h:79
std::vector< UInt_t > GetScaler() const
Definition TScaler.h:70
UInt_t GetAddress() const
Definition TScaler.h:66
void Clear(Option_t *opt="") override
Definition TScaler.cxx:228