source: trunk/MagicSoft/Mars/msignal/MExtractTimeAndChargeDigitalFilter.cc@ 5845

Last change on this file since 5845 was 5833, checked in by gaug, 21 years ago
*** empty log message ***
File size: 33.0 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Hendrik Bartko, 09/2004 <mailto:hbartko@mppmu.mpg.de>
19! Author(s): Markus Gaug, 05/2004 <mailto:markus@ifae.es>
20! Author(s): Diego Tescaro, 05/2004 <mailto:tescaro@pd.infn.it>
21!
22! Copyright: MAGIC Software Development, 2000-2004
23!
24!
25\* ======================================================================== */
26//////////////////////////////////////////////////////////////////////////////
27//
28// MExtractTimeAndChargeDigitalFilter
29//
30// Hendrik has promised to write more documentation
31//
32//
33// The following variables have to be set by the derived class and
34// do not have defaults:
35// - fNumHiGainSamples
36// - fNumLoGainSamples
37// - fSqrtHiGainSamples
38// - fSqrtLoGainSamples
39//
40// Input Containers:
41// MRawEvtData
42// MRawRunHeader
43// MPedestalCam
44//
45// Output Containers:
46// MArrivalTimeCam
47// MExtractedSignalCam
48//
49//////////////////////////////////////////////////////////////////////////////
50#include "MExtractTimeAndChargeDigitalFilter.h"
51
52#include <errno.h>
53#include <fstream>
54
55#include <TFile.h>
56#include <TH1F.h>
57#include <TH2F.h>
58#include <TString.h>
59#include <TMatrix.h>
60
61#include "MLog.h"
62#include "MLogManip.h"
63
64#include "MPedestalPix.h"
65
66ClassImp(MExtractTimeAndChargeDigitalFilter);
67
68using namespace std;
69
70const Byte_t MExtractTimeAndChargeDigitalFilter::fgHiGainFirst = 0;
71const Byte_t MExtractTimeAndChargeDigitalFilter::fgHiGainLast = 14;
72const Byte_t MExtractTimeAndChargeDigitalFilter::fgLoGainFirst = 3;
73const Byte_t MExtractTimeAndChargeDigitalFilter::fgLoGainLast = 14;
74const Int_t MExtractTimeAndChargeDigitalFilter::fgWindowSizeHiGain = 6;
75const Int_t MExtractTimeAndChargeDigitalFilter::fgWindowSizeLoGain = 6;
76const Int_t MExtractTimeAndChargeDigitalFilter::fgBinningResolutionHiGain = 10;
77const Int_t MExtractTimeAndChargeDigitalFilter::fgBinningResolutionLoGain = 10;
78const Int_t MExtractTimeAndChargeDigitalFilter::fgSignalStartBinHiGain = 4;
79const Int_t MExtractTimeAndChargeDigitalFilter::fgSignalStartBinLoGain = 4;
80const TString MExtractTimeAndChargeDigitalFilter::fgNameWeightsFile = "msignal/cosmics_weights.dat";
81const Float_t MExtractTimeAndChargeDigitalFilter::fgOffsetLoGain = 1.8; // 5 ns
82// --------------------------------------------------------------------------
83//
84// Default constructor.
85//
86// Calls:
87// - SetWindowSize();
88// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
89// - SetBinningResolution();
90//
91// Sets all weights to 1.
92//
93MExtractTimeAndChargeDigitalFilter::MExtractTimeAndChargeDigitalFilter(const char *name, const char *title)
94 : fWeightsSet(kFALSE), fRandomIter(0)
95{
96 fName = name ? name : "MExtractTimeAndChargeDigitalFilter";
97 fTitle = title ? title : "Digital Filter";
98
99 SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
100 SetWindowSize();
101 SetBinningResolution();
102 SetSignalStartBin();
103
104 SetNameWeightsFile();
105 SetOffsetLoGain(fgOffsetLoGain);
106}
107
108// ---------------------------------------------------------------------------------------
109//
110// Checks:
111// - if a window is bigger than the one defined by the ranges, set it to the available range
112//
113// Sets:
114// - fNumHiGainSamples to: (Float_t)fWindowSizeHiGain
115// - fNumLoGainSamples to: (Float_t)fWindowSizeLoGain
116//
117void MExtractTimeAndChargeDigitalFilter::SetWindowSize(Int_t windowh, Int_t windowl)
118{
119
120 if (windowh != fgWindowSizeHiGain)
121 *fLog << warn << GetDescriptor()
122 << ": ATTENTION!!! If you are not Hendrik Bartko, do NOT use a different window size than the default." << endl;
123 if (windowl != fgWindowSizeLoGain)
124 *fLog << warn << GetDescriptor()
125 << ": ATTENTION!!! If you are not Hendrik Bartko, do NOT use a different window size than the default" << endl;
126
127 fWindowSizeHiGain = windowh;
128 fWindowSizeLoGain = windowl;
129
130 const Int_t availhirange = (Int_t)(fHiGainLast-fHiGainFirst+1);
131
132 if (fWindowSizeHiGain > availhirange)
133 {
134 // Please simplify this!
135 *fLog << warn << GetDescriptor()
136 << Form("%s%2i%s%2i%s%2i%s",": Hi Gain window size: ",fWindowSizeHiGain,
137 " is bigger than available range: [",(int)fHiGainFirst,",",(int)fHiGainLast,"]") << endl;
138 fHiGainLast = fHiGainFirst + fWindowSizeHiGain;
139 *fLog << warn << GetDescriptor()
140 << ": Will set the upper range to: " << (int)fHiGainLast << endl;
141 }
142
143 if (fWindowSizeHiGain < 2)
144 {
145 fWindowSizeHiGain = 2;
146 *fLog << warn << GetDescriptor() << ": High Gain window size set to two samples" << endl;
147 }
148
149 if (fLoGainLast != 0 && fWindowSizeLoGain != 0)
150 {
151 const Int_t availlorange = (Int_t)(fLoGainLast-fLoGainFirst+1);
152
153 if (fWindowSizeLoGain > availlorange)
154 {
155 // Please simplify this!
156 *fLog << warn << GetDescriptor()
157 << Form("%s%2i%s%2i%s%2i%s",": Lo Gain window size: ",fWindowSizeLoGain,
158 " is bigger than available range: [",(int)fLoGainFirst,",",(int)fLoGainLast,"]") << endl;
159 fLoGainLast = fLoGainFirst + fWindowSizeLoGain;
160 *fLog << warn << GetDescriptor()
161 << ": Will set the upper range to: " << (int)fLoGainLast << endl;
162 }
163
164 if (fWindowSizeLoGain<2)
165 {
166 fWindowSizeLoGain = 2;
167 *fLog << warn << GetDescriptor() << ": Low Gain window size set to two samples" << endl;
168 }
169 }
170 //
171 // We need here the effective number of samples which is about 2.5 in the case of a window
172 // size of 6. The exact numbers have to be found still.
173 //
174 fNumHiGainSamples = (Float_t)fWindowSizeHiGain/2.4;
175 fNumLoGainSamples = (Float_t)fWindowSizeLoGain/2.4;
176 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
177 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
178
179}
180
181// --------------------------------------------------------------------------
182//
183// InitArrays
184//
185// Gets called in the ReInit() and initialized the arrays
186//
187Bool_t MExtractTimeAndChargeDigitalFilter::InitArrays()
188{
189
190 Int_t range = (Int_t)(fHiGainLast - fHiGainFirst + 1 + fHiLoLast);
191
192 fHiGainSignal.Set(range);
193
194 range = (Int_t)(fLoGainLast - fLoGainFirst + 1);
195
196 fLoGainSignal.Set(range);
197
198 if (!fWeightsSet)
199 if (!ReadWeightsFile(fNameWeightsFile))
200 return kFALSE;
201
202 fTimeShiftHiGain = (Float_t)fHiGainFirst + 0.5 + 1./fBinningResolutionHiGain;
203 fTimeShiftLoGain = 0.5 + 1./fBinningResolutionLoGain;
204 //
205 // We need here the effective number of samples which is about 2.5 in the case of a window
206 // size of 6. The exact numbers have to be found still.
207 //
208 fNumHiGainSamples = (Float_t)fWindowSizeHiGain/2.4;
209 fNumLoGainSamples = (Float_t)fWindowSizeLoGain/2.4;
210 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
211 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
212
213 return kTRUE;
214}
215
216void MExtractTimeAndChargeDigitalFilter::CalcBinningResArrays()
217{
218
219 fArrBinningResHiGain.Set(fWindowSizeHiGain);
220 fArrBinningResHalfHiGain.Set(fWindowSizeHiGain);
221
222 for (int i=0; i<fWindowSizeHiGain; i++)
223 {
224 fArrBinningResHiGain[i] = fBinningResolutionHiGain*i;
225 fArrBinningResHalfHiGain[i] = fArrBinningResHiGain[i] + fBinningResolutionHalfHiGain;
226 }
227
228 fArrBinningResLoGain.Set(fWindowSizeLoGain);
229 fArrBinningResHalfLoGain.Set(fWindowSizeLoGain);
230
231 for (int i=0; i<fWindowSizeLoGain; i++)
232 {
233 fArrBinningResLoGain[i] = fBinningResolutionLoGain*i;
234 fArrBinningResHalfLoGain[i] = fArrBinningResLoGain[i] + fBinningResolutionHalfLoGain;
235 }
236}
237
238// --------------------------------------------------------------------------
239//
240// Apply the digital filter algorithm to the high-gain slices.
241//
242void MExtractTimeAndChargeDigitalFilter::FindTimeAndChargeHiGain(Byte_t *ptr, Byte_t *logain, Float_t &sum, Float_t &dsum,
243 Float_t &time, Float_t &dtime,
244 Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
245{
246
247 Int_t range = fHiGainLast - fHiGainFirst + 1;
248 Int_t maxpos = 0;
249
250 const Byte_t *end = ptr + range;
251 Byte_t *p = ptr;
252 Byte_t max = 0;
253
254 //
255 // Preparations for the pedestal subtraction (with AB-noise correction)
256 //
257 const Float_t pedes = ped.GetPedestal();
258 const Float_t ABoffs = ped.GetPedestalABoffset();
259
260 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
261
262 //
263 // Check for saturation in all other slices
264 //
265 Int_t ids = fHiGainFirst;
266 Float_t *sample = fHiGainSignal.GetArray();
267 while (p<end)
268 {
269 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
270
271 if (*p > max)
272 {
273 max = *p;
274 maxpos = p-ptr;
275 }
276
277 if (*p++ >= fSaturationLimit)
278 if (!sat)
279 sat = ids-4;
280 }
281
282 if (fHiLoLast != 0)
283 {
284
285 end = logain + fHiLoLast;
286
287 while (logain<end)
288 {
289
290 *sample++ = (Float_t)*logain - pedmean[(ids++ + abflag) & 0x1];
291
292 if (*logain > max)
293 {
294 max = *logain;
295 maxpos = range;
296 }
297
298 if (*logain++ >= fSaturationLimit)
299 if (!sat)
300 sat = ids-4;
301
302 range++;
303 }
304 }
305
306 //
307 // allow no saturated slice
308 //
309 if (sat > 0)
310 return;
311
312 //
313 // Slide with a window of size fWindowSizeHiGain over the sample
314 // and multiply the entries with the corresponding weights
315 //
316 if (IsNoiseCalculation())
317 {
318 if (fRandomIter == fBinningResolutionHiGain)
319 fRandomIter = 0;
320 for (Int_t ids=0; ids < fWindowSizeHiGain; ids++)
321 {
322 const Int_t idx = fArrBinningResHiGain[ids] + fRandomIter;
323 sum += fAmpWeightsHiGain [idx]*fHiGainSignal[ids];
324 }
325 fRandomIter++;
326 return;
327 }
328
329 dtime = 1.0;
330 const Int_t uplim = range-fWindowSizeHiGain+1;
331
332 if (maxpos > uplim)
333 {
334 time = uplim+1;
335 return;
336 }
337
338 maxpos -= 3;
339
340 if (maxpos < 0)
341 {
342 time = 0.;
343 return;
344 }
345
346 Float_t time_sum = 0.;
347 Float_t fmax = 0.;
348 Float_t ftime_max = 0.;
349 Int_t max_p = 0;
350
351 //
352 // Perform first only 3 temptative calculations around the maximum!
353 //
354 for (Int_t i=maxpos; i<maxpos+fWindowSizeHiGain+3;i++)
355 {
356 sum = 0.;
357 time_sum = 0.;
358
359 //
360 // Slide with a window of size fWindowSizeHiGain over the sample
361 // and multiply the entries with the corresponding weights
362 //
363 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
364 {
365 const Int_t idx = fArrBinningResHalfHiGain[sample];
366 const Float_t pex = fHiGainSignal[sample+i];
367 sum += fAmpWeightsHiGain [idx]*pex;
368 time_sum += fTimeWeightsHiGain[idx]*pex;
369 }
370
371 if (sum>fmax)
372 {
373 fmax = sum;
374 ftime_max = time_sum;
375 max_p = i;
376 }
377 } /* for (Int_t i=maxpos;i<maxpos+fWindowSizeHiGain+3;i++) */
378
379 time = 0;
380 if (fmax==0)
381 return;
382
383 ftime_max /= fmax;
384 Int_t t_iter = Int_t(ftime_max*fBinningResolutionHiGain);
385 Int_t sample_iter = 0;
386
387 while ( t_iter > fBinningResolutionHalfHiGain-1 || t_iter < -fBinningResolutionHalfHiGain )
388 {
389 if (t_iter > fBinningResolutionHalfHiGain-1)
390 {
391 t_iter -= fBinningResolutionHiGain;
392 max_p--;
393 sample_iter--;
394 }
395 if (t_iter < -fBinningResolutionHalfHiGain)
396 {
397 t_iter += fBinningResolutionHiGain;
398 max_p++;
399 sample_iter++;
400 }
401 }
402
403 sum = 0.;
404 time_sum = 0.;
405 //
406 // Slide with a window of size fWindowSizeHiGain over the sample
407 // and multiply the entries with the corresponding weights
408 //
409 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
410 {
411 const Int_t idx = fArrBinningResHalfHiGain[sample] + t_iter;
412 const Int_t ids = max_p + sample;
413 const Float_t pex = ids < 0 ? 0. : ( ids >= range ? 0. : fHiGainSignal[ids]);
414 sum += fAmpWeightsHiGain [idx]*pex;
415 time_sum += fTimeWeightsHiGain[idx]*pex;
416 }
417
418 if (sum == 0)
419 return;
420
421 time = max_p + fTimeShiftHiGain /* this shifts the time to the start of the rising edge */
422 - ((Float_t)t_iter)/fBinningResolutionHiGain - time_sum/sum;
423}
424
425// --------------------------------------------------------------------------
426//
427// Apply the digital filter algorithm to the low-gain slices.
428//
429void MExtractTimeAndChargeDigitalFilter::FindTimeAndChargeLoGain(Byte_t *ptr, Float_t &sum, Float_t &dsum,
430 Float_t &time, Float_t &dtime,
431 Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
432{
433
434 const Int_t range = fLoGainLast - fLoGainFirst + 1;
435
436 const Byte_t *end = ptr + range;
437 Byte_t *p = ptr;
438 //
439 // Prepare the low-gain pedestal
440 //
441 const Float_t pedes = ped.GetPedestal();
442 const Float_t ABoffs = ped.GetPedestalABoffset();
443
444 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
445
446 //
447 // Check for saturation in all other slices
448 //
449 Float_t *sample = fLoGainSignal.GetArray();
450 Int_t ids = fLoGainFirst;
451 while (p<end)
452 {
453 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
454
455 if (*p++ >= fSaturationLimit)
456 sat++;
457 }
458
459 //
460 // Slide with a window of size fWindowSizeLoGain over the sample
461 // and multiply the entries with the corresponding weights
462 //
463 if (IsNoiseCalculation())
464 {
465 if (fRandomIter == fBinningResolutionLoGain)
466 fRandomIter = 0;
467 for (Int_t ids=0; ids < fWindowSizeLoGain; ids++)
468 {
469 const Int_t idx = fArrBinningResLoGain[ids] + fRandomIter;
470 sum += fAmpWeightsLoGain [idx]*fLoGainSignal[ids];
471 }
472 return;
473 }
474
475 Float_t time_sum = 0.;
476 Float_t fmax = 0.;
477 Float_t ftime_max = 0.;
478 Int_t max_p = 0;
479
480 //
481 // Calculate the sum of the first fWindowSize slices
482 //
483 for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++)
484 {
485 sum = 0.;
486 time_sum = 0.;
487
488 //
489 // Slide with a window of size fWindowSizeLoGain over the sample
490 // and multiply the entries with the corresponding weights
491 //
492 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
493 {
494 const Int_t idx = fArrBinningResHalfLoGain[sample];
495 const Float_t pex = fLoGainSignal[sample+i];
496 sum += fAmpWeightsLoGain [idx]*pex;
497 time_sum += fTimeWeightsLoGain[idx]*pex;
498 }
499
500 if (sum>fmax)
501 {
502 fmax = sum;
503 ftime_max = time_sum;
504 max_p = i;
505 }
506 } /* for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++) */
507
508 time = 0;
509 if (fmax==0)
510 return;
511
512 ftime_max /= fmax;
513 Int_t t_iter = Int_t(ftime_max*fBinningResolutionLoGain);
514 Int_t sample_iter = 0;
515
516 while ( t_iter > fBinningResolutionHalfLoGain-1 || t_iter < -fBinningResolutionHalfLoGain )
517 {
518 if (t_iter > fBinningResolutionHalfLoGain-1)
519 {
520 t_iter -= fBinningResolutionLoGain;
521 max_p--;
522 sample_iter--;
523 }
524 if (t_iter < -fBinningResolutionHalfLoGain)
525 {
526 t_iter += fBinningResolutionLoGain;
527 max_p++;
528 sample_iter++;
529 }
530 }
531
532 sum = 0.;
533 time_sum = 0.;
534
535 //
536 // Slide with a window of size fWindowSizeLoGain over the sample
537 // and multiply the entries with the corresponding weights
538 //
539 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
540 {
541 const Int_t idx = fArrBinningResHalfLoGain[sample] + t_iter;
542 const Int_t ids = max_p + sample;
543 const Float_t pex = ids < 0 ? 0. : ( ids >= range ? 0. : fLoGainSignal[ids]);
544 sum += fAmpWeightsLoGain [idx]*pex;
545 time_sum += fTimeWeightsLoGain[idx]*pex;
546 }
547
548 if (sum == 0)
549 return;
550
551 time = max_p + fTimeShiftLoGain + (Float_t)fLoGainFirst /* this shifts the time to the start of the rising edge */
552 - ((Float_t)t_iter)/fBinningResolutionLoGain - time_sum/sum;
553}
554
555// --------------------------------------------------------------------------
556//
557// Read the setup from a TEnv, eg:
558// MJPedestal.MExtractor.WindowSizeHiGain: 6
559// MJPedestal.MExtractor.WindowSizeLoGain: 6
560// MJPedestal.MExtractor.BinningResolutionHiGain: 10
561// MJPedestal.MExtractor.BinningResolutionLoGain: 10
562// MJPedestal.MExtractor.WeightsFile: filename
563//
564Int_t MExtractTimeAndChargeDigitalFilter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
565{
566
567 Byte_t hw = fWindowSizeHiGain;
568 Byte_t lw = fWindowSizeLoGain;
569 Bool_t rc = kFALSE;
570
571 if (IsEnvDefined(env, prefix, "WindowSizeHiGain", print))
572 {
573 hw = GetEnvValue(env, prefix, "WindowSizeHiGain", hw);
574 rc = kTRUE;
575 }
576 if (IsEnvDefined(env, prefix, "WindowSizeLoGain", print))
577 {
578 lw = GetEnvValue(env, prefix, "WindowSizeLoGain", lw);
579 rc = kTRUE;
580 }
581
582 if (rc)
583 SetWindowSize(hw, lw);
584
585 Bool_t rc2 = kFALSE;
586 Int_t brh = fBinningResolutionHiGain;
587 Int_t brl = fBinningResolutionLoGain;
588
589 if (IsEnvDefined(env, prefix, "BinningResolutionHiGain", print))
590 {
591 brh = GetEnvValue(env, prefix, brh);
592 rc2 = kTRUE;
593 }
594 if (IsEnvDefined(env, prefix, "BinningResolutionLoGain", print))
595 {
596 brl = GetEnvValue(env, prefix, brl);
597 rc2 = kTRUE;
598 }
599
600 if (rc2)
601 {
602 SetBinningResolution(brh, brl);
603 rc = kTRUE;
604 }
605
606 if (IsEnvDefined(env, prefix, "WeightsFile", print))
607 {
608 if (!ReadWeightsFile(GetEnvValue(env, prefix, "WeightsFile", "")))
609 return kERROR;
610 rc = kTRUE;
611 }
612
613 return MExtractTimeAndCharge::ReadEnv(env, prefix, print) ? kTRUE : rc;
614}
615
616//----------------------------------------------------------------------------
617//
618// Read a pre-defined weights file into the class.
619// This is mandatory for the extraction
620//
621// If filenname is empty, then all weights will be set to 1.
622//
623Bool_t MExtractTimeAndChargeDigitalFilter::ReadWeightsFile(TString filename)
624{
625
626 // This is a fix for TEnv files edited with windows editors
627 filename.ReplaceAll("\015", "");
628
629 SetNameWeightsFile(filename);
630
631 fAmpWeightsHiGain .Set(fBinningResolutionHiGain*fWindowSizeHiGain);
632 fAmpWeightsLoGain .Set(fBinningResolutionLoGain*fWindowSizeLoGain);
633 fTimeWeightsHiGain.Set(fBinningResolutionHiGain*fWindowSizeHiGain);
634 fTimeWeightsLoGain.Set(fBinningResolutionLoGain*fWindowSizeLoGain);
635
636 if (fNameWeightsFile.IsNull())
637 {
638 fAmpWeightsHiGain.Reset(1);
639 fTimeWeightsHiGain.Reset(1);
640 fAmpWeightsLoGain.Reset(1);
641 fTimeWeightsLoGain.Reset(1);
642 return kTRUE;
643 }
644
645 ifstream fin(filename.Data());
646 if (!fin)
647 {
648 *fLog << err << GetDescriptor() << ": ERROR - Cannot open file " << filename << ": ";
649 *fLog << strerror(errno) << endl;
650 return kFALSE;
651 }
652
653 *fLog << inf << "Reading weights file " << filename << "..." << flush;
654
655 Int_t len = 0;
656 Int_t cnt = 0;
657 Int_t line = 0;
658 Bool_t hi = kFALSE;
659 Bool_t lo = kFALSE;
660
661 TString str;
662
663 while (1)
664 {
665 str.ReadLine(fin);
666 if (!fin)
667 break;
668
669 line++;
670
671 if (str.Contains("# High Gain Weights:"))
672 {
673 if (hi)
674 {
675 *fLog << err << "ERROR - 'High Gain Weights' found twice in line #" << line << "." << endl;
676 return kFALSE;
677 }
678
679 if (2!=sscanf(str.Data(), "# High Gain Weights:%2i %2i", &fWindowSizeHiGain, &fBinningResolutionHiGain))
680 {
681 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
682 *fLog << str << endl;
683 return kFALSE;
684 }
685
686 len = fBinningResolutionHiGain*fWindowSizeHiGain;
687 fAmpWeightsHiGain .Set(len);
688 fTimeWeightsHiGain.Set(len);
689 hi = kTRUE;
690 continue;
691 }
692
693 if (str.Contains("# Low Gain Weights:"))
694 {
695 if (lo)
696 {
697 *fLog << err << "ERROR - 'Lo Gain Weights' found twice in line #" << line << "." << endl;
698 return kFALSE;
699 }
700
701 if (2!=sscanf(str.Data(),"# Low Gain Weights:%2i %2i", &fWindowSizeLoGain, &fBinningResolutionLoGain))
702 {
703 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
704 *fLog << str << endl;
705 return kFALSE;
706 }
707
708 len = fBinningResolutionLoGain*fWindowSizeLoGain;
709 fAmpWeightsLoGain .Set(len);
710 fTimeWeightsLoGain.Set(len);
711 lo = kTRUE;
712 continue;
713 }
714
715 // Handle lines with comments
716 if (str.Contains("#"))
717 continue;
718
719 // Nothing found so far
720 if (len == 0)
721 continue;
722
723 if (2!=sscanf(str.Data(), "%f %f",
724 lo ? &fAmpWeightsLoGain [cnt] : &fAmpWeightsHiGain [cnt],
725 lo ? &fTimeWeightsLoGain[cnt] : &fTimeWeightsHiGain[cnt]))
726 {
727 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
728 *fLog << str << endl;
729 return kFALSE;
730 }
731
732 if (++cnt == len)
733 {
734 len = 0;
735 cnt = 0;
736 }
737 }
738
739 if (cnt != len)
740 {
741 *fLog << err << "Size mismatch in weights file " << filename << endl;
742 return kFALSE;
743 }
744
745 if (!hi)
746 {
747 *fLog << err << "No correct header found in weights file " << filename << endl;
748 return kFALSE;
749 }
750
751 *fLog << "done." << endl;
752
753 *fLog << inf << " File contains " << fWindowSizeHiGain << " hi-gain slices ";
754 *fLog << "with a resolution of " << fBinningResolutionHiGain << endl;
755
756 *fLog << inf << " File contains " << fWindowSizeLoGain << " lo-gain slices ";
757 *fLog << "with a resolution of " << fBinningResolutionLoGain << endl;
758
759 CalcBinningResArrays();
760
761 fWeightsSet = kTRUE;
762
763 return kTRUE;
764}
765
766//----------------------------------------------------------------------------
767//
768// Create the weights file
769// Beware that the shape-histogram has to contain the pulse starting at bin 1
770//
771Bool_t MExtractTimeAndChargeDigitalFilter::WriteWeightsFile(TString filename, TH1F *shapehi, TH2F *autocorrhi,
772 TH1F *shapelo, TH2F *autocorrlo )
773{
774
775 const Int_t nbinshi = shapehi->GetNbinsX();
776 Float_t binwidth = shapehi->GetBinWidth(1);
777
778 TH1F *derivativehi = new TH1F(Form("%s%s",shapehi->GetName(),"_der"),
779 Form("%s%s",shapehi->GetTitle()," derivative"),
780 nbinshi,
781 shapehi->GetBinLowEdge(1),
782 shapehi->GetBinLowEdge(nbinshi)+binwidth);
783
784 //
785 // Calculate the derivative of shapehi
786 //
787 for (Int_t i = 1; i<nbinshi+1;i++)
788 {
789 derivativehi->SetBinContent(i,
790 ((shapehi->GetBinContent(i+1)-shapehi->GetBinContent(i-1))/2./binwidth));
791 derivativehi->SetBinError(i,
792 (sqrt(shapehi->GetBinError(i+1)*shapehi->GetBinError(i+1)
793 +shapehi->GetBinError(i-1)*shapehi->GetBinError(i-1))/2./binwidth));
794 }
795
796 //
797 // normalize the shapehi, such that the integral for fWindowSize slices is one!
798 //
799 Float_t sum = 0;
800 Int_t lasttemp = fBinningResolutionHiGain * (fSignalStartBinHiGain + fWindowSizeHiGain);
801 lasttemp = lasttemp > nbinshi ? nbinshi : lasttemp;
802
803 for (Int_t i=fBinningResolutionHiGain*fSignalStartBinHiGain; i<lasttemp; i++) {
804 sum += shapehi->GetBinContent(i);
805 }
806 sum /= fBinningResolutionHiGain;
807
808 shapehi->Scale(1./sum);
809 derivativehi->Scale(1./sum);
810
811 //
812 // read in the noise auto-correlation function:
813 //
814 TMatrix Bhi(fWindowSizeHiGain,fWindowSizeHiGain);
815
816 for (Int_t i=0; i<fWindowSizeHiGain; i++){
817 for (Int_t j=0; j<fWindowSizeHiGain; j++){
818 Bhi[i][j]=autocorrhi->GetBinContent(i+1,j+1); //+fSignalStartBinHiGain +fSignalStartBinHiGain
819 }
820 }
821 Bhi.Invert();
822
823 const Int_t nsizehi = fWindowSizeHiGain*fBinningResolutionHiGain;
824 fAmpWeightsHiGain.Set(nsizehi);
825 fTimeWeightsHiGain.Set(nsizehi);
826
827 //
828 // Loop over relative time in one BinningResolution interval
829 //
830 Int_t start = fBinningResolutionHiGain*(fSignalStartBinHiGain + 1);
831
832 for (Int_t i = -fBinningResolutionHalfHiGain+1; i<=fBinningResolutionHalfHiGain; i++)
833 {
834
835 TMatrix g(fWindowSizeHiGain,1);
836 TMatrix gT(1,fWindowSizeHiGain);
837 TMatrix d(fWindowSizeHiGain,1);
838 TMatrix dT(1,fWindowSizeHiGain);
839
840 for (Int_t count=0; count < fWindowSizeHiGain; count++){
841
842 g[count][0]=shapehi->GetBinContent(start
843 +fBinningResolutionHiGain*count+i);
844 gT[0][count]=shapehi->GetBinContent(start
845 +fBinningResolutionHiGain*count+i);
846 d[count][0]=derivativehi->GetBinContent(start
847 +fBinningResolutionHiGain*count+i);
848 dT[0][count]=derivativehi->GetBinContent(start
849 +fBinningResolutionHiGain*count+i);
850 }
851
852 TMatrix m_denom = (gT*(Bhi*g))*(dT*(Bhi*d)) - (dT*(Bhi*g))*(dT*(Bhi*g));
853 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
854
855 TMatrix m_first = dT*(Bhi*d); // ROOT thinks, m_first is still a matrix
856 Float_t first = m_first[0][0]/denom;
857
858 TMatrix m_last = gT*(Bhi*d); // ROOT thinks, m_last is still a matrix
859 Float_t last = m_last[0][0]/denom;
860
861 TMatrix m1 = gT*Bhi;
862 m1 *= first;
863
864 TMatrix m2 = dT*Bhi;
865 m2 *=last;
866
867 TMatrix w_amp = m1 - m2;
868
869 TMatrix m_first1 = gT*(Bhi*g);
870 Float_t first1 = m_first1[0][0]/denom;
871
872 TMatrix m_last1 = gT*(Bhi*d);
873 Float_t last1 = m_last1 [0][0]/denom;
874
875 TMatrix m11 = dT*Bhi;
876 m11 *=first1;
877
878 TMatrix m21 = gT*Bhi;
879 m21 *=last1;
880
881 TMatrix w_time= m11 - m21;
882
883 for (Int_t count=0; count < fWindowSizeHiGain; count++)
884 {
885 const Int_t idx = i+fBinningResolutionHalfHiGain+fBinningResolutionHiGain*count-1;
886 fAmpWeightsHiGain [idx] = w_amp [0][count];
887 fTimeWeightsHiGain[idx] = w_time[0][count];
888 }
889
890 } // end loop over i
891
892 //
893 // Low Gain histograms
894 //
895 TH1F *derivativelo = NULL;
896 if (shapelo)
897 {
898 const Int_t nbinslo = shapelo->GetNbinsX();
899 binwidth = shapelo->GetBinWidth(1);
900
901 derivativelo = new TH1F(Form("%s%s",shapelo->GetName(),"_der"),
902 Form("%s%s",shapelo->GetTitle()," derivative"),
903 nbinslo,
904 shapelo->GetBinLowEdge(1),
905 shapelo->GetBinLowEdge(nbinslo)+binwidth);
906
907 //
908 // Calculate the derivative of shapelo
909 //
910 for (Int_t i = 1; i<nbinslo+1;i++)
911 {
912 derivativelo->SetBinContent(i,
913 ((shapelo->GetBinContent(i+1)-shapelo->GetBinContent(i-1))/2./binwidth));
914 derivativelo->SetBinError(i,
915 (sqrt(shapelo->GetBinError(i+1)*shapelo->GetBinError(i+1)
916 +shapelo->GetBinError(i-1)*shapelo->GetBinError(i-1))/2./binwidth));
917 }
918
919 //
920 // normalize the shapelo, such that the integral for fWindowSize slices is one!
921 //
922 sum = 0;
923 lasttemp = fBinningResolutionLoGain * (fSignalStartBinLoGain + fWindowSizeLoGain);
924 lasttemp = lasttemp > nbinslo ? nbinslo : lasttemp;
925
926 for (Int_t i=fBinningResolutionLoGain*fSignalStartBinLoGain; i<lasttemp; i++)
927 sum += shapelo->GetBinContent(i);
928
929 sum /= fBinningResolutionLoGain;
930
931 shapelo->Scale(1./sum);
932 derivativelo->Scale(1./sum);
933
934 //
935 // read in the noise auto-correlation function:
936 //
937 TMatrix Blo(fWindowSizeLoGain,fWindowSizeLoGain);
938
939 for (Int_t i=0; i<fWindowSizeLoGain; i++){
940 for (Int_t j=0; j<fWindowSizeLoGain; j++){
941 Blo[i][j]=autocorrlo->GetBinContent(i+1+fSignalStartBinLoGain,j+1+fSignalStartBinLoGain);
942 }
943 }
944 Blo.Invert();
945
946 const Int_t nsizelo = fWindowSizeLoGain*fBinningResolutionLoGain;
947 fAmpWeightsLoGain.Set(nsizelo);
948 fTimeWeightsLoGain.Set(nsizelo);
949
950 //
951 // Loop over relative time in one BinningResolution interval
952 //
953 Int_t start = fBinningResolutionLoGain*fSignalStartBinLoGain + fBinningResolutionHalfLoGain;
954
955 for (Int_t i = -fBinningResolutionHalfLoGain+1; i<=fBinningResolutionHalfLoGain; i++)
956 {
957
958 TMatrix g(fWindowSizeLoGain,1);
959 TMatrix gT(1,fWindowSizeLoGain);
960 TMatrix d(fWindowSizeLoGain,1);
961 TMatrix dT(1,fWindowSizeLoGain);
962
963 for (Int_t count=0; count < fWindowSizeLoGain; count++){
964
965 g[count][0] = shapelo->GetBinContent(start
966 +fBinningResolutionLoGain*count+i);
967 gT[0][count]= shapelo->GetBinContent(start
968 +fBinningResolutionLoGain*count+i);
969 d[count][0] = derivativelo->GetBinContent(start
970 +fBinningResolutionLoGain*count+i);
971 dT[0][count]= derivativelo->GetBinContent(start
972 +fBinningResolutionLoGain*count+i);
973 }
974
975 TMatrix m_denom = (gT*(Blo*g))*(dT*(Blo*d)) - (dT*(Blo*g))*(dT*(Blo*g));
976 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
977
978 TMatrix m_first = dT*(Blo*d); // ROOT thinks, m_first is still a matrix
979 Float_t first = m_first[0][0]/denom;
980
981 TMatrix m_last = gT*(Blo*d); // ROOT thinks, m_last is still a matrix
982 Float_t last = m_last[0][0]/denom;
983
984 TMatrix m1 = gT*Blo;
985 m1 *= first;
986
987 TMatrix m2 = dT*Blo;
988 m2 *=last;
989
990 TMatrix w_amp = m1 - m2;
991
992 TMatrix m_first1 = gT*(Blo*g);
993 Float_t first1 = m_first1[0][0]/denom;
994
995 TMatrix m_last1 = gT*(Blo*d);
996 Float_t last1 = m_last1 [0][0]/denom;
997
998 TMatrix m11 = dT*Blo;
999 m11 *=first1;
1000
1001 TMatrix m21 = gT*Blo;
1002 m21 *=last1;
1003
1004 TMatrix w_time= m11 - m21;
1005
1006 for (Int_t count=0; count < fWindowSizeLoGain; count++)
1007 {
1008 const Int_t idx = i+fBinningResolutionHalfLoGain+fBinningResolutionLoGain*count-1;
1009 fAmpWeightsLoGain [idx] = w_amp [0][count];
1010 fTimeWeightsLoGain[idx] = w_time[0][count];
1011 }
1012
1013 } // end loop over i
1014 }
1015
1016 ofstream fn(filename.Data());
1017
1018 fn << "# High Gain Weights: " << fWindowSizeHiGain << " " << fBinningResolutionHiGain << endl;
1019 fn << "# (Amplitude) (Time) " << endl;
1020
1021 for (Int_t i=0; i<nsizehi; i++)
1022 fn << "\t" << fAmpWeightsHiGain[i] << "\t" << fTimeWeightsHiGain[i] << endl;
1023
1024 fn << "# Low Gain Weights: " << fWindowSizeLoGain << " " << fBinningResolutionLoGain << endl;
1025 fn << "# (Amplitude) (Time) " << endl;
1026
1027 for (Int_t i=0; i<nsizehi; i++)
1028 fn << "\t" << fAmpWeightsLoGain[i] << "\t" << fTimeWeightsLoGain[i] << endl;
1029
1030 delete derivativehi;
1031 if (derivativelo)
1032 delete derivativelo;
1033
1034 return kTRUE;
1035}
1036
1037void MExtractTimeAndChargeDigitalFilter::Print(Option_t *o) const
1038{
1039 if (IsA()==Class())
1040 *fLog << GetDescriptor() << ":" << endl;
1041
1042 MExtractTimeAndCharge::Print(o);
1043 *fLog << " Time Shift HiGain: " << fTimeShiftHiGain << " LoGain: " << fTimeShiftLoGain << endl;
1044 *fLog << " Window Size HiGain: " << fWindowSizeHiGain << " LoGain: " << fWindowSizeLoGain << endl;
1045 *fLog << " Binning Res HiGain: " << fBinningResolutionHiGain << " LoGain: " << fBinningResolutionHiGain << endl;
1046 *fLog << " Weights File: " << fNameWeightsFile << endl;
1047
1048 TString opt(o);
1049 if (!opt.Contains("weights"))
1050 return;
1051
1052 *fLog << endl;
1053 *fLog << inf << "Using the following weights: " << endl;
1054 *fLog << "Hi-Gain:" << endl;
1055 for (Int_t i=0; i<fBinningResolutionHiGain*fWindowSizeHiGain; i++)
1056 *fLog << " " << fAmpWeightsHiGain[i] << " \t " << fTimeWeightsHiGain[i] << endl;
1057
1058 *fLog << "Lo-Gain:" << endl;
1059 for (Int_t i=0; i<fBinningResolutionLoGain*fWindowSizeLoGain; i++)
1060 *fLog << " " << fAmpWeightsLoGain[i] << " \t " << fTimeWeightsLoGain[i] << endl;
1061}
Note: See TracBrowser for help on using the repository browser.