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

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