source: trunk/Mars/mcore/DrsCalib.h@ 14896

Last change on this file since 14896 was 14865, checked in by tbretz, 12 years ago
Added a new function to treat the spikes; improved treatment for steps when they coincide with spikes and the rms is large
File size: 33.9 KB
Line 
1#ifndef MARS_DrsCalib
2#define MARS_DrsCalib
3
4#include <math.h> // fabs
5#include <errno.h> // errno
6
7#include "fits.h"
8#include "ofits.h"
9
10#ifdef __MARS__
11#include "MTime.h"
12#endif
13
14class DrsCalibrate
15{
16protected:
17 uint64_t fNumEntries;
18
19 size_t fNumSamples;
20 size_t fNumChannels;
21
22 std::vector<int64_t> fSum;
23 std::vector<int64_t> fSum2;
24
25public:
26 DrsCalibrate() : fNumEntries(0), fNumSamples(0), fNumChannels(0) { }
27 void Reset()
28 {
29 fNumEntries = 0;
30 fNumSamples = 0;
31 fNumChannels = 0;
32
33 fSum.clear();
34 fSum2.clear();
35 }
36
37 void InitSize(uint16_t channels, uint16_t samples)
38 {
39 fNumChannels = channels;
40 fNumSamples = samples;
41
42 fSum.resize(samples*channels);
43 fSum2.resize(samples*channels);
44 }
45
46 void AddRel(const int16_t *val, const int16_t *start)
47 {
48 for (size_t ch=0; ch<fNumChannels; ch++)
49 {
50 const int16_t spos = start[ch];
51 if (spos<0)
52 continue;
53
54 const size_t pos = ch*1024;
55
56 for (size_t i=0; i<1024; i++)
57 {
58 // Value is relative to trigger
59 // Abs is corresponding index relative to DRS pipeline
60 const size_t rel = pos + i;
61 const size_t abs = pos + (spos+i)%1024;
62
63 const int64_t v = val[rel];
64
65 fSum[abs] += v;
66 fSum2[abs] += v*v;
67 }
68 }
69
70 fNumEntries++;
71 }
72
73 void AddRel(const int16_t *val, const int16_t *start,
74 const int32_t *offset, const uint32_t scale)
75 {
76 for (size_t ch=0; ch<fNumChannels; ch++)
77 {
78 const int16_t spos = start[ch];
79 if (spos<0)
80 continue;
81
82 const size_t pos = ch*fNumSamples;
83
84 for (size_t i=0; i<fNumSamples; i++)
85 {
86 // Value is relative to trigger
87 // Offset is relative to DRS pipeline
88 // Abs is corresponding index relative to DRS pipeline
89 const size_t rel = pos + i;
90 const size_t abs = pos + (spos+i)%fNumSamples;
91
92 const int64_t v = int64_t(val[rel])*scale-offset[abs];
93
94 fSum[abs] += v;
95 fSum2[abs] += v*v;
96 }
97 }
98
99 fNumEntries++;
100 }
101 /*
102 void AddAbs(const int16_t *val, const int16_t *start,
103 const int32_t *offset, const uint32_t scale)
104 {
105 for (size_t ch=0; ch<fNumChannels; ch++)
106 {
107 const int16_t spos = start[ch];
108 if (spos<0)
109 continue;
110
111 const size_t pos = ch*fNumSamples;
112
113 for (size_t i=0; i<fNumSamples; i++)
114 {
115 // Value is relative to trigger
116 // Offset is relative to DRS pipeline
117 // Abs is corresponding index relative to DRS pipeline
118 const size_t rel = pos + i;
119 const size_t abs = pos + (spos+i)%fNumSamples;
120
121 const int64_t v = int64_t(val[rel])*scale-offset[abs];
122
123 fSum[rel] += v;
124 fSum2[rel] += v*v;
125 }
126 }
127
128 fNumEntries++;
129 }*/
130
131 void AddAbs(const int16_t *val, const int16_t *start,
132 const int32_t *offset, const uint32_t scale)
133 {
134 // 1440 without tm, 1600 with tm
135 for (size_t ch=0; ch<fNumChannels; ch++)
136 {
137 const int16_t spos = start[ch];
138 if (spos<0)
139 continue;
140
141 const size_t pos = ch*fNumSamples;
142 const size_t drs = ch>1439 ? ((ch-1440)*9+8)*1024 : ch*1024;
143
144 for (size_t i=0; i<fNumSamples; i++)
145 {
146 // Value is relative to trigger
147 // Offset is relative to DRS pipeline
148 // Abs is corresponding index relative to DRS pipeline
149 const size_t rel = pos + i;
150 const size_t abs = drs + (spos+i)%1024;
151
152 const int64_t v = int64_t(val[rel])*scale-offset[abs];
153
154 fSum[rel] += v;
155 fSum2[rel] += v*v;
156 }
157 }
158
159 fNumEntries++;
160 }
161
162
163 static void ApplyCh(float *vec, const int16_t *val, int16_t start, uint32_t roi,
164 const int32_t *offset, const uint32_t scaleabs,
165 const int64_t *gain, const uint64_t scalegain)
166 {
167 if (start<0)
168 {
169 memset(vec, 0, roi);
170 return;
171 }
172
173 for (size_t i=0; i<roi; i++)
174 {
175 // Value is relative to trigger
176 // Offset is relative to DRS pipeline
177 // Abs is corresponding index relative to DRS pipeline
178 const size_t abs = (start+i)%1024;
179
180 const int64_t v =
181 + int64_t(val[i])*scaleabs-offset[abs]
182 ;
183
184 const int64_t div = gain[abs];
185 vec[i] = div==0 ? 0 : double(v)*scalegain/div;
186 }
187 }
188
189 static void ApplyCh(float *vec, const int16_t *val, int16_t start, uint32_t roi,
190 const int32_t *offset, const uint32_t scaleabs,
191 const int64_t *gain, const uint64_t scalegain,
192 const int64_t *trgoff, const uint64_t scalerel)
193 {
194 if (start<0)
195 {
196 memset(vec, 0, roi);
197 return;
198 }
199
200 for (size_t i=0; i<roi; i++)
201 {
202 // Value is relative to trigger
203 // Offset is relative to DRS pipeline
204 // Abs is corresponding index relative to DRS pipeline
205 const size_t abs = (start+i)%1024;
206
207 const int64_t v =
208 + (int64_t(val[i])*scaleabs-offset[abs])*scalerel
209 - trgoff[i]
210 ;
211
212 const int64_t div = gain[abs]*scalerel;
213 vec[i] = div==0 ? 0 : double(v)*scalegain/div;
214 }
215 }
216
217 static double FindStep(const size_t ch0, const float *vec, int16_t roi, const int16_t pos, const uint16_t *map)
218 {
219 // We have about 1% of all cases which are not ahndled here,
220 // because the baseline jumps up just before the readout window
221 // and down just after it. In this cases we could determine the jump
222 // from the board time difference or throw the event away.
223 if (pos==0 || pos>=roi)
224 return 0;
225
226 double step = 0; // before
227 double rms = 0; // before
228 int cnt = 0;
229
230 // Exclude TM channel
231 for (int p=0; p<8; p++)
232 {
233 const size_t hw = ch0+p;
234 const size_t sw = map[hw]*roi + pos;
235
236 const double diff = vec[sw]-vec[sw-1];
237
238 step += diff;
239 rms += (vec[sw]-vec[sw-1])*(vec[sw]-vec[sw-1]);
240
241 cnt++;
242 }
243
244 return cnt==0 ? 0 : step/cnt;
245 }
246
247 static void SubtractStep(const size_t ch0, const double avg, float *vec, int16_t roi, int32_t pos, const uint16_t *map)
248 {
249 if (pos==0 || pos>=roi)
250 return;
251
252 const int begin = avg>0 ? pos : 0;
253 const int end = avg>0 ? roi : pos;
254
255 const double sub = fabs(avg);
256
257 for (int p=0; p<9; p++)
258 {
259 for (int j=begin; j<end; j++)
260 {
261 const size_t hw = ch0+p;
262 const size_t sw = map[hw]*roi + j;
263
264 vec[sw] -= sub;
265 }
266 }
267 }
268
269 struct Step
270 {
271 Step() : avg(0), rms(0), pos(0), cnt(0) { }
272 double avg;
273 double rms;
274 double pos;
275 uint16_t cnt;
276
277 static bool sort(const Step &s, const Step &r) { return s.avg<r.avg; }
278 };
279
280 static Step AverageSteps(const std::vector<Step>::iterator beg, const std::vector<Step>::iterator end)
281 {
282 Step rc;
283 for (auto it=beg; it!=end; it++)
284 {
285 rc.pos += it->pos;
286 rc.avg += it->avg;
287 rc.rms += it->avg*it->avg;
288 }
289
290 rc.cnt = end-beg;
291
292 rc.pos /= rc.cnt;
293 rc.avg /= rc.cnt;
294 rc.rms /= rc.cnt;
295
296 rc.rms = sqrt(rc.rms-rc.avg*rc.avg);
297
298 return rc;
299 }
300
301
302 static Step CorrectStep(float *vec, uint16_t nch, uint16_t roi,
303 const int16_t *prev, const int16_t *start,
304 const int16_t offset, const uint16_t *map)
305 {
306
307 std::vector<Step> list;
308
309 // Fill steps into array
310 // Exclude broken pixels?
311 // Remove maximum and minimum patches (4max and 4min)?
312 for (size_t ch=0; ch<nch; ch += 9)
313 {
314 if (prev[ch]<0 || start[ch]<0)
315 continue;
316
317 const int16_t dist = (prev[ch]-start[ch]+1024+offset)%1024;
318 const double step = FindStep(ch, vec, roi, dist, map);
319 if (step==0)
320 continue;
321
322 Step rc;
323 rc.pos = dist;
324 rc.avg = step;
325 list.push_back(rc);
326 }
327
328 if (list.size()==0)
329 return Step();
330
331 Step rc = AverageSteps(list.begin(), list.begin()+list.size());;
332
333 if (rc.avg==0)
334 return Step();
335
336 // std::cout << " A0=" << rc.avg << " rms=" << rc.rms << std::endl;
337 if (rc.rms>5)
338 {
339 sort(list.begin(), list.end(), Step::sort);
340
341 //for (auto it=list.begin(); it!=list.end(); it++)
342 // std::cout << " " << it->avg << std::endl;
343
344 const Int_t skip = list.size()/10;
345 rc = AverageSteps(list.begin()+skip, list.begin()+list.size()-skip);
346
347 // std::cout << " A1=" << rc.avg << " rms=" << rc.rms << std::endl;
348 }
349
350 for (size_t ch=0; ch<nch; ch += 9)
351 {
352 const int16_t dist = (prev[ch]-start[ch]+1024+offset)%1024;
353 SubtractStep(ch, rc.avg, vec, roi, dist, map);
354 }
355
356 return rc;
357 }
358
359 static void RemoveSpikes(float *vec, uint32_t roi)
360 {
361 if (roi<4)
362 return;
363
364 for (size_t ch=0; ch<1440; ch++)
365 {
366 float *p = vec + ch*roi;
367
368 for (size_t i=1; i<roi-2; i++)
369 {
370 if (p[i]-p[i-1]>25 && p[i]-p[i+1]>25)
371 {
372 p[i] = (p[i-1]+p[i+1])/2;
373 }
374
375 if (p[i]-p[i-1]>22 && fabs(p[i]-p[i+1])<4 && p[i+1]-p[i+2]>22)
376 {
377 p[i] = (p[i-1]+p[i+2])/2;
378 p[i+1] = p[i];
379 }
380 }
381 }
382 }
383
384 static void RemoveSpikes2(float *vec, uint32_t roi)//from Werner
385 {
386 if (roi<4)
387 return;
388
389 for (size_t ch=0; ch<1440; ch++)
390 {
391 float *p = vec + ch*roi;
392
393 std::vector<float> Ameas(p, p+roi);
394
395 std::vector<float> diff(roi);
396 for (size_t i=1; i<roi-1; i++)
397 diff[i] = (p[i-1] + p[i+1])/2 - p[i];
398
399 //std::vector<float> N1mean(roi);
400 //for (size_t i=1; i<roi-1; i++)
401 // N1mean[i] = (p[i-1] + p[i+1])/2;
402
403 const float fract = 0.8;
404
405 for (size_t i=0; i<roi-3; i++)
406 {
407 if (diff[i]<5)
408 continue;
409
410 if (Ameas[i+2] - (Ameas[i] + Ameas[i+3])/2 > 10)
411 {
412 p[i+1]= (Ameas[i+3] - Ameas[i])/3 + Ameas[i];
413 p[i+2]= 2*(Ameas[i+3] - Ameas[i])/3 + Ameas[i];
414
415 i += 3;
416
417 continue;
418 }
419
420 if ( (diff[i+1]<-diff[i]*fract*2) && (diff[i+2]>10) )
421 {
422 p[i+1] = (Ameas[i]+Ameas[i+2])/2;
423 diff[i+2] = (p[i+1] + Ameas[i+3])/2 - Ameas[i+2];
424
425 i += 2;
426 }
427
428 // const float x = Ameas[i] - N1mean[i];
429 // if (x > -5.)
430 // continue;
431
432 // if (Ameas[i+2] - (Ameas[i] + Ameas[i+3])/2. > 10.)
433 // {
434 // p[i+1]= (Ameas[i+3] - Ameas[i])/3 + Ameas[i];
435 // p[i+2]= 2*(Ameas[i+3] - Ameas[i])/3 + Ameas[i];
436 // i += 3;
437 // continue;
438 // }
439
440 // const float xp = Ameas[i+1] - N1mean[i+1];
441 // const float xpp = Ameas[i+2] - N1mean[i+2];
442
443 // if ( (xp > -2.*x*fract) && (xpp < -10.) )
444 // {
445 // p[i+1] = N1mean[i+1];
446 // N1mean[i+2] = Ameas[i+1] - Ameas[i+3]/2;
447 //
448 // i += 2;
449 // }
450 }
451 }
452 }
453
454 static void RemoveSpikes3(float *vec, uint32_t roi)//from Werner
455 {
456 const float SingleCandidateTHR = -10.;
457 const float DoubleCandidateTHR = -5.;
458
459 const std::vector<float> src(vec, vec+roi);
460
461 std::vector<float> diff(roi);
462 for (size_t i=1; i<roi-1; i++)
463 diff[i] = src[i] - (src[i-1] + src[i+1])/2;
464
465 // find the spike and replace it by mean value of neighbors
466 for (unsigned int i=1; i<roi-3; i++)
467 {
468 // Speed up (no leading edge)
469 if (diff[i]>=DoubleCandidateTHR)
470 continue;
471
472 //bool checkDouble = false;
473
474 // a single spike candidate
475 if (diff[i]<SingleCandidateTHR)
476 {
477 // check consistency with a single channel spike
478 if (diff[i+1] > -1.6*diff[i])
479 {
480 vec[i+1] = (src[i] + src[i+2]) / 2;
481
482 i += 2;
483
484 /*** NEW ***/
485 continue;
486 /*** NEW ***/
487 }
488 /*
489 else
490 {
491 // do nothing - not really a single spike,
492 // but check if it is a double
493 checkDouble = true;
494 }*/
495 }
496
497 // a double spike candidate
498 //if (diff[i]>DoubleCandidateTHR || checkDouble == 1)
499 {
500 // check the consistency with a double spike
501 if ((diff[i+1] > -DoubleCandidateTHR) &&
502 (diff[i+2] > -DoubleCandidateTHR))
503 {
504 vec[i+1] = (src[i+3] - src[i])/3 + src[i];
505 vec[i+2] = 2*(src[i+3] - src[i])/3 + src[i];
506
507 //vec[i] = (src[i-1] + src[i+2]) / 2.;
508 //vec[i+1] = (src[i-1] + src[i+2]) / 2.;
509
510 //do not care about the next sample it was the spike
511 i += 3;
512 }
513 }
514 }
515 }
516
517 static void SlidingAverage(float *const vec, const uint32_t roi, const uint16_t w)
518 {
519 if (w==0)
520 return;
521
522 if (w>roi)
523 return;
524
525 for (float *pix=vec; pix<vec+1440*roi; pix += roi)
526 {
527 for (float *ptr=pix; ptr<pix+roi-w; ptr++)
528 {
529 for (float *p=ptr+1; p<ptr+w; p++)
530 *ptr += *p;
531 *ptr /= w;
532 }
533 }
534 }
535
536 std::pair<std::vector<double>,std::vector<double> > GetSampleStats() const
537 {
538 if (fNumEntries==0)
539 return make_pair(std::vector<double>(),std::vector<double>());
540
541 std::vector<double> mean(fSum.size());
542 std::vector<double> error(fSum.size());
543
544 std::vector<int64_t>::const_iterator it = fSum.begin();
545 std::vector<int64_t>::const_iterator i2 = fSum2.begin();
546 std::vector<double>::iterator im = mean.begin();
547 std::vector<double>::iterator ie = error.begin();
548
549 while (it!=fSum.end())
550 {
551 *im = /*cnt<fResult.size() ? fResult[cnt] :*/ double(*it)/fNumEntries;
552 *ie = sqrt(double(*i2*int64_t(fNumEntries) - *it * *it))/fNumEntries;
553
554 im++;
555 ie++;
556 it++;
557 i2++;
558 }
559
560
561 /*
562 valarray<double> ...
563
564 mean /= fNumEntries;
565 error = sqrt(error/fNumEntries - mean*mean);
566 */
567
568 return make_pair(mean, error);
569 }
570
571 void GetSampleStats(float *ptr, float scale) const
572 {
573 const size_t sz = fNumSamples*fNumChannels;
574
575 if (fNumEntries==0)
576 {
577 memset(ptr, 0, sizeof(float)*sz*2);
578 return;
579 }
580
581 std::vector<int64_t>::const_iterator it = fSum.begin();
582 std::vector<int64_t>::const_iterator i2 = fSum2.begin();
583
584 while (it!=fSum.end())
585 {
586 *ptr = scale*double(*it)/fNumEntries;
587 *(ptr+sz) = scale*sqrt(double(*i2*int64_t(fNumEntries) - *it * *it))/fNumEntries;
588
589 ptr++;
590 it++;
591 i2++;
592 }
593 }
594
595 static double GetPixelStats(float *ptr, const float *data, uint16_t roi)
596 {
597 if (roi==0)
598 return -1;
599
600 const int beg = roi>10 ? 10 : 0;
601
602 double max = 0;
603 for (int i=0; i<1440; i++)
604 {
605 const float *vec = data+i*roi;
606
607 int pos = beg;
608 double sum = vec[beg];
609 double sum2 = vec[beg]*vec[beg];
610 for (int j=beg+1; j<roi; j++)
611 {
612 sum += vec[j];
613 sum2 += vec[j]*vec[j];
614
615 if (vec[j]>vec[pos])
616 pos = j;
617 }
618 sum /= roi-beg;
619 sum2 /= roi-beg;
620
621 if (vec[pos]>0)
622 max = vec[pos];
623
624 *(ptr+0*1440+i) = sum;
625 *(ptr+1*1440+i) = sqrt(sum2 - sum * sum);
626 *(ptr+2*1440+i) = vec[pos];
627 *(ptr+3*1440+i) = pos;
628 }
629
630 return max;
631 }
632
633 static void GetPixelMax(float *max, const float *data, uint16_t roi, int32_t first, int32_t last)
634 {
635 if (roi==0 || first<0 || last<0 || first>=roi || last>=roi || last<first)
636 return;
637
638 for (int i=0; i<1440; i++)
639 {
640 const float *beg = data+i*roi+first;
641 const float *end = data+i*roi+last;
642
643 const float *pmax = beg;
644
645 for (const float *ptr=beg+1; ptr<=end; ptr++)
646 if (*ptr>*pmax)
647 pmax = ptr;
648
649 max[i] = *pmax;
650 }
651 }
652
653 const std::vector<int64_t> &GetSum() const { return fSum; }
654
655 uint64_t GetNumEntries() const { return fNumEntries; }
656};
657
658class DrsCalibrateTime
659{
660public:
661 uint64_t fNumEntries;
662
663 size_t fNumSamples;
664 size_t fNumChannels;
665
666 std::vector<std::pair<double, double>> fStat;
667
668public:
669 DrsCalibrateTime() : fNumEntries(0), fNumSamples(0), fNumChannels(0)
670 {
671 InitSize(160, 1024);
672 }
673
674 DrsCalibrateTime(const DrsCalibrateTime &p) : fNumEntries(p.fNumEntries), fNumSamples(p.fNumSamples), fNumChannels(p.fNumChannels), fStat(p.fStat)
675 {
676 }
677 virtual ~DrsCalibrateTime()
678 {
679 }
680
681 double Sum(uint32_t i) const { return fStat[i].first; }
682 double W(uint32_t i) const { return fStat[i].second; }
683
684 virtual void InitSize(uint16_t channels, uint16_t samples)
685 {
686 fNumChannels = channels;
687 fNumSamples = samples;
688
689 fNumEntries = 0;
690
691 fStat.clear();
692
693 fStat.resize(samples*channels);
694 }
695
696 void AddT(const float *val, const int16_t *start, signed char edge=0)
697 {
698 if (fNumSamples!=1024 || fNumChannels!=160)
699 return;
700
701 // Rising or falling edge detection has the advantage that
702 // we are much less sensitive to baseline shifts
703
704 for (size_t ch=0; ch<160; ch++)
705 {
706 const size_t tm = ch*9+8;
707
708 const int16_t spos = start[tm];
709 if (spos<0)
710 continue;
711
712 const size_t pos = ch*1024;
713
714 double p_prev = 0;
715 int32_t i_prev = -1;
716
717 for (size_t i=0; i<1024-1; i++)
718 {
719 const size_t rel = tm*1024 + i;
720
721 const float &v0 = val[rel]; //-avg;
722 const float &v1 = val[rel+1];//-avg;
723
724 // Is rising edge?
725 if (edge>0 && v0>0)
726 continue;
727
728 // Is falling edge?
729 if (edge<0 && v0<0)
730 continue;
731
732 // Has sign changed?
733 if ((v0<0 && v1<0) || (v0>0 && v1>0))
734 continue;
735
736 const double p = v0==v1 ? 0.5 : v0/(v0-v1);
737
738 const double l = i+p - (i_prev+p_prev);
739
740 if (i_prev>=0)
741 {
742 const double w0 = 1-p_prev;
743 const double w1 = p;
744
745 fStat[pos+(spos+i_prev)%1024].first += w0*l;
746 fStat[pos+(spos+i_prev)%1024].second += w0;
747
748 for (size_t k=i_prev+1; k<i; k++)
749 {
750 fStat[pos+(spos+k)%1024].first += l;
751 fStat[pos+(spos+k)%1024].second += 1;
752 }
753
754 fStat[pos+(spos+i)%1024].first += w1*l;
755 fStat[pos+(spos+i)%1024].second += w1;
756 }
757
758 p_prev = p;
759 i_prev = i;
760 }
761 }
762 fNumEntries++;
763 }
764
765 void FillEmptyBins()
766 {
767 for (int ch=0; ch<160; ch++)
768 {
769 const auto beg = fStat.begin() + ch*1024;
770 const auto end = beg + 1024;
771
772 double avg = 0;
773 uint32_t num = 0;
774 for (auto it=beg; it!=end; it++)
775 {
776 if (it->second<fNumEntries-0.5)
777 continue;
778
779 avg += it->first / it->second;
780 num++;
781 }
782 avg /= num;
783
784 for (auto it=beg; it!=end; it++)
785 {
786 if (it->second>=fNumEntries-0.5)
787 continue;
788
789 // {
790 // result[i+1].first = *is2;
791 // result[i+1].second = *iw2;
792 // }
793 // else
794 // {
795 it->first = avg*fNumEntries;
796 it->second = fNumEntries;
797 // }
798 }
799 }
800 }
801
802 DrsCalibrateTime GetComplete() const
803 {
804 DrsCalibrateTime rc(*this);
805 rc.FillEmptyBins();
806 return rc;
807 }
808
809 void CalcResult()
810 {
811 for (int ch=0; ch<160; ch++)
812 {
813 const auto beg = fStat.begin() + ch*1024;
814 const auto end = beg + 1024;
815
816 // Calculate mean
817 double s = 0;
818 double w = 0;
819
820 for (auto it=beg; it!=end; it++)
821 {
822 s += it->first;
823 w += it->second;
824 }
825
826 s /= w;
827
828 double sumw = 0;
829 double sumv = 0;
830 int n = 0;
831
832 // Sums about many values are numerically less stable than
833 // just sums over less. So we do the exercise from both sides.
834 for (auto it=beg; it!=end-512; it++, n++)
835 {
836 const double valv = it->first;
837 const double valw = it->second;
838
839 it->first = sumv>0 ? n*(1-s*sumw/sumv) :0;
840
841 sumv += valv;
842 sumw += valw;
843 }
844
845 sumw = 0;
846 sumv = 0;
847 n = 1;
848
849 for (auto it=end-1; it!=beg-1+512; it--, n++)
850 {
851 const double valv = it->first;
852 const double valw = it->second;
853
854 sumv += valv;
855 sumw += valw;
856
857 it->first = sumv>0 ? n*(s*sumw/sumv-1) : 0;
858 }
859 }
860 }
861
862 DrsCalibrateTime GetResult() const
863 {
864 DrsCalibrateTime rc(*this);
865 rc.CalcResult();
866 return rc;
867 }
868
869 double Offset(uint32_t ch, double pos) const
870 {
871 const auto p = fStat.begin() + ch*1024;
872
873 const uint32_t f = floor(pos);
874
875 const double v0 = p[f].first;
876 const double v1 = p[(f+1)%1024].first;
877
878 return v0 + fmod(pos, 1)*(v1-v0);
879 }
880
881 double Calib(uint32_t ch, double pos) const
882 {
883 return pos-Offset(ch, pos);
884 }
885};
886
887struct DrsCalibration
888{
889 std::vector<int32_t> fOffset;
890 std::vector<int64_t> fGain;
891 std::vector<int64_t> fTrgOff;
892
893 uint64_t fNumOffset;
894 uint64_t fNumGain;
895 uint64_t fNumTrgOff;
896
897 uint32_t fStep;
898 uint16_t fRoi; // Region of interest for trgoff
899 uint16_t fNumTm; // Number of time marker channels in trgoff
900
901// uint16_t fDAC[8];
902
903 DrsCalibration() :
904 fOffset (1440*1024, 0),
905 fGain (1440*1024, 4096),
906 fTrgOff (1600*1024, 0),
907 fNumOffset(1),
908 fNumGain(2000),
909 fNumTrgOff(1),
910 fStep(0)
911 {
912 }
913
914 void Clear()
915 {
916 // Default gain:
917 // 0.575*[45590]*2.5V / 2^16 = 0.99999 V
918 fOffset.assign(1440*1024, 0);
919 fGain.assign (1440*1024, 4096);
920 fTrgOff.assign(1600*1024, 0);
921
922 fNumOffset = 1;
923 fNumGain = 2000;
924 fNumTrgOff = 1;
925
926 fStep = 0;
927 }
928
929 std::string ReadFitsImp(const std::string &str, std::vector<float> &vec)
930 {
931 std::fits file(str);
932 if (!file)
933 {
934 std::ostringstream msg;
935 msg << "Could not open file '" << str << "': " << strerror(errno);
936 return msg.str();
937 }
938
939 if (file.GetStr("TELESCOP")!="FACT")
940 {
941 std::ostringstream msg;
942 msg << "Reading '" << str << "' failed: Not a valid FACT file (TELESCOP not FACT in header)";
943 return msg.str();
944 }
945
946 if (!file.HasKey("STEP"))
947 {
948 std::ostringstream msg;
949 msg << "Reading '" << str << "' failed: Is not a DRS calib file (STEP not found in header)";
950 return msg.str();
951 }
952
953 if (file.GetNumRows()!=1)
954 {
955 std::ostringstream msg;
956 msg << "Reading '" << str << "' failed: Number of rows in table is not 1.";
957 return msg.str();
958 }
959
960 fStep = file.GetUInt("STEP");
961 fNumOffset = file.GetUInt("NBOFFSET");
962 fNumGain = file.GetUInt("NBGAIN");
963 fNumTrgOff = file.GetUInt("NBTRGOFF");
964 fRoi = file.GetUInt("NROI");
965 fNumTm = file.HasKey("NTM") ? file.GetUInt("NTM") : 0;
966/*
967 fDAC[0] = file.GetUInt("DAC_A");
968 fDAC[1] = file.GetUInt("DAC_B");
969 fDAC[4] = file.GetUInt("DAC_C");
970*/
971 vec.resize(1440*1024*4 + (1440+fNumTm)*fRoi*2 + 4);
972
973 float *base = vec.data();
974
975 reinterpret_cast<uint32_t*>(base)[0] = fRoi;
976
977 file.SetPtrAddress("RunNumberBaseline", base+1, 1);
978 file.SetPtrAddress("RunNumberGain", base+2, 1);
979 file.SetPtrAddress("RunNumberTriggerOffset", base+3, 1);
980 file.SetPtrAddress("BaselineMean", base+4+0*1024*1440, 1024*1440);
981 file.SetPtrAddress("BaselineRms", base+4+1*1024*1440, 1024*1440);
982 file.SetPtrAddress("GainMean", base+4+2*1024*1440, 1024*1440);
983 file.SetPtrAddress("GainRms", base+4+3*1024*1440, 1024*1440);
984 file.SetPtrAddress("TriggerOffsetMean", base+4+4*1024*1440, fRoi*1440);
985 file.SetPtrAddress("TriggerOffsetRms", base+4+4*1024*1440+fRoi*1440, fRoi*1440);
986 if (fNumTm>0)
987 {
988 file.SetPtrAddress("TriggerOffsetTMMean", base+4+4*1024*1440+ 2*fRoi*1440, fRoi*fNumTm);
989 file.SetPtrAddress("TriggerOffsetTMRms", base+4+4*1024*1440+ 2*fRoi*1440+ fRoi*fNumTm, fRoi*fNumTm);
990 }
991
992 if (!file.GetNextRow())
993 {
994 std::ostringstream msg;
995 msg << "Reading data from " << str << " failed.";
996 return msg.str();
997 }
998/*
999 fDAC[2] = fDAC[1];
1000 fDAC[4] = fDAC[1];
1001
1002 fDAC[5] = fDAC[4];
1003 fDAC[6] = fDAC[4];
1004 fDAC[7] = fDAC[4];
1005*/
1006 fOffset.resize(1024*1440);
1007 fGain.resize(1024*1440);
1008
1009 fTrgOff.resize(fRoi*(1440+fNumTm));
1010
1011 // Convert back to ADC counts: 256/125 = 4096/2000
1012 // Convert back to sum (mean * num_entries)
1013 for (int i=0; i<1024*1440; i++)
1014 {
1015 fOffset[i] = fNumOffset *256*base[i+1024*1440*0+4]/125;
1016 fGain[i] = fNumOffset*fNumGain*256*base[i+1024*1440*2+4]/125;
1017 }
1018
1019 for (int i=0; i<fRoi*1440; i++)
1020 fTrgOff[i] = fNumOffset*fNumTrgOff*256*base[i+1024*1440*4+4]/125;
1021
1022 for (int i=0; i<fRoi*fNumTm; i++)
1023 fTrgOff[i+1440*fRoi] = fNumOffset*fNumTrgOff*256*base[i+1024*1440*4+2*fRoi*1440+4]/125;
1024
1025
1026 // DAC: 0..2.5V == 0..65535
1027 // V-mV: 1000
1028 //fNumGain *= 2500*50000;
1029 //for (int i=0; i<1024*1440; i++)
1030 // fGain[i] *= 65536;
1031 if (fStep==0)
1032 {
1033 for (int i=0; i<1024*1440; i++)
1034 fGain[i] = fNumOffset*4096;
1035 }
1036 else
1037 {
1038 fNumGain *= 1953125;
1039 for (int i=0; i<1024*1440; i++)
1040 fGain[i] *= 1024;
1041 }
1042
1043 // Now mark the stored DRS data as "officially valid"
1044 // However, this is not thread safe. It only ensures that
1045 // this data is not used before it is completely and correctly
1046 // read.
1047 fStep++;
1048
1049 return std::string();
1050 }
1051
1052 std::string WriteFitsImp(const std::string &filename, const std::vector<float> &vec) const
1053 {
1054 const size_t n = 1440*1024*4 + 1440*fRoi*2 + fNumTm*fRoi*2 + 3;
1055
1056 std::ofits file(filename.c_str());
1057 if (!file)
1058 {
1059 std::ostringstream msg;
1060 msg << "Could not open file '" << filename << "': " << strerror(errno);
1061 return msg.str();
1062 }
1063
1064 file.AddColumnInt("RunNumberBaseline");
1065 file.AddColumnInt("RunNumberGain");
1066 file.AddColumnInt("RunNumberTriggerOffset");
1067
1068 file.AddColumnFloat(1024*1440, "BaselineMean", "mV");
1069 file.AddColumnFloat(1024*1440, "BaselineRms", "mV");
1070 file.AddColumnFloat(1024*1440, "GainMean", "mV");
1071 file.AddColumnFloat(1024*1440, "GainRms", "mV");
1072 file.AddColumnFloat(fRoi*1440, "TriggerOffsetMean", "mV");
1073 file.AddColumnFloat(fRoi*1440, "TriggerOffsetRms", "mV");
1074 file.AddColumnFloat(fRoi*fNumTm, "TriggerOffsetTMMean", "mV");
1075 file.AddColumnFloat(fRoi*fNumTm, "TriggerOffsetTMRms", "mV");
1076
1077#ifdef __MARS__
1078 const MTime now(-1);
1079 file.SetStr( "TELESCOP", "FACT", "Telescope that acquired this data");
1080 file.SetStr( "PACKAGE", "MARS", "Package name");
1081 file.SetStr( "VERSION", "1.0", "Package description");
1082 //file.SetStr( "CREATOR", "root", "Program that wrote this file");
1083 file.SetFloat("EXTREL", 1.0, "Release Number");
1084 file.SetStr( "COMPILED", __DATE__" "__TIME__, "Compile time");
1085 //file.SetStr( "REVISION", REVISION, "SVN revision");
1086 file.SetStr( "ORIGIN", "FACT", "Institution that wrote the file");
1087 file.SetStr( "DATE", now.GetStringFmt("%Y-%m-%dT%H:%M:%S").Data(), "File creation date");
1088 file.SetInt( "NIGHT", now.GetNightAsInt(), "Night as int");
1089 file.SetStr( "TIMESYS", "UTC", "Time system");
1090 file.SetStr( "TIMEUNIT", "d", "Time given in days w.r.t. to MJDREF");
1091 file.SetInt( "MJDREF", 40587, "MJD to UNIX time (seconds since 1970/1/1)");
1092#else
1093 DataWriteFits2::WriteDefaultKeys(file);
1094#endif
1095 file.SetInt("STEP", fStep, "");
1096
1097 file.SetInt("ADCRANGE", 2000, "Dynamic range of the ADC in mV");
1098 file.SetInt("DACRANGE", 2500, "Dynamic range of the DAC in mV");
1099 file.SetInt("ADC", 12, "Resolution of ADC in bits");
1100 file.SetInt("DAC", 16, "Resolution of DAC in bits");
1101 file.SetInt("NPIX", 1440, "Number of channels in the camera");
1102 file.SetInt("NTM", fNumTm, "Number of time marker channels");
1103 file.SetInt("NROI", fRoi, "Region of interest");
1104
1105 file.SetInt("NBOFFSET", fNumOffset, "Num of entries for offset calibration");
1106 file.SetInt("NBGAIN", fNumGain/1953125, "Num of entries for gain calibration");
1107 file.SetInt("NBTRGOFF", fNumTrgOff, "Num of entries for trigger offset calibration");
1108
1109 // file.WriteKeyNT("DAC_A", fData.fDAC[0], "Level of DAC 0 in DAC counts") ||
1110 // file.WriteKeyNT("DAC_B", fData.fDAC[1], "Leval of DAC 1-3 in DAC counts") ||
1111 // file.WriteKeyNT("DAC_C", fData.fDAC[4], "Leval of DAC 4-7 in DAC counts") ||
1112
1113 file.WriteTableHeader("DrsCalibration");
1114
1115 if (!file.WriteRow(vec.data()+1, n*sizeof(float)))
1116 {
1117 std::ostringstream msg;
1118 msg << "Writing data to " << filename << " failed.";
1119 return msg.str();
1120 }
1121
1122 return std::string();
1123 }
1124
1125
1126 std::string ReadFitsImp(const std::string &str)
1127 {
1128 std::vector<float> vec;
1129 return ReadFitsImp(str, vec);
1130 }
1131
1132 bool IsValid() { return fStep>2; }
1133
1134 bool Apply(float *vec, const int16_t *val, const int16_t *start, uint32_t roi)
1135 {
1136 if (roi!=fRoi)
1137 {
1138 for (size_t ch=0; ch<1440; ch++)
1139 {
1140 const size_t pos = ch*roi;
1141 const size_t drs = ch*1024;
1142
1143 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1144 fOffset.data()+drs, fNumOffset,
1145 fGain.data() +drs, fNumGain);
1146 }
1147
1148 return false;
1149 }
1150
1151 for (size_t ch=0; ch<1440; ch++)
1152 {
1153 const size_t pos = ch*fRoi;
1154 const size_t drs = ch*1024;
1155
1156 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1157 fOffset.data()+drs, fNumOffset,
1158 fGain.data() +drs, fNumGain,
1159 fTrgOff.data()+pos, fNumTrgOff);
1160 }
1161
1162 for (size_t ch=0; ch<fNumTm; ch++)
1163 {
1164 const size_t pos = (ch+1440)*fRoi;
1165 const size_t drs = (ch*9+8)*1024;
1166
1167 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1168 fOffset.data()+drs, fNumOffset,
1169 fGain.data() +drs, fNumGain,
1170 fTrgOff.data()+pos, fNumTrgOff);
1171 }
1172
1173 return true;
1174 }
1175};
1176
1177#endif
Note: See TracBrowser for help on using the repository browser.