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

Last change on this file since 14631 was 14619, checked in by tbretz, 12 years ago
Fixed a bug in GetPixelStats... the first value was added twice
File size: 30.7 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> &list, uint32_t n)
281 {
282 Step rc;
283 for (auto it=list.begin(); it!=list.begin()+n; it++)
284 {
285 rc.pos += it->pos;
286 rc.avg += it->avg;
287 rc.rms += it->avg*it->avg;
288 }
289
290 rc.cnt = n;
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, list.size());;
332
333 if (rc.avg==0)
334 return Step();
335
336 if (rc.rms>5)
337 {
338 partial_sort(list.begin(), list.begin()+list.size()/2,
339 list.end(), Step::sort);
340
341 rc = AverageSteps(list, list.size()/2);
342 }
343
344 for (size_t ch=0; ch<nch; ch += 9)
345 {
346 const int16_t dist = (prev[ch]-start[ch]+1024+offset)%1024;
347 SubtractStep(ch, rc.avg, vec, roi, dist, map);
348 }
349
350 return rc;
351 }
352
353 static void RemoveSpikes(float *vec, uint32_t roi)
354 {
355 if (roi<4)
356 return;
357
358 for (size_t ch=0; ch<1440; ch++)
359 {
360 float *p = vec + ch*roi;
361
362 for (size_t i=1; i<roi-2; i++)
363 {
364 if (p[i]-p[i-1]>25 && p[i]-p[i+1]>25)
365 {
366 p[i] = (p[i-1]+p[i+1])/2;
367 }
368
369 if (p[i]-p[i-1]>22 && fabs(p[i]-p[i+1])<4 && p[i+1]-p[i+2]>22)
370 {
371 p[i] = (p[i-1]+p[i+2])/2;
372 p[i+1] = p[i];
373 }
374 }
375 }
376 }
377
378 static void RemoveSpikes2(float *vec, uint32_t roi)//from Werner
379 {
380 if (roi<4)
381 return;
382
383 for (size_t ch=0; ch<1440; ch++)
384 {
385 float *p = vec + ch*roi;
386
387 std::vector<float> Ameas(p, p+roi);
388
389 std::vector<float> N1mean(roi);
390 for (size_t i=2; i<roi-2; i++)
391 {
392 N1mean[i] = (p[i-1] + p[i+1])/2.;
393 }
394
395 const float fract = 0.8;
396 float x, xp, xpp;
397
398 for (size_t i=0; i<roi-3; i++)
399 {
400 x = Ameas[i] - N1mean[i];
401 if ( x > -5.)
402 continue;
403
404 xp = Ameas[i+1] - N1mean[i+1];
405 xpp = Ameas[i+2] - N1mean[i+2];
406
407 if(Ameas[i+2] - (Ameas[i] + Ameas[i+3])/2. > 10.)
408 {
409 p[i+1]=(Ameas[i] + Ameas[i+3])/2.;
410 p[i+2]=(Ameas[i] + Ameas[i+3])/2.;
411
412 i += 3;
413
414 continue;
415 }
416 if ( (xp > -2.*x*fract) && (xpp < -10.) )
417 {
418 p[i+1] = N1mean[i+1];
419 N1mean[i+2] = (Ameas[i+1] - Ameas[i+3] / 2.);
420
421 i += 2;
422 }
423 }
424 }
425 }
426
427 static void SlidingAverage(float *const vec, const uint32_t roi, const uint16_t w)
428 {
429 if (w==0)
430 return;
431
432 if (w>roi)
433 return;
434
435 for (float *pix=vec; pix<vec+1440*roi; pix += roi)
436 {
437 for (float *ptr=pix; ptr<pix+roi-w; ptr++)
438 {
439 for (float *p=ptr+1; p<ptr+w; p++)
440 *ptr += *p;
441 *ptr /= w;
442 }
443 }
444 }
445
446 std::pair<std::vector<double>,std::vector<double> > GetSampleStats() const
447 {
448 if (fNumEntries==0)
449 return make_pair(std::vector<double>(),std::vector<double>());
450
451 std::vector<double> mean(fSum.size());
452 std::vector<double> error(fSum.size());
453
454 std::vector<int64_t>::const_iterator it = fSum.begin();
455 std::vector<int64_t>::const_iterator i2 = fSum2.begin();
456 std::vector<double>::iterator im = mean.begin();
457 std::vector<double>::iterator ie = error.begin();
458
459 while (it!=fSum.end())
460 {
461 *im = /*cnt<fResult.size() ? fResult[cnt] :*/ double(*it)/fNumEntries;
462 *ie = sqrt(double(*i2*int64_t(fNumEntries) - *it * *it))/fNumEntries;
463
464 im++;
465 ie++;
466 it++;
467 i2++;
468 }
469
470
471 /*
472 valarray<double> ...
473
474 mean /= fNumEntries;
475 error = sqrt(error/fNumEntries - mean*mean);
476 */
477
478 return make_pair(mean, error);
479 }
480
481 void GetSampleStats(float *ptr, float scale) const
482 {
483 const size_t sz = fNumSamples*fNumChannels;
484
485 if (fNumEntries==0)
486 {
487 memset(ptr, 0, sizeof(float)*sz*2);
488 return;
489 }
490
491 std::vector<int64_t>::const_iterator it = fSum.begin();
492 std::vector<int64_t>::const_iterator i2 = fSum2.begin();
493
494 while (it!=fSum.end())
495 {
496 *ptr = scale*double(*it)/fNumEntries;
497 *(ptr+sz) = scale*sqrt(double(*i2*int64_t(fNumEntries) - *it * *it))/fNumEntries;
498
499 ptr++;
500 it++;
501 i2++;
502 }
503 }
504
505 static double GetPixelStats(float *ptr, const float *data, uint16_t roi)
506 {
507 if (roi==0)
508 return -1;
509
510 const int beg = roi>10 ? 10 : 0;
511
512 double max = 0;
513 for (int i=0; i<1440; i++)
514 {
515 const float *vec = data+i*roi;
516
517 int pos = beg;
518 double sum = vec[beg];
519 double sum2 = vec[beg]*vec[beg];
520 for (int j=beg+1; j<roi; j++)
521 {
522 sum += vec[j];
523 sum2 += vec[j]*vec[j];
524
525 if (vec[j]>vec[pos])
526 pos = j;
527 }
528 sum /= roi-beg;
529 sum2 /= roi-beg;
530
531 if (vec[pos]>0)
532 max = vec[pos];
533
534 *(ptr+0*1440+i) = sum;
535 *(ptr+1*1440+i) = sqrt(sum2 - sum * sum);
536 *(ptr+2*1440+i) = vec[pos];
537 *(ptr+3*1440+i) = pos;
538 }
539
540 return max;
541 }
542
543 static void GetPixelMax(float *max, const float *data, uint16_t roi, int32_t first, int32_t last)
544 {
545 if (roi==0 || first<0 || last<0 || first>=roi || last>=roi || last<first)
546 return;
547
548 for (int i=0; i<1440; i++)
549 {
550 const float *beg = data+i*roi+first;
551 const float *end = data+i*roi+last;
552
553 const float *pmax = beg;
554
555 for (const float *ptr=beg+1; ptr<=end; ptr++)
556 if (*ptr>*pmax)
557 pmax = ptr;
558
559 max[i] = *pmax;
560 }
561 }
562
563 const std::vector<int64_t> &GetSum() const { return fSum; }
564
565 uint64_t GetNumEntries() const { return fNumEntries; }
566};
567
568class DrsCalibrateTime
569{
570public:
571 uint64_t fNumEntries;
572
573 size_t fNumSamples;
574 size_t fNumChannels;
575
576 std::vector<std::pair<double, double>> fStat;
577
578public:
579 DrsCalibrateTime() : fNumEntries(0), fNumSamples(0), fNumChannels(0)
580 {
581 InitSize(160, 1024);
582 }
583
584 DrsCalibrateTime(const DrsCalibrateTime &p) : fNumEntries(p.fNumEntries), fNumSamples(p.fNumSamples), fNumChannels(p.fNumChannels), fStat(p.fStat)
585 {
586 }
587 virtual ~DrsCalibrateTime()
588 {
589 }
590
591 double Sum(uint32_t i) const { return fStat[i].first; }
592 double W(uint32_t i) const { return fStat[i].second; }
593
594 virtual void InitSize(uint16_t channels, uint16_t samples)
595 {
596 fNumChannels = channels;
597 fNumSamples = samples;
598
599 fNumEntries = 0;
600
601 fStat.clear();
602
603 fStat.resize(samples*channels);
604 }
605
606 void AddT(const float *val, const int16_t *start, signed char edge=0)
607 {
608 if (fNumSamples!=1024 || fNumChannels!=160)
609 return;
610
611 // Rising or falling edge detection has the advantage that
612 // we are much less sensitive to baseline shifts
613
614 for (size_t ch=0; ch<160; ch++)
615 {
616 const size_t tm = ch*9+8;
617
618 const int16_t spos = start[tm];
619 if (spos<0)
620 continue;
621
622 const size_t pos = ch*1024;
623
624 double p_prev = 0;
625 int32_t i_prev = -1;
626
627 for (size_t i=0; i<1024-1; i++)
628 {
629 const size_t rel = tm*1024 + i;
630
631 const float &v0 = val[rel]; //-avg;
632 const float &v1 = val[rel+1];//-avg;
633
634 // Is rising edge?
635 if (edge>0 && v0>0)
636 continue;
637
638 // Is falling edge?
639 if (edge<0 && v0<0)
640 continue;
641
642 // Has sign changed?
643 if ((v0<0 && v1<0) || (v0>0 && v1>0))
644 continue;
645
646 const double p = v0==v1 ? 0.5 : v0/(v0-v1);
647
648 const double l = i+p - (i_prev+p_prev);
649
650 if (i_prev>=0)
651 {
652 const double w0 = 1-p_prev;
653 const double w1 = p;
654
655 fStat[pos+(spos+i_prev)%1024].first += w0*l;
656 fStat[pos+(spos+i_prev)%1024].second += w0;
657
658 for (size_t k=i_prev+1; k<i; k++)
659 {
660 fStat[pos+(spos+k)%1024].first += l;
661 fStat[pos+(spos+k)%1024].second += 1;
662 }
663
664 fStat[pos+(spos+i)%1024].first += w1*l;
665 fStat[pos+(spos+i)%1024].second += w1;
666 }
667
668 p_prev = p;
669 i_prev = i;
670 }
671 }
672 fNumEntries++;
673 }
674
675 void FillEmptyBins()
676 {
677 for (int ch=0; ch<160; ch++)
678 {
679 const auto beg = fStat.begin() + ch*1024;
680 const auto end = beg + 1024;
681
682 double avg = 0;
683 uint32_t num = 0;
684 for (auto it=beg; it!=end; it++)
685 {
686 if (it->second<fNumEntries-0.5)
687 continue;
688
689 avg += it->first / it->second;
690 num++;
691 }
692 avg /= num;
693
694 for (auto it=beg; it!=end; it++)
695 {
696 if (it->second>=fNumEntries-0.5)
697 continue;
698
699 // {
700 // result[i+1].first = *is2;
701 // result[i+1].second = *iw2;
702 // }
703 // else
704 // {
705 it->first = avg*fNumEntries;
706 it->second = fNumEntries;
707 // }
708 }
709 }
710 }
711
712 DrsCalibrateTime GetComplete() const
713 {
714 DrsCalibrateTime rc(*this);
715 rc.FillEmptyBins();
716 return rc;
717 }
718
719 void CalcResult()
720 {
721 for (int ch=0; ch<160; ch++)
722 {
723 const auto beg = fStat.begin() + ch*1024;
724 const auto end = beg + 1024;
725
726 // Calculate mean
727 double s = 0;
728 double w = 0;
729
730 for (auto it=beg; it!=end; it++)
731 {
732 s += it->first;
733 w += it->second;
734 }
735
736 s /= w;
737
738 double sumw = 0;
739 double sumv = 0;
740 int n = 0;
741
742 // Sums about many values are numerically less stable than
743 // just sums over less. So we do the exercise from both sides.
744 for (auto it=beg; it!=end-512; it++, n++)
745 {
746 const double valv = it->first;
747 const double valw = it->second;
748
749 it->first = sumv>0 ? n*(1-s*sumw/sumv) :0;
750
751 sumv += valv;
752 sumw += valw;
753 }
754
755 sumw = 0;
756 sumv = 0;
757 n = 1;
758
759 for (auto it=end-1; it!=beg-1+512; it--, n++)
760 {
761 const double valv = it->first;
762 const double valw = it->second;
763
764 sumv += valv;
765 sumw += valw;
766
767 it->first = sumv>0 ? n*(s*sumw/sumv-1) : 0;
768 }
769 }
770 }
771
772 DrsCalibrateTime GetResult() const
773 {
774 DrsCalibrateTime rc(*this);
775 rc.CalcResult();
776 return rc;
777 }
778
779 double Offset(uint32_t ch, double pos) const
780 {
781 const auto p = fStat.begin() + ch*1024;
782
783 const uint32_t f = floor(pos);
784
785 const double v0 = p[f].first;
786 const double v1 = p[(f+1)%1024].first;
787
788 return v0 + fmod(pos, 1)*(v1-v0);
789 }
790
791 double Calib(uint32_t ch, double pos) const
792 {
793 return pos-Offset(ch, pos);
794 }
795};
796
797struct DrsCalibration
798{
799 std::vector<int32_t> fOffset;
800 std::vector<int64_t> fGain;
801 std::vector<int64_t> fTrgOff;
802
803 uint64_t fNumOffset;
804 uint64_t fNumGain;
805 uint64_t fNumTrgOff;
806
807 uint32_t fStep;
808 uint16_t fRoi; // Region of interest for trgoff
809 uint16_t fNumTm; // Number of time marker channels in trgoff
810
811// uint16_t fDAC[8];
812
813 DrsCalibration() :
814 fOffset (1440*1024, 0),
815 fGain (1440*1024, 4096),
816 fTrgOff (1600*1024, 0),
817 fNumOffset(1),
818 fNumGain(2000),
819 fNumTrgOff(1),
820 fStep(0)
821 {
822 }
823
824 void Clear()
825 {
826 // Default gain:
827 // 0.575*[45590]*2.5V / 2^16 = 0.99999 V
828 fOffset.assign(1440*1024, 0);
829 fGain.assign (1440*1024, 4096);
830 fTrgOff.assign(1600*1024, 0);
831
832 fNumOffset = 1;
833 fNumGain = 2000;
834 fNumTrgOff = 1;
835
836 fStep = 0;
837 }
838
839 std::string ReadFitsImp(const std::string &str, std::vector<float> &vec)
840 {
841 std::fits file(str);
842 if (!file)
843 {
844 std::ostringstream msg;
845 msg << "Could not open file '" << str << "': " << strerror(errno);
846 return msg.str();
847 }
848
849 if (file.GetStr("TELESCOP")!="FACT")
850 {
851 std::ostringstream msg;
852 msg << "Reading '" << str << "' failed: Not a valid FACT file (TELESCOP not FACT in header)";
853 return msg.str();
854 }
855
856 if (!file.HasKey("STEP"))
857 {
858 std::ostringstream msg;
859 msg << "Reading '" << str << "' failed: Is not a DRS calib file (STEP not found in header)";
860 return msg.str();
861 }
862
863 if (file.GetNumRows()!=1)
864 {
865 std::ostringstream msg;
866 msg << "Reading '" << str << "' failed: Number of rows in table is not 1.";
867 return msg.str();
868 }
869
870 fStep = file.GetUInt("STEP");
871 fNumOffset = file.GetUInt("NBOFFSET");
872 fNumGain = file.GetUInt("NBGAIN");
873 fNumTrgOff = file.GetUInt("NBTRGOFF");
874 fRoi = file.GetUInt("NROI");
875 fNumTm = file.HasKey("NTM") ? file.GetUInt("NTM") : 0;
876/*
877 fDAC[0] = file.GetUInt("DAC_A");
878 fDAC[1] = file.GetUInt("DAC_B");
879 fDAC[4] = file.GetUInt("DAC_C");
880*/
881 vec.resize(1440*1024*4 + (1440+fNumTm)*fRoi*2 + 4);
882
883 float *base = vec.data();
884
885 reinterpret_cast<uint32_t*>(base)[0] = fRoi;
886
887 file.SetPtrAddress("RunNumberBaseline", base+1, 1);
888 file.SetPtrAddress("RunNumberGain", base+2, 1);
889 file.SetPtrAddress("RunNumberTriggerOffset", base+3, 1);
890 file.SetPtrAddress("BaselineMean", base+4+0*1024*1440, 1024*1440);
891 file.SetPtrAddress("BaselineRms", base+4+1*1024*1440, 1024*1440);
892 file.SetPtrAddress("GainMean", base+4+2*1024*1440, 1024*1440);
893 file.SetPtrAddress("GainRms", base+4+3*1024*1440, 1024*1440);
894 file.SetPtrAddress("TriggerOffsetMean", base+4+4*1024*1440, fRoi*1440);
895 file.SetPtrAddress("TriggerOffsetRms", base+4+4*1024*1440+fRoi*1440, fRoi*1440);
896 if (fNumTm>0)
897 {
898 file.SetPtrAddress("TriggerOffsetTMMean", base+4+4*1024*1440+ 2*fRoi*1440, fRoi*fNumTm);
899 file.SetPtrAddress("TriggerOffsetTMRms", base+4+4*1024*1440+ 2*fRoi*1440+ fRoi*fNumTm, fRoi*fNumTm);
900 }
901
902 if (!file.GetNextRow())
903 {
904 std::ostringstream msg;
905 msg << "Reading data from " << str << " failed.";
906 return msg.str();
907 }
908/*
909 fDAC[2] = fDAC[1];
910 fDAC[4] = fDAC[1];
911
912 fDAC[5] = fDAC[4];
913 fDAC[6] = fDAC[4];
914 fDAC[7] = fDAC[4];
915*/
916 fOffset.resize(1024*1440);
917 fGain.resize(1024*1440);
918
919 fTrgOff.resize(fRoi*(1440+fNumTm));
920
921 // Convert back to ADC counts: 256/125 = 4096/2000
922 // Convert back to sum (mean * num_entries)
923 for (int i=0; i<1024*1440; i++)
924 {
925 fOffset[i] = fNumOffset *256*base[i+1024*1440*0+4]/125;
926 fGain[i] = fNumOffset*fNumGain*256*base[i+1024*1440*2+4]/125;
927 }
928
929 for (int i=0; i<fRoi*1440; i++)
930 fTrgOff[i] = fNumOffset*fNumTrgOff*256*base[i+1024*1440*4+4]/125;
931
932 for (int i=0; i<fRoi*fNumTm; i++)
933 fTrgOff[i+1440*fRoi] = fNumOffset*fNumTrgOff*256*base[i+1024*1440*4+2*fRoi*1440+4]/125;
934
935
936 // DAC: 0..2.5V == 0..65535
937 // V-mV: 1000
938 //fNumGain *= 2500*50000;
939 //for (int i=0; i<1024*1440; i++)
940 // fGain[i] *= 65536;
941 if (fStep==0)
942 {
943 for (int i=0; i<1024*1440; i++)
944 fGain[i] = fNumOffset*4096;
945 }
946 else
947 {
948 fNumGain *= 1953125;
949 for (int i=0; i<1024*1440; i++)
950 fGain[i] *= 1024;
951 }
952
953 // Now mark the stored DRS data as "officially valid"
954 // However, this is not thread safe. It only ensures that
955 // this data is not used before it is completely and correctly
956 // read.
957 fStep++;
958
959 return std::string();
960 }
961
962 std::string WriteFitsImp(const std::string &filename, const std::vector<float> &vec) const
963 {
964 const size_t n = 1440*1024*4 + 1440*fRoi*2 + fNumTm*fRoi*2 + 3;
965
966 std::ofits file(filename.c_str());
967 if (!file)
968 {
969 std::ostringstream msg;
970 msg << "Could not open file '" << filename << "': " << strerror(errno);
971 return msg.str();
972 }
973
974 file.AddColumnInt("RunNumberBaseline");
975 file.AddColumnInt("RunNumberGain");
976 file.AddColumnInt("RunNumberTriggerOffset");
977
978 file.AddColumnFloat(1024*1440, "BaselineMean", "mV");
979 file.AddColumnFloat(1024*1440, "BaselineRms", "mV");
980 file.AddColumnFloat(1024*1440, "GainMean", "mV");
981 file.AddColumnFloat(1024*1440, "GainRms", "mV");
982 file.AddColumnFloat(fRoi*1440, "TriggerOffsetMean", "mV");
983 file.AddColumnFloat(fRoi*1440, "TriggerOffsetRms", "mV");
984 file.AddColumnFloat(fRoi*fNumTm, "TriggerOffsetTMMean", "mV");
985 file.AddColumnFloat(fRoi*fNumTm, "TriggerOffsetTMRms", "mV");
986
987#ifdef __MARS__
988 const MTime now(-1);
989 file.SetStr( "TELESCOP", "FACT", "Telescope that acquired this data");
990 file.SetStr( "PACKAGE", "MARS", "Package name");
991 file.SetStr( "VERSION", "1.0", "Package description");
992 //file.SetStr( "CREATOR", "root", "Program that wrote this file");
993 file.SetFloat("EXTREL", 1.0, "Release Number");
994 file.SetStr( "COMPILED", __DATE__" "__TIME__, "Compile time");
995 //file.SetStr( "REVISION", REVISION, "SVN revision");
996 file.SetStr( "ORIGIN", "FACT", "Institution that wrote the file");
997 file.SetStr( "DATE", now.GetStringFmt("%Y-%m-%dT%H:%M:%S").Data(), "File creation date");
998 file.SetInt( "NIGHT", now.GetNightAsInt(), "Night as int");
999 file.SetStr( "TIMESYS", "UTC", "Time system");
1000 file.SetStr( "TIMEUNIT", "d", "Time given in days w.r.t. to MJDREF");
1001 file.SetInt( "MJDREF", 40587, "MJD to UNIX time (seconds since 1970/1/1)");
1002#else
1003 DataWriteFits2::WriteDefaultKeys(file);
1004#endif
1005 file.SetInt("STEP", fStep, "");
1006
1007 file.SetInt("ADCRANGE", 2000, "Dynamic range of the ADC in mV");
1008 file.SetInt("DACRANGE", 2500, "Dynamic range of the DAC in mV");
1009 file.SetInt("ADC", 12, "Resolution of ADC in bits");
1010 file.SetInt("DAC", 16, "Resolution of DAC in bits");
1011 file.SetInt("NPIX", 1440, "Number of channels in the camera");
1012 file.SetInt("NTM", fNumTm, "Number of time marker channels");
1013 file.SetInt("NROI", fRoi, "Region of interest");
1014
1015 file.SetInt("NBOFFSET", fNumOffset, "Num of entries for offset calibration");
1016 file.SetInt("NBGAIN", fNumGain/1953125, "Num of entries for gain calibration");
1017 file.SetInt("NBTRGOFF", fNumTrgOff, "Num of entries for trigger offset calibration");
1018
1019 // file.WriteKeyNT("DAC_A", fData.fDAC[0], "Level of DAC 0 in DAC counts") ||
1020 // file.WriteKeyNT("DAC_B", fData.fDAC[1], "Leval of DAC 1-3 in DAC counts") ||
1021 // file.WriteKeyNT("DAC_C", fData.fDAC[4], "Leval of DAC 4-7 in DAC counts") ||
1022
1023 file.WriteTableHeader("DrsCalibration");
1024
1025 if (!file.WriteRow(vec.data()+1, n*sizeof(float)))
1026 {
1027 std::ostringstream msg;
1028 msg << "Writing data to " << filename << " failed.";
1029 return msg.str();
1030 }
1031
1032 return std::string();
1033 }
1034
1035
1036 std::string ReadFitsImp(const std::string &str)
1037 {
1038 std::vector<float> vec;
1039 return ReadFitsImp(str, vec);
1040 }
1041
1042 bool IsValid() { return fStep>2; }
1043
1044 bool Apply(float *vec, const int16_t *val, const int16_t *start, uint32_t roi)
1045 {
1046 if (roi!=fRoi)
1047 {
1048 for (size_t ch=0; ch<1440; ch++)
1049 {
1050 const size_t pos = ch*roi;
1051 const size_t drs = ch*1024;
1052
1053 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1054 fOffset.data()+drs, fNumOffset,
1055 fGain.data() +drs, fNumGain);
1056 }
1057
1058 return false;
1059 }
1060
1061 for (size_t ch=0; ch<1440; ch++)
1062 {
1063 const size_t pos = ch*fRoi;
1064 const size_t drs = ch*1024;
1065
1066 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1067 fOffset.data()+drs, fNumOffset,
1068 fGain.data() +drs, fNumGain,
1069 fTrgOff.data()+pos, fNumTrgOff);
1070 }
1071
1072 for (size_t ch=0; ch<fNumTm; ch++)
1073 {
1074 const size_t pos = (ch+1440)*fRoi;
1075 const size_t drs = (ch*9+8)*1024;
1076
1077 DrsCalibrate::ApplyCh(vec+pos, val+pos, start[ch], roi,
1078 fOffset.data()+drs, fNumOffset,
1079 fGain.data() +drs, fNumGain,
1080 fTrgOff.data()+pos, fNumTrgOff);
1081 }
1082
1083 return true;
1084 }
1085};
1086
1087#endif
Note: See TracBrowser for help on using the repository browser.