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

Last change on this file since 5980 was 5980, checked in by gaug, 20 years ago
*** empty log message ***
File size: 32.8 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
249 const Byte_t *end = ptr + range;
250 Byte_t *p = ptr;
251
252 //
253 // Preparations for the pedestal subtraction (with AB-noise correction)
254 //
255 const Float_t pedes = ped.GetPedestal();
256 const Float_t ABoffs = ped.GetPedestalABoffset();
257
258 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
259
260 fMaxBinContent = 0;
261 //
262 // Check for saturation in all other slices
263 //
264 Int_t ids = fHiGainFirst;
265 Float_t *sample = fHiGainSignal.GetArray();
266 while (p<end)
267 {
268 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
269
270 if (*p > fMaxBinContent)
271 fMaxBinContent = *p;
272
273 if (*p++ >= fSaturationLimit)
274 if (!sat)
275 sat = ids-4;
276 }
277
278 if (fHiLoLast != 0)
279 {
280
281 end = logain + fHiLoLast;
282
283 while (logain<end)
284 {
285
286 *sample++ = (Float_t)*logain - pedmean[(ids++ + abflag) & 0x1];
287
288 if (*logain > fMaxBinContent)
289 fMaxBinContent = *logain;
290
291 if (*logain++ >= fSaturationLimit)
292 if (!sat)
293 sat = ids-4;
294
295 range++;
296 }
297 }
298
299 //
300 // allow no saturated slice
301 //
302 if (sat > 0)
303 return;
304
305 //
306 // Slide with a window of size fWindowSizeHiGain over the sample
307 // and multiply the entries with the corresponding weights
308 //
309 if (IsNoiseCalculation())
310 {
311 if (fRandomIter == fBinningResolutionHiGain)
312 fRandomIter = 0;
313 for (Int_t ids=0; ids < fWindowSizeHiGain; ids++)
314 {
315 const Int_t idx = fArrBinningResHiGain[ids] + fRandomIter;
316 sum += fAmpWeightsHiGain [idx]*fHiGainSignal[ids];
317 }
318 fRandomIter++;
319 return;
320 }
321
322 Float_t time_sum = 0.;
323 Float_t fmax = 0.;
324 Float_t ftime_max = 0.;
325 Int_t max_p = 0;
326
327 //
328 // Calculate the sum of the first fWindowSize slices
329 //
330 for (Int_t i=0;i<range-fWindowSizeHiGain+1;i++)
331 {
332 sum = 0.;
333 time_sum = 0.;
334
335 //
336 // Slide with a window of size fWindowSizeHiGain over the sample
337 // and multiply the entries with the corresponding weights
338 //
339 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
340 {
341 const Int_t idx = fBinningResolutionHiGain*sample+fBinningResolutionHalfHiGain;
342 const Float_t pex = fHiGainSignal[sample+i];
343 sum += fAmpWeightsHiGain [idx]*pex;
344 time_sum += fTimeWeightsHiGain[idx]*pex;
345 }
346
347 if (sum>fmax)
348 {
349 fmax = sum;
350 ftime_max = time_sum;
351 max_p = i;
352 }
353 } /* for (Int_t i=0;i<range-fWindowSizeHiGain+1;i++) */
354
355 if (fmax==0)
356 return;
357
358 ftime_max /= fmax;
359 Int_t t_iter = Int_t(ftime_max*fBinningResolutionHiGain);
360 Int_t sample_iter = 0;
361
362 while ( t_iter > fBinningResolutionHalfHiGain-1 || t_iter < -fBinningResolutionHalfHiGain )
363 {
364 if (t_iter > fBinningResolutionHalfHiGain-1)
365 {
366 t_iter -= fBinningResolutionHiGain;
367 max_p--;
368 sample_iter--;
369 }
370 if (t_iter < -fBinningResolutionHalfHiGain)
371 {
372 t_iter += fBinningResolutionHiGain;
373 max_p++;
374 sample_iter++;
375 }
376 }
377
378 sum = 0.;
379 time_sum = 0.;
380 //
381 // Slide with a window of size fWindowSizeHiGain over the sample
382 // and multiply the entries with the corresponding weights
383 //
384 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
385 {
386 const Int_t idx = fArrBinningResHalfHiGain[sample] + t_iter;
387 const Int_t ids = max_p + sample;
388 const Float_t pex = ids < 0 ? 0. : ( ids >= range ? 0. : fHiGainSignal[ids]);
389 sum += fAmpWeightsHiGain [idx]*pex;
390 time_sum += fTimeWeightsHiGain[idx]*pex;
391 }
392
393 if (sum == 0)
394 return;
395
396 time = max_p + fTimeShiftHiGain /* this shifts the time to the start of the rising edge */
397 - ((Float_t)t_iter)/fBinningResolutionHiGain;
398
399 const Float_t timefineadjust = time_sum/sum;
400
401 if (timefineadjust < 2./fBinningResolutionHiGain)
402 time -= timefineadjust;
403
404}
405
406// --------------------------------------------------------------------------
407//
408// Apply the digital filter algorithm to the low-gain slices.
409//
410void MExtractTimeAndChargeDigitalFilter::FindTimeAndChargeLoGain(Byte_t *ptr, Float_t &sum, Float_t &dsum,
411 Float_t &time, Float_t &dtime,
412 Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
413{
414
415 const Int_t range = fLoGainLast - fLoGainFirst + 1;
416
417 const Byte_t *end = ptr + range;
418 Byte_t *p = ptr;
419 //
420 // Prepare the low-gain pedestal
421 //
422 const Float_t pedes = ped.GetPedestal();
423 const Float_t ABoffs = ped.GetPedestalABoffset();
424
425 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
426
427 //
428 // Check for saturation in all other slices
429 //
430 Float_t *sample = fLoGainSignal.GetArray();
431 Int_t ids = fLoGainFirst;
432 while (p<end)
433 {
434 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
435
436 if (*p++ >= fSaturationLimit)
437 sat++;
438 }
439
440 //
441 // Slide with a window of size fWindowSizeLoGain over the sample
442 // and multiply the entries with the corresponding weights
443 //
444 if (IsNoiseCalculation())
445 {
446 if (fRandomIter == fBinningResolutionLoGain)
447 fRandomIter = 0;
448 for (Int_t ids=0; ids < fWindowSizeLoGain; ids++)
449 {
450 const Int_t idx = fArrBinningResLoGain[ids] + fRandomIter;
451 sum += fAmpWeightsLoGain [idx]*fLoGainSignal[ids];
452 }
453 return;
454 }
455
456 Float_t time_sum = 0.;
457 Float_t fmax = 0.;
458 Float_t ftime_max = 0.;
459 Int_t max_p = 0;
460
461 //
462 // Calculate the sum of the first fWindowSize slices
463 //
464 for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++)
465 {
466 sum = 0.;
467 time_sum = 0.;
468
469 //
470 // Slide with a window of size fWindowSizeLoGain over the sample
471 // and multiply the entries with the corresponding weights
472 //
473 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
474 {
475 const Int_t idx = fArrBinningResHalfLoGain[sample];
476 const Float_t pex = fLoGainSignal[sample+i];
477 sum += fAmpWeightsLoGain [idx]*pex;
478 time_sum += fTimeWeightsLoGain[idx]*pex;
479 }
480
481 if (sum>fmax)
482 {
483 fmax = sum;
484 ftime_max = time_sum;
485 max_p = i;
486 }
487 } /* for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++) */
488
489 time = 0;
490 if (fmax==0)
491 return;
492
493 ftime_max /= fmax;
494 Int_t t_iter = Int_t(ftime_max*fBinningResolutionLoGain);
495 Int_t sample_iter = 0;
496
497 while ( t_iter > fBinningResolutionHalfLoGain-1 || t_iter < -fBinningResolutionHalfLoGain )
498 {
499 if (t_iter > fBinningResolutionHalfLoGain-1)
500 {
501 t_iter -= fBinningResolutionLoGain;
502 max_p--;
503 sample_iter--;
504 }
505 if (t_iter < -fBinningResolutionHalfLoGain)
506 {
507 t_iter += fBinningResolutionLoGain;
508 max_p++;
509 sample_iter++;
510 }
511 }
512
513 sum = 0.;
514 time_sum = 0.;
515
516 //
517 // Slide with a window of size fWindowSizeLoGain over the sample
518 // and multiply the entries with the corresponding weights
519 //
520 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
521 {
522 const Int_t idx = fArrBinningResHalfLoGain[sample] + t_iter;
523 const Int_t ids = max_p + sample;
524 const Float_t pex = ids < 0 ? 0. : ( ids >= range ? 0. : fLoGainSignal[ids]);
525 sum += fAmpWeightsLoGain [idx]*pex;
526 time_sum += fTimeWeightsLoGain[idx]*pex;
527 }
528
529 if (sum == 0)
530 return;
531
532 time = max_p + fTimeShiftLoGain + (Float_t)fLoGainFirst /* this shifts the time to the start of the rising edge */
533 - ((Float_t)t_iter)/fBinningResolutionLoGain - time_sum/sum;
534}
535
536// --------------------------------------------------------------------------
537//
538// Read the setup from a TEnv, eg:
539// MJPedestal.MExtractor.WindowSizeHiGain: 6
540// MJPedestal.MExtractor.WindowSizeLoGain: 6
541// MJPedestal.MExtractor.BinningResolutionHiGain: 10
542// MJPedestal.MExtractor.BinningResolutionLoGain: 10
543// MJPedestal.MExtractor.WeightsFile: filename
544//
545Int_t MExtractTimeAndChargeDigitalFilter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
546{
547
548 Byte_t hw = fWindowSizeHiGain;
549 Byte_t lw = fWindowSizeLoGain;
550 Bool_t rc = kFALSE;
551
552 if (IsEnvDefined(env, prefix, "WindowSizeHiGain", print))
553 {
554 hw = GetEnvValue(env, prefix, "WindowSizeHiGain", hw);
555 rc = kTRUE;
556 }
557 if (IsEnvDefined(env, prefix, "WindowSizeLoGain", print))
558 {
559 lw = GetEnvValue(env, prefix, "WindowSizeLoGain", lw);
560 rc = kTRUE;
561 }
562
563 if (rc)
564 SetWindowSize(hw, lw);
565
566 Bool_t rc2 = kFALSE;
567 Int_t brh = fBinningResolutionHiGain;
568 Int_t brl = fBinningResolutionLoGain;
569
570 if (IsEnvDefined(env, prefix, "BinningResolutionHiGain", print))
571 {
572 brh = GetEnvValue(env, prefix, brh);
573 rc2 = kTRUE;
574 }
575 if (IsEnvDefined(env, prefix, "BinningResolutionLoGain", print))
576 {
577 brl = GetEnvValue(env, prefix, brl);
578 rc2 = kTRUE;
579 }
580
581 if (rc2)
582 {
583 SetBinningResolution(brh, brl);
584 rc = kTRUE;
585 }
586
587 if (IsEnvDefined(env, prefix, "WeightsFile", print))
588 {
589 if (!ReadWeightsFile(GetEnvValue(env, prefix, "WeightsFile", "")))
590 return kERROR;
591 rc = kTRUE;
592 }
593
594 return MExtractTimeAndCharge::ReadEnv(env, prefix, print) ? kTRUE : rc;
595}
596
597//----------------------------------------------------------------------------
598//
599// Read a pre-defined weights file into the class.
600// This is mandatory for the extraction
601//
602// If filenname is empty, then all weights will be set to 1.
603//
604Bool_t MExtractTimeAndChargeDigitalFilter::ReadWeightsFile(TString filename)
605{
606
607 // This is a fix for TEnv files edited with windows editors
608 filename.ReplaceAll("\015", "");
609
610 SetNameWeightsFile(filename);
611
612 fAmpWeightsHiGain .Set(fBinningResolutionHiGain*fWindowSizeHiGain);
613 fAmpWeightsLoGain .Set(fBinningResolutionLoGain*fWindowSizeLoGain);
614 fTimeWeightsHiGain.Set(fBinningResolutionHiGain*fWindowSizeHiGain);
615 fTimeWeightsLoGain.Set(fBinningResolutionLoGain*fWindowSizeLoGain);
616
617 if (fNameWeightsFile.IsNull())
618 {
619 fAmpWeightsHiGain.Reset(1);
620 fTimeWeightsHiGain.Reset(1);
621 fAmpWeightsLoGain.Reset(1);
622 fTimeWeightsLoGain.Reset(1);
623 return kTRUE;
624 }
625
626 ifstream fin(filename.Data());
627 if (!fin)
628 {
629 *fLog << err << GetDescriptor() << ": ERROR - Cannot open file " << filename << ": ";
630 *fLog << strerror(errno) << endl;
631 return kFALSE;
632 }
633
634 *fLog << inf << "Reading weights file " << filename << "..." << flush;
635
636 Int_t len = 0;
637 Int_t cnt = 0;
638 Int_t line = 0;
639 Bool_t hi = kFALSE;
640 Bool_t lo = kFALSE;
641
642 TString str;
643
644 while (1)
645 {
646 str.ReadLine(fin);
647 if (!fin)
648 break;
649
650 line++;
651
652 if (str.Contains("# High Gain Weights:"))
653 {
654 if (hi)
655 {
656 *fLog << err << "ERROR - 'High Gain Weights' found twice in line #" << line << "." << endl;
657 return kFALSE;
658 }
659
660 if (2!=sscanf(str.Data(), "# High Gain Weights:%2i %2i", &fWindowSizeHiGain, &fBinningResolutionHiGain))
661 {
662 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
663 *fLog << str << endl;
664 return kFALSE;
665 }
666
667 len = fBinningResolutionHiGain*fWindowSizeHiGain;
668 fAmpWeightsHiGain .Set(len);
669 fTimeWeightsHiGain.Set(len);
670 hi = kTRUE;
671 continue;
672 }
673
674 if (str.Contains("# Low Gain Weights:"))
675 {
676 if (lo)
677 {
678 *fLog << err << "ERROR - 'Lo Gain Weights' found twice in line #" << line << "." << endl;
679 return kFALSE;
680 }
681
682 if (2!=sscanf(str.Data(),"# Low Gain Weights:%2i %2i", &fWindowSizeLoGain, &fBinningResolutionLoGain))
683 {
684 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
685 *fLog << str << endl;
686 return kFALSE;
687 }
688
689 len = fBinningResolutionLoGain*fWindowSizeLoGain;
690 fAmpWeightsLoGain .Set(len);
691 fTimeWeightsLoGain.Set(len);
692 lo = kTRUE;
693 continue;
694 }
695
696 // Handle lines with comments
697 if (str.Contains("#"))
698 continue;
699
700 // Nothing found so far
701 if (len == 0)
702 continue;
703
704 if (2!=sscanf(str.Data(), "%f %f",
705 lo ? &fAmpWeightsLoGain [cnt] : &fAmpWeightsHiGain [cnt],
706 lo ? &fTimeWeightsLoGain[cnt] : &fTimeWeightsHiGain[cnt]))
707 {
708 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
709 *fLog << str << endl;
710 return kFALSE;
711 }
712
713 if (++cnt == len)
714 {
715 len = 0;
716 cnt = 0;
717 }
718 }
719
720 if (cnt != len)
721 {
722 *fLog << err << "Size mismatch in weights file " << filename << endl;
723 return kFALSE;
724 }
725
726 if (!hi)
727 {
728 *fLog << err << "No correct header found in weights file " << filename << endl;
729 return kFALSE;
730 }
731
732 *fLog << "done." << endl;
733
734 *fLog << inf << " File contains " << fWindowSizeHiGain << " hi-gain slices ";
735 *fLog << "with a resolution of " << fBinningResolutionHiGain << endl;
736
737 *fLog << inf << " File contains " << fWindowSizeLoGain << " lo-gain slices ";
738 *fLog << "with a resolution of " << fBinningResolutionLoGain << endl;
739
740 CalcBinningResArrays();
741
742 fWeightsSet = kTRUE;
743
744 return kTRUE;
745}
746
747//----------------------------------------------------------------------------
748//
749// Create the weights file
750// Beware that the shape-histogram has to contain the pulse starting at bin 1
751//
752Bool_t MExtractTimeAndChargeDigitalFilter::WriteWeightsFile(TString filename, TH1F *shapehi, TH2F *autocorrhi,
753 TH1F *shapelo, TH2F *autocorrlo )
754{
755
756 const Int_t nbinshi = shapehi->GetNbinsX();
757 Float_t binwidth = shapehi->GetBinWidth(1);
758
759 TH1F *derivativehi = new TH1F(Form("%s%s",shapehi->GetName(),"_der"),
760 Form("%s%s",shapehi->GetTitle()," derivative"),
761 nbinshi,
762 shapehi->GetBinLowEdge(1),
763 shapehi->GetBinLowEdge(nbinshi)+binwidth);
764
765 //
766 // Calculate the derivative of shapehi
767 //
768 for (Int_t i = 1; i<nbinshi+1;i++)
769 {
770 derivativehi->SetBinContent(i,
771 ((shapehi->GetBinContent(i+1)-shapehi->GetBinContent(i-1))/2./binwidth));
772 derivativehi->SetBinError(i,
773 (sqrt(shapehi->GetBinError(i+1)*shapehi->GetBinError(i+1)
774 +shapehi->GetBinError(i-1)*shapehi->GetBinError(i-1))/2./binwidth));
775 }
776
777 //
778 // normalize the shapehi, such that the integral for fWindowSize slices is one!
779 //
780 Float_t sum = 0;
781 Int_t lasttemp = fBinningResolutionHiGain * (fSignalStartBinHiGain + fWindowSizeHiGain);
782 lasttemp = lasttemp > nbinshi ? nbinshi : lasttemp;
783
784 for (Int_t i=fBinningResolutionHiGain*fSignalStartBinHiGain; i<lasttemp; i++) {
785 sum += shapehi->GetBinContent(i);
786 }
787 sum /= fBinningResolutionHiGain;
788
789 shapehi->Scale(1./sum);
790 derivativehi->Scale(1./sum);
791
792 //
793 // read in the noise auto-correlation function:
794 //
795 TMatrix Bhi(fWindowSizeHiGain,fWindowSizeHiGain);
796
797 for (Int_t i=0; i<fWindowSizeHiGain; i++){
798 for (Int_t j=0; j<fWindowSizeHiGain; j++){
799 Bhi[i][j]=autocorrhi->GetBinContent(i+1,j+1); //+fSignalStartBinHiGain +fSignalStartBinHiGain
800 }
801 }
802 Bhi.Invert();
803
804 const Int_t nsizehi = fWindowSizeHiGain*fBinningResolutionHiGain;
805 fAmpWeightsHiGain.Set(nsizehi);
806 fTimeWeightsHiGain.Set(nsizehi);
807
808 //
809 // Loop over relative time in one BinningResolution interval
810 //
811 Int_t start = fBinningResolutionHiGain*(fSignalStartBinHiGain + 1);
812
813 for (Int_t i = -fBinningResolutionHalfHiGain+1; i<=fBinningResolutionHalfHiGain; i++)
814 {
815
816 TMatrix g(fWindowSizeHiGain,1);
817 TMatrix gT(1,fWindowSizeHiGain);
818 TMatrix d(fWindowSizeHiGain,1);
819 TMatrix dT(1,fWindowSizeHiGain);
820
821 for (Int_t count=0; count < fWindowSizeHiGain; count++){
822
823 g[count][0]=shapehi->GetBinContent(start
824 +fBinningResolutionHiGain*count+i);
825 gT[0][count]=shapehi->GetBinContent(start
826 +fBinningResolutionHiGain*count+i);
827 d[count][0]=derivativehi->GetBinContent(start
828 +fBinningResolutionHiGain*count+i);
829 dT[0][count]=derivativehi->GetBinContent(start
830 +fBinningResolutionHiGain*count+i);
831 }
832
833 TMatrix m_denom = (gT*(Bhi*g))*(dT*(Bhi*d)) - (dT*(Bhi*g))*(dT*(Bhi*g));
834 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
835
836 TMatrix m_first = dT*(Bhi*d); // ROOT thinks, m_first is still a matrix
837 Float_t first = m_first[0][0]/denom;
838
839 TMatrix m_last = gT*(Bhi*d); // ROOT thinks, m_last is still a matrix
840 Float_t last = m_last[0][0]/denom;
841
842 TMatrix m1 = gT*Bhi;
843 m1 *= first;
844
845 TMatrix m2 = dT*Bhi;
846 m2 *=last;
847
848 TMatrix w_amp = m1 - m2;
849
850 TMatrix m_first1 = gT*(Bhi*g);
851 Float_t first1 = m_first1[0][0]/denom;
852
853 TMatrix m_last1 = gT*(Bhi*d);
854 Float_t last1 = m_last1 [0][0]/denom;
855
856 TMatrix m11 = dT*Bhi;
857 m11 *=first1;
858
859 TMatrix m21 = gT*Bhi;
860 m21 *=last1;
861
862 TMatrix w_time= m11 - m21;
863
864 for (Int_t count=0; count < fWindowSizeHiGain; count++)
865 {
866 const Int_t idx = i+fBinningResolutionHalfHiGain+fBinningResolutionHiGain*count-1;
867 fAmpWeightsHiGain [idx] = w_amp [0][count];
868 fTimeWeightsHiGain[idx] = w_time[0][count];
869 }
870
871 } // end loop over i
872
873 //
874 // Low Gain histograms
875 //
876 TH1F *derivativelo = NULL;
877 if (shapelo)
878 {
879 const Int_t nbinslo = shapelo->GetNbinsX();
880 binwidth = shapelo->GetBinWidth(1);
881
882 derivativelo = new TH1F(Form("%s%s",shapelo->GetName(),"_der"),
883 Form("%s%s",shapelo->GetTitle()," derivative"),
884 nbinslo,
885 shapelo->GetBinLowEdge(1),
886 shapelo->GetBinLowEdge(nbinslo)+binwidth);
887
888 //
889 // Calculate the derivative of shapelo
890 //
891 for (Int_t i = 1; i<nbinslo+1;i++)
892 {
893 derivativelo->SetBinContent(i,
894 ((shapelo->GetBinContent(i+1)-shapelo->GetBinContent(i-1))/2./binwidth));
895 derivativelo->SetBinError(i,
896 (sqrt(shapelo->GetBinError(i+1)*shapelo->GetBinError(i+1)
897 +shapelo->GetBinError(i-1)*shapelo->GetBinError(i-1))/2./binwidth));
898 }
899
900 //
901 // normalize the shapelo, such that the integral for fWindowSize slices is one!
902 //
903 sum = 0;
904 lasttemp = fBinningResolutionLoGain * (fSignalStartBinLoGain + fWindowSizeLoGain);
905 lasttemp = lasttemp > nbinslo ? nbinslo : lasttemp;
906
907 for (Int_t i=fBinningResolutionLoGain*fSignalStartBinLoGain; i<lasttemp; i++)
908 sum += shapelo->GetBinContent(i);
909
910 sum /= fBinningResolutionLoGain;
911
912 shapelo->Scale(1./sum);
913 derivativelo->Scale(1./sum);
914
915 //
916 // read in the noise auto-correlation function:
917 //
918 TMatrix Blo(fWindowSizeLoGain,fWindowSizeLoGain);
919
920 for (Int_t i=0; i<fWindowSizeLoGain; i++){
921 for (Int_t j=0; j<fWindowSizeLoGain; j++){
922 Blo[i][j]=autocorrlo->GetBinContent(i+1+fSignalStartBinLoGain,j+1+fSignalStartBinLoGain);
923 }
924 }
925 Blo.Invert();
926
927 const Int_t nsizelo = fWindowSizeLoGain*fBinningResolutionLoGain;
928 fAmpWeightsLoGain.Set(nsizelo);
929 fTimeWeightsLoGain.Set(nsizelo);
930
931 //
932 // Loop over relative time in one BinningResolution interval
933 //
934 Int_t start = fBinningResolutionLoGain*fSignalStartBinLoGain + fBinningResolutionHalfLoGain;
935
936 for (Int_t i = -fBinningResolutionHalfLoGain+1; i<=fBinningResolutionHalfLoGain; i++)
937 {
938
939 TMatrix g(fWindowSizeLoGain,1);
940 TMatrix gT(1,fWindowSizeLoGain);
941 TMatrix d(fWindowSizeLoGain,1);
942 TMatrix dT(1,fWindowSizeLoGain);
943
944 for (Int_t count=0; count < fWindowSizeLoGain; count++){
945
946 g[count][0] = shapelo->GetBinContent(start
947 +fBinningResolutionLoGain*count+i);
948 gT[0][count]= shapelo->GetBinContent(start
949 +fBinningResolutionLoGain*count+i);
950 d[count][0] = derivativelo->GetBinContent(start
951 +fBinningResolutionLoGain*count+i);
952 dT[0][count]= derivativelo->GetBinContent(start
953 +fBinningResolutionLoGain*count+i);
954 }
955
956 TMatrix m_denom = (gT*(Blo*g))*(dT*(Blo*d)) - (dT*(Blo*g))*(dT*(Blo*g));
957 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
958
959 TMatrix m_first = dT*(Blo*d); // ROOT thinks, m_first is still a matrix
960 Float_t first = m_first[0][0]/denom;
961
962 TMatrix m_last = gT*(Blo*d); // ROOT thinks, m_last is still a matrix
963 Float_t last = m_last[0][0]/denom;
964
965 TMatrix m1 = gT*Blo;
966 m1 *= first;
967
968 TMatrix m2 = dT*Blo;
969 m2 *=last;
970
971 TMatrix w_amp = m1 - m2;
972
973 TMatrix m_first1 = gT*(Blo*g);
974 Float_t first1 = m_first1[0][0]/denom;
975
976 TMatrix m_last1 = gT*(Blo*d);
977 Float_t last1 = m_last1 [0][0]/denom;
978
979 TMatrix m11 = dT*Blo;
980 m11 *=first1;
981
982 TMatrix m21 = gT*Blo;
983 m21 *=last1;
984
985 TMatrix w_time= m11 - m21;
986
987 for (Int_t count=0; count < fWindowSizeLoGain; count++)
988 {
989 const Int_t idx = i+fBinningResolutionHalfLoGain+fBinningResolutionLoGain*count-1;
990 fAmpWeightsLoGain [idx] = w_amp [0][count];
991 fTimeWeightsLoGain[idx] = w_time[0][count];
992 }
993
994 } // end loop over i
995 }
996
997 ofstream fn(filename.Data());
998
999 fn << "# High Gain Weights: " << fWindowSizeHiGain << " " << fBinningResolutionHiGain << endl;
1000 fn << "# (Amplitude) (Time) " << endl;
1001
1002 for (Int_t i=0; i<nsizehi; i++)
1003 fn << "\t" << fAmpWeightsHiGain[i] << "\t" << fTimeWeightsHiGain[i] << endl;
1004
1005 fn << "# Low Gain Weights: " << fWindowSizeLoGain << " " << fBinningResolutionLoGain << endl;
1006 fn << "# (Amplitude) (Time) " << endl;
1007
1008 for (Int_t i=0; i<nsizehi; i++)
1009 fn << "\t" << fAmpWeightsLoGain[i] << "\t" << fTimeWeightsLoGain[i] << endl;
1010
1011 delete derivativehi;
1012 if (derivativelo)
1013 delete derivativelo;
1014
1015 return kTRUE;
1016}
1017
1018void MExtractTimeAndChargeDigitalFilter::Print(Option_t *o) const
1019{
1020 if (IsA()==Class())
1021 *fLog << GetDescriptor() << ":" << endl;
1022
1023 MExtractTimeAndCharge::Print(o);
1024 *fLog << " Time Shift HiGain: " << fTimeShiftHiGain << " LoGain: " << fTimeShiftLoGain << endl;
1025 *fLog << " Window Size HiGain: " << fWindowSizeHiGain << " LoGain: " << fWindowSizeLoGain << endl;
1026 *fLog << " Binning Res HiGain: " << fBinningResolutionHiGain << " LoGain: " << fBinningResolutionHiGain << endl;
1027 *fLog << " Weights File: " << fNameWeightsFile << endl;
1028
1029 TString opt(o);
1030 if (!opt.Contains("weights"))
1031 return;
1032
1033 *fLog << endl;
1034 *fLog << inf << "Using the following weights: " << endl;
1035 *fLog << "Hi-Gain:" << endl;
1036 for (Int_t i=0; i<fBinningResolutionHiGain*fWindowSizeHiGain; i++)
1037 *fLog << " " << fAmpWeightsHiGain[i] << " \t " << fTimeWeightsHiGain[i] << endl;
1038
1039 *fLog << "Lo-Gain:" << endl;
1040 for (Int_t i=0; i<fBinningResolutionLoGain*fWindowSizeLoGain; i++)
1041 *fLog << " " << fAmpWeightsLoGain[i] << " \t " << fTimeWeightsLoGain[i] << endl;
1042}
Note: See TracBrowser for help on using the repository browser.