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

Last change on this file since 7523 was 7430, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 41.7 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// The following variables have to be set by the derived class and
33// do not have defaults:
34// - fNumHiGainSamples
35// - fNumLoGainSamples
36// - fSqrtHiGainSamples
37// - fSqrtLoGainSamples
38//
39// The reading of automatic weights files (color, type) can be switched
40// off using EnableAutomaticWeights(kFALSE).
41//
42// An empty name or "-" as the weights file name is a synonym for
43// setting all weights to 1
44//
45// Input Containers:
46// MRawEvtData
47// MRawRunHeader
48// MPedestalCam
49// [MCalibrationPattern]
50//
51// Output Containers:
52// MArrivalTimeCam
53// MExtractedSignalCam
54//
55//////////////////////////////////////////////////////////////////////////////
56#include "MExtractTimeAndChargeDigitalFilter.h"
57
58#include <errno.h>
59#include <fstream>
60
61#include <TH1.h>
62#include <TH2.h>
63#include <TMatrix.h>
64
65#include "MLog.h"
66#include "MLogManip.h"
67
68#include "MParList.h"
69
70#include "MRawRunHeader.h"
71#include "MCalibrationPattern.h"
72#include "MExtractedSignalCam.h"
73
74#include "MPedestalPix.h"
75
76ClassImp(MExtractTimeAndChargeDigitalFilter);
77
78using namespace std;
79
80const Byte_t MExtractTimeAndChargeDigitalFilter::fgHiGainFirst = 0;
81const Byte_t MExtractTimeAndChargeDigitalFilter::fgHiGainLast = 15;
82const Byte_t MExtractTimeAndChargeDigitalFilter::fgLoGainFirst = 1;
83const Byte_t MExtractTimeAndChargeDigitalFilter::fgLoGainLast = 14;
84const Int_t MExtractTimeAndChargeDigitalFilter::fgWindowSizeHiGain = 4;
85const Int_t MExtractTimeAndChargeDigitalFilter::fgWindowSizeLoGain = 6;
86const Int_t MExtractTimeAndChargeDigitalFilter::fgBinningResolutionHiGain = 10;
87const Int_t MExtractTimeAndChargeDigitalFilter::fgBinningResolutionLoGain = 10;
88const Int_t MExtractTimeAndChargeDigitalFilter::fgSignalStartBinHiGain = 4;
89const Int_t MExtractTimeAndChargeDigitalFilter::fgSignalStartBinLoGain = 4;
90const Float_t MExtractTimeAndChargeDigitalFilter::fgOffsetLoGain = 1.7; // 5 ns
91const Float_t MExtractTimeAndChargeDigitalFilter::fgLoGainStartShift = -1.8;
92
93// --------------------------------------------------------------------------
94//
95// Default constructor.
96//
97// Calls:
98// - SetWindowSize();
99// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
100// - SetBinningResolution();
101//
102// Sets all weights to 1.
103//
104MExtractTimeAndChargeDigitalFilter::MExtractTimeAndChargeDigitalFilter(const char *name, const char *title)
105 : fTimeShiftHiGain(0.), fTimeShiftLoGain(0.), fAutomaticWeights(kTRUE), fRandomIter(0)
106{
107 fName = name ? name : "MExtractTimeAndChargeDigitalFilter";
108 fTitle = title ? title : "Digital Filter";
109
110 SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
111 SetWindowSize();
112 SetBinningResolution();
113 SetSignalStartBin();
114
115 SetOffsetLoGain(fgOffsetLoGain);
116 SetLoGainStartShift(fgLoGainStartShift);
117}
118
119// ---------------------------------------------------------------------------------------
120//
121// Checks:
122// - if a window is bigger than the one defined by the ranges, set it to the available range
123//
124// Sets:
125// - fNumHiGainSamples to: (Float_t)fWindowSizeHiGain
126// - fNumLoGainSamples to: (Float_t)fWindowSizeLoGain
127//
128void MExtractTimeAndChargeDigitalFilter::SetWindowSize(Int_t windowh, Int_t windowl)
129{
130
131 if (windowh != fgWindowSizeHiGain)
132 *fLog << warn << GetDescriptor()
133 << ": ATTENTION!!! If you are not Hendrik Bartko, do NOT use a different window size than the default." << endl;
134 if (windowl != fgWindowSizeLoGain)
135 *fLog << warn << GetDescriptor()
136 << ": ATTENTION!!! If you are not Hendrik Bartko, do NOT use a different window size than the default" << endl;
137
138 fWindowSizeHiGain = windowh;
139 fWindowSizeLoGain = windowl;
140
141 const Int_t availhirange = (Int_t)(fHiGainLast-fHiGainFirst+1);
142
143 if (fWindowSizeHiGain > availhirange)
144 {
145 *fLog << warn << GetDescriptor() << ": Hi Gain window size: " << Form("%2i",fWindowSizeHiGain);
146 *fLog << " is bigger than available range: [" << Form("%2i", (int)fHiGainFirst);
147 *fLog << "," << Form("%21", (int)fHiGainLast) << "]" << endl;
148
149 fHiGainLast = fHiGainFirst + fWindowSizeHiGain;
150
151 *fLog << warn << GetDescriptor() << ": Will set the upper range to: " << (int)fHiGainLast << endl;
152 }
153
154 if (fWindowSizeHiGain < 2)
155 {
156 fWindowSizeHiGain = 2;
157 *fLog << warn << GetDescriptor() << ": High Gain window size set to two samples" << endl;
158 }
159
160 if (fLoGainLast != 0 && fWindowSizeLoGain != 0)
161 {
162 const Int_t availlorange = (Int_t)(fLoGainLast-fLoGainFirst+1);
163
164 if (fWindowSizeLoGain > availlorange)
165 {
166 *fLog << warn << GetDescriptor() << ": Lo Gain window size: " << Form("%2i",fWindowSizeLoGain);
167 *fLog << " is bigger than available range: [" << Form("%2i", (int)fLoGainFirst);
168 *fLog << "," << Form("%21", (int)fLoGainLast) << "]" << endl;
169
170 fLoGainLast = fLoGainFirst + fWindowSizeLoGain;
171
172 *fLog << warn << GetDescriptor() << ": Will set the upper range to: " << (int)fLoGainLast << endl;
173 }
174
175 if (fWindowSizeLoGain<2)
176 {
177 fWindowSizeLoGain = 2;
178 *fLog << warn << GetDescriptor() << ": Low Gain window size set to two samples" << endl;
179 }
180 }
181 //
182 // We need here the effective number of samples which is about 2.5 in the case of a window
183 // size of 6. The exact numbers have to be found still.
184 //
185 fNumHiGainSamples = (Float_t)fWindowSizeHiGain/2.4;
186 fNumLoGainSamples = (Float_t)fWindowSizeLoGain/2.4;
187 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
188 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
189
190}
191
192// --------------------------------------------------------------------------
193//
194// Executing MExtractTimeAndCharge::PreProcess and searching for
195// MCalibrationPattern
196//
197Int_t MExtractTimeAndChargeDigitalFilter::PreProcess(MParList *pList)
198{
199 if (!MExtractTimeAndCharge::PreProcess(pList))
200 return kFALSE;
201
202 fCalibPattern = (MCalibrationPattern*)pList->FindObject("MCalibrationPattern");
203 return kTRUE;
204}
205
206// --------------------------------------------------------------------------
207//
208// The weights are determined using GetAutimaticWeights().
209//
210// kFALSE is returned if it returned an error.
211// kTRUE is returned if no new weights were set.
212//
213// If new weights are set
214// fTimeShiftHiGain
215// fTimeShiftLoGain
216// fNumHiGainSamples
217// fNumLoGainSamples
218// fSqrtHiGainSamples
219// fSqrtLoGainSamples
220// and
221// fSignals->SetUsedFADCSlices(...)
222// is updated accordingly.
223//
224Bool_t MExtractTimeAndChargeDigitalFilter::GetWeights()
225{
226 switch (GetAutomaticWeights())
227 {
228 case kERROR: // An error occured
229 return kFALSE;
230 case kFALSE: // No new weights set
231 return kTRUE;
232 }
233
234 //
235 // shift the times back to the right reference (start counting from 0)
236 //
237 // The high-gain has the extraction offset (fHiGainFirst) already included
238 // here for speed reason. The low-gain has a changing extraction offset,
239 // so it will be added at every event (in FindTimeAndChargeLoGain)
240 fTimeShiftHiGain = 0.5 + 1./fBinningResolutionHiGain + fHiGainFirst;
241 fTimeShiftLoGain = 0.5 + 1./fBinningResolutionLoGain;
242
243 //
244 // We need here the effective number of samples which is about 2.5 in the case of a window
245 // size of 6. The exact numbers have to be found still.
246 //
247 fNumHiGainSamples = (Float_t)fWindowSizeHiGain/2.4;
248 fNumLoGainSamples = (Float_t)fWindowSizeLoGain/2.4;
249 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
250 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
251
252 // From MExtractTimeAndCharge::ReInit
253 if (fSignals)
254 fSignals->SetUsedFADCSlices(fHiGainFirst, fHiGainLast+fHiLoLast, fNumHiGainSamples,
255 fLoGainFirst, fLoGainLast, fNumLoGainSamples);
256 return kTRUE;
257}
258
259// --------------------------------------------------------------------------
260//
261// InitArrays
262//
263// Gets called in the ReInit() and initialized the arrays
264//
265Bool_t MExtractTimeAndChargeDigitalFilter::InitArrays()
266{
267 if (!fRunHeader)
268 return kFALSE;
269
270 const Int_t rangehi = (Int_t)(fHiGainLast - fHiGainFirst + 1 + fHiLoLast);
271 const Int_t rangelo = (Int_t)(fLoGainLast - fLoGainFirst + 1);
272
273 fHiGainSignal.Set(rangehi);
274 fLoGainSignal.Set(rangelo);
275
276 return GetWeights();
277}
278
279// --------------------------------------------------------------------------
280//
281// Check if reading a new weights file is necessary because the calibration
282// pattern has changed. (Cannot be done in ReInit, because at this time
283// the calibration pattern is not available.
284// Then process the event.
285//
286Int_t MExtractTimeAndChargeDigitalFilter::Process()
287{
288 // Change Weights if the calibration patter changes
289 if (!GetWeights())
290 return kERROR;
291
292 // Process event
293 return MExtractTimeAndCharge::Process();
294}
295
296void MExtractTimeAndChargeDigitalFilter::CalcBinningResArrays()
297{
298
299 fArrBinningResHiGain.Set(fWindowSizeHiGain);
300 fArrBinningResHalfHiGain.Set(fWindowSizeHiGain);
301
302 for (int i=0; i<fWindowSizeHiGain; i++)
303 {
304 fArrBinningResHiGain[i] = fBinningResolutionHiGain*i;
305 fArrBinningResHalfHiGain[i] = fArrBinningResHiGain[i] + fBinningResolutionHalfHiGain;
306 }
307
308 fArrBinningResLoGain.Set(fWindowSizeLoGain);
309 fArrBinningResHalfLoGain.Set(fWindowSizeLoGain);
310
311 for (int i=0; i<fWindowSizeLoGain; i++)
312 {
313 fArrBinningResLoGain[i] = fBinningResolutionLoGain*i;
314 fArrBinningResHalfLoGain[i] = fArrBinningResLoGain[i] + fBinningResolutionHalfLoGain;
315 }
316}
317
318// --------------------------------------------------------------------------
319//
320// Apply the digital filter algorithm to the high-gain slices.
321//
322void MExtractTimeAndChargeDigitalFilter::FindTimeAndChargeHiGain(Byte_t *ptr, Byte_t *logain, Float_t &sum, Float_t &dsum,
323 Float_t &time, Float_t &dtime,
324 Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
325{
326
327 Int_t range = fHiGainLast - fHiGainFirst + 1;
328
329 const Byte_t *end = ptr + range;
330 Byte_t *p = ptr;
331 Byte_t maxpos = 0;
332
333 //
334 // Preparations for the pedestal subtraction (with AB-noise correction)
335 //
336 const Float_t pedes = ped.GetPedestal();
337 const Float_t ABoffs = ped.GetPedestalABoffset();
338
339 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
340
341 range += fHiLoLast;
342 fMaxBinContent = 0;
343 //
344 // Check for saturation in all other slices
345 //
346 Int_t ids = fHiGainFirst;
347 Float_t *sample = fHiGainSignal.GetArray();
348 while (p<end)
349 {
350 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
351
352 if (*p > fMaxBinContent)
353 {
354 maxpos = p-ptr;
355 if (maxpos > 1 && maxpos < (range - fWindowSizeHiGain + 1))
356 fMaxBinContent = *p;
357 }
358
359 if (*p++ >= fSaturationLimit)
360 if (!sat)
361 sat = ids-3;
362 }
363
364 //
365 // Slide with a window of size fWindowSizeHiGain over the sample
366 // and multiply the entries with the corresponding weights
367 //
368 if (IsNoiseCalculation())
369 {
370 if (fRandomIter == fBinningResolutionHiGain)
371 fRandomIter = 0;
372 for (Int_t ids=0; ids < fWindowSizeHiGain; ids++)
373 {
374 const Int_t idx = fArrBinningResHiGain[ids] + fRandomIter;
375 sum += fAmpWeightsHiGain [idx]*fHiGainSignal[ids];
376 }
377 fRandomIter++;
378 return;
379 }
380
381 if (fHiLoLast != 0)
382 {
383
384 end = logain + fHiLoLast;
385
386 while (logain<end)
387 {
388
389 *sample++ = (Float_t)*logain - pedmean[(ids++ + abflag) & 0x1];
390
391 if (*logain++ >= fSaturationLimit)
392 if (!sat)
393 sat = ids-3;
394 }
395 }
396
397 //
398 // allow no saturated slice
399 //
400 // if (sat > 0)
401 // return;
402
403 Float_t time_sum = 0.;
404 Float_t fmax = -9999.;
405 Float_t ftime_max = 0.;
406 Int_t max_p = 0;
407
408 //
409 // Calculate the sum of the first fWindowSize slices
410 //
411 for (Int_t i=0;i<range-fWindowSizeHiGain+1;i++)
412 {
413 sum = 0.;
414 time_sum = 0.;
415
416 //
417 // Slide with a window of size fWindowSizeHiGain over the sample
418 // and multiply the entries with the corresponding weights
419 //
420 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
421 {
422 const Int_t idx = fBinningResolutionHiGain*sample+fBinningResolutionHalfHiGain;
423 const Float_t pex = fHiGainSignal[sample+i];
424 sum += fAmpWeightsHiGain [idx]*pex;
425 time_sum += fTimeWeightsHiGain[idx]*pex;
426 }
427
428 if (sum>fmax)
429 {
430 fmax = sum;
431 ftime_max = time_sum;
432 max_p = i;
433 }
434 } /* for (Int_t i=0;i<range-fWindowSizeHiGain+1;i++) */
435
436 if (fmax==0)
437 return;
438
439 ftime_max /= fmax;
440 Int_t t_iter = Int_t(ftime_max*fBinningResolutionHiGain);
441 Int_t sample_iter = 0;
442
443 while ( t_iter > fBinningResolutionHalfHiGain-1 || t_iter < -fBinningResolutionHalfHiGain )
444 {
445 if (t_iter > fBinningResolutionHalfHiGain-1)
446 {
447 t_iter -= fBinningResolutionHiGain;
448 max_p--;
449 sample_iter--;
450 }
451 if (t_iter < -fBinningResolutionHalfHiGain)
452 {
453 t_iter += fBinningResolutionHiGain;
454 max_p++;
455 sample_iter++;
456 }
457 }
458
459 sum = 0.;
460 time_sum = 0.;
461 if (max_p < 0)
462 max_p = 0;
463 if (max_p+fWindowSizeHiGain > range)
464 max_p = range-fWindowSizeHiGain;
465 //
466 // Slide with a window of size fWindowSizeHiGain over the sample
467 // and multiply the entries with the corresponding weights
468 //
469 for (Int_t sample=0; sample < fWindowSizeHiGain; sample++)
470 {
471 const Int_t idx = fArrBinningResHalfHiGain[sample] + t_iter;
472 const Int_t ids = max_p + sample;
473 const Float_t pex = fHiGainSignal[ids];
474 sum += fAmpWeightsHiGain [idx]*pex;
475 time_sum += fTimeWeightsHiGain[idx]*pex;
476 }
477
478 if (sum == 0)
479 return;
480
481 // here, the first high-gain slice is already included in the fTimeShiftHiGain
482 time = fTimeShiftHiGain + max_p - Float_t(t_iter)/fBinningResolutionHiGain;
483
484 const Float_t timefineadjust = time_sum/sum;
485
486 if (TMath::Abs(timefineadjust) < 4./fBinningResolutionHiGain)
487 time -= timefineadjust;
488
489}
490
491// --------------------------------------------------------------------------
492//
493// Apply the digital filter algorithm to the low-gain slices.
494//
495void MExtractTimeAndChargeDigitalFilter::FindTimeAndChargeLoGain(Byte_t *ptr, Float_t &sum, Float_t &dsum,
496 Float_t &time, Float_t &dtime,
497 Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
498{
499
500 const Int_t range = fLoGainLast - fLoGainFirst + 1;
501
502 const Byte_t *end = ptr + range;
503 Byte_t *p = ptr;
504 //
505 // Prepare the low-gain pedestal
506 //
507 const Float_t pedes = ped.GetPedestal();
508 const Float_t ABoffs = ped.GetPedestalABoffset();
509
510 const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
511
512 //
513 // Check for saturation in all other slices
514 //
515 Float_t *sample = fLoGainSignal.GetArray();
516 Int_t ids = fLoGainFirst;
517 while (p<end)
518 {
519 *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
520
521 if (*p++ >= fSaturationLimit)
522 sat++;
523 }
524
525 //
526 // Slide with a window of size fWindowSizeLoGain over the sample
527 // and multiply the entries with the corresponding weights
528 //
529 if (IsNoiseCalculation())
530 {
531 if (fRandomIter == fBinningResolutionLoGain)
532 fRandomIter = 0;
533 for (Int_t ids=0; ids < fWindowSizeLoGain; ids++)
534 {
535 const Int_t idx = fArrBinningResLoGain[ids] + fRandomIter;
536 sum += fAmpWeightsLoGain [idx]*fLoGainSignal[ids];
537 }
538 return;
539 }
540
541 Float_t time_sum = 0.;
542 Float_t fmax = 0.;
543 Float_t ftime_max = 0.;
544 Int_t max_p = 0;
545
546 //
547 // Calculate the sum of the first fWindowSize slices
548 //
549 for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++)
550 {
551 sum = 0.;
552 time_sum = 0.;
553
554 //
555 // Slide with a window of size fWindowSizeLoGain over the sample
556 // and multiply the entries with the corresponding weights
557 //
558 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
559 {
560 const Int_t idx = fArrBinningResHalfLoGain[sample];
561 const Float_t pex = fLoGainSignal[sample+i];
562 sum += fAmpWeightsLoGain [idx]*pex;
563 time_sum += fTimeWeightsLoGain[idx]*pex;
564 }
565
566 if (sum>fmax)
567 {
568 fmax = sum;
569 ftime_max = time_sum;
570 max_p = i;
571 }
572 } /* for (Int_t i=0;i<range-fWindowSizeLoGain+1;i++) */
573
574 time = 0;
575 if (fmax==0)
576 return;
577
578 ftime_max /= fmax;
579 Int_t t_iter = Int_t(ftime_max*fBinningResolutionLoGain);
580 Int_t sample_iter = 0;
581
582 while ( t_iter > fBinningResolutionHalfLoGain-1 || t_iter < -fBinningResolutionHalfLoGain )
583 {
584 if (t_iter > fBinningResolutionHalfLoGain-1)
585 {
586 t_iter -= fBinningResolutionLoGain;
587 max_p--;
588 sample_iter--;
589 }
590 if (t_iter < -fBinningResolutionHalfLoGain)
591 {
592 t_iter += fBinningResolutionLoGain;
593 max_p++;
594 sample_iter++;
595 }
596 }
597
598 sum = 0.;
599 time_sum = 0.;
600
601 //
602 // Slide with a window of size fWindowSizeLoGain over the sample
603 // and multiply the entries with the corresponding weights
604 //
605 for (Int_t sample=0; sample < fWindowSizeLoGain; sample++)
606 {
607 const Int_t idx = fArrBinningResHalfLoGain[sample] + t_iter;
608 const Int_t ids = max_p + sample;
609 const Float_t pex = ids < 0 ? 0. : ( ids >= range ? 0. : fLoGainSignal[ids]);
610 sum += fAmpWeightsLoGain [idx]*pex;
611 time_sum += fTimeWeightsLoGain[idx]*pex;
612 }
613
614 if (sum == 0)
615 return;
616
617 time = max_p + fTimeShiftLoGain + (Float_t)fLoGainFirst /* this shifts the time to the start of the rising edge */
618 - ((Float_t)t_iter)/fBinningResolutionLoGain;
619
620 const Float_t timefineadjust = time_sum/sum;
621
622 if (TMath::Abs(timefineadjust) < 4./fBinningResolutionLoGain)
623 time -= timefineadjust;
624
625}
626
627// --------------------------------------------------------------------------
628//
629// Read the setup from a TEnv, eg:
630// MJPedestal.MExtractor.WindowSizeHiGain: 6
631// MJPedestal.MExtractor.WindowSizeLoGain: 6
632// MJPedestal.MExtractor.BinningResolutionHiGain: 10
633// MJPedestal.MExtractor.BinningResolutionLoGain: 10
634// MJPedestal.MExtractor.WeightsFile: filename
635// MJPedestal.MExtractor.AutomaticWeights: off
636//
637Int_t MExtractTimeAndChargeDigitalFilter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
638{
639
640 Byte_t hw = fWindowSizeHiGain;
641 Byte_t lw = fWindowSizeLoGain;
642 Bool_t rc = kFALSE;
643
644 if (IsEnvDefined(env, prefix, "AutomaticWeights", print))
645 {
646 EnableAutomaticWeights(GetEnvValue(env, prefix, "AutomaticWeights", fAutomaticWeights));
647 rc = kTRUE;
648 }
649
650 if (IsEnvDefined(env, prefix, "WindowSizeHiGain", print))
651 {
652 hw = GetEnvValue(env, prefix, "WindowSizeHiGain", hw);
653 rc = kTRUE;
654 }
655 if (IsEnvDefined(env, prefix, "WindowSizeLoGain", print))
656 {
657 lw = GetEnvValue(env, prefix, "WindowSizeLoGain", lw);
658 rc = kTRUE;
659 }
660
661 if (rc)
662 SetWindowSize(hw, lw);
663
664 Bool_t rc2 = kFALSE;
665 Int_t brh = fBinningResolutionHiGain;
666 Int_t brl = fBinningResolutionLoGain;
667
668 if (IsEnvDefined(env, prefix, "BinningResolutionHiGain", print))
669 {
670 brh = GetEnvValue(env, prefix, brh);
671 rc2 = kTRUE;
672 }
673 if (IsEnvDefined(env, prefix, "BinningResolutionLoGain", print))
674 {
675 brl = GetEnvValue(env, prefix, brl);
676 rc2 = kTRUE;
677 }
678
679 if (rc2)
680 {
681 SetBinningResolution(brh, brl);
682 rc = kTRUE;
683 }
684
685 if (IsEnvDefined(env, prefix, "WeightsFile", print))
686 {
687 SetNameWeightsFile(GetEnvValue(env, prefix, "WeightsFile", ""));
688 rc = kTRUE;
689 }
690
691 return MExtractTimeAndCharge::ReadEnv(env, prefix, print) ? kTRUE : rc;
692}
693
694//----------------------------------------------------------------------------
695//
696// If automatic weights are requested, no default weights (name.IsNull())
697// are requested, fRunHeader is available and fRunHeader->IsMonteCarloRun()
698// is true prepend "MC_" in from of the name.
699//
700// return poth+name;
701//
702TString MExtractTimeAndChargeDigitalFilter::CompileWeightFileName(TString path, const TString &name) const
703{
704 if (fAutomaticWeights && !name.IsNull() && fRunHeader && fRunHeader->IsMonteCarloRun())
705 path += "MC_";
706
707 path += name;
708
709 return path;
710}
711
712//----------------------------------------------------------------------------
713//
714// Read a pre-defined weights file into the class.
715// This is mandatory for the extraction
716//
717// If filenname is empty, then all weights will be set to 1.
718//
719// Returns:
720// kTRUE: new weights set
721// kFALSE: no weights set
722// kERROR: error
723//
724Int_t MExtractTimeAndChargeDigitalFilter::ReadWeightsFile(TString filename, TString path)
725{
726 if (filename.IsNull())
727 {
728 fAmpWeightsHiGain .Set(fBinningResolutionHiGain*fWindowSizeHiGain);
729 fAmpWeightsLoGain .Set(fBinningResolutionLoGain*fWindowSizeLoGain);
730 fTimeWeightsHiGain.Set(fBinningResolutionHiGain*fWindowSizeHiGain);
731 fTimeWeightsLoGain.Set(fBinningResolutionLoGain*fWindowSizeLoGain);
732
733 fAmpWeightsHiGain.Reset(1);
734 fTimeWeightsHiGain.Reset(1);
735 fAmpWeightsLoGain.Reset(1);
736 fTimeWeightsLoGain.Reset(1);
737 return kTRUE;
738 }
739
740 // Add "MC_" in front of the filename if necessary
741 filename = CompileWeightFileName(path, filename);
742
743 //filename = MJob::ExpandPath(filename);
744
745 if (fNameWeightsFileSet==filename)
746 return kFALSE; // No file read
747
748 ifstream fin(filename.Data());
749 if (!fin)
750 {
751 *fLog << err << GetDescriptor() << ": ERROR - Cannot open file " << filename << ": ";
752 *fLog << strerror(errno) << endl;
753 return kERROR;
754 }
755
756 *fLog << all << GetDescriptor() << ": Reading weights in " << filename << "..." << flush;
757
758 Int_t len = 0;
759 Int_t cnt = 0;
760 Int_t line = 0;
761 Bool_t hi = kFALSE;
762 Bool_t lo = kFALSE;
763
764 TString str;
765
766 while (1)
767 {
768 str.ReadLine(fin);
769 if (!fin)
770 break;
771
772 line++;
773
774 if (str.Contains("# High Gain Weights:"))
775 {
776 if (hi)
777 {
778 *fLog << err << "ERROR - 'High Gain Weights' found twice in line #" << line << "." << endl;
779 return kERROR;
780 }
781
782 if (2!=sscanf(str.Data(), "# High Gain Weights:%2i %2i", &fWindowSizeHiGain, &fBinningResolutionHiGain))
783 {
784 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
785 *fLog << str << endl;
786 return kERROR;
787 }
788
789 len = fBinningResolutionHiGain*fWindowSizeHiGain;
790 fAmpWeightsHiGain .Set(len);
791 fTimeWeightsHiGain.Set(len);
792 hi = kTRUE;
793 continue;
794 }
795
796 if (str.Contains("# Low Gain Weights:"))
797 {
798 if (lo)
799 {
800 *fLog << err << "ERROR - 'Lo Gain Weights' found twice in line #" << line << "." << endl;
801 return kERROR;
802 }
803
804 if (2!=sscanf(str.Data(),"# Low Gain Weights:%2i %2i", &fWindowSizeLoGain, &fBinningResolutionLoGain))
805 {
806 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
807 *fLog << str << endl;
808 return kERROR;
809 }
810
811 len = fBinningResolutionLoGain*fWindowSizeLoGain;
812 fAmpWeightsLoGain .Set(len);
813 fTimeWeightsLoGain.Set(len);
814 lo = kTRUE;
815 continue;
816 }
817
818 // Handle lines with comments
819 if (str.Contains("#"))
820 continue;
821
822 // Nothing found so far
823 if (len == 0)
824 continue;
825
826 if (2!=sscanf(str.Data(), "%f %f",
827 lo ? &fAmpWeightsLoGain [cnt] : &fAmpWeightsHiGain [cnt],
828 lo ? &fTimeWeightsLoGain[cnt] : &fTimeWeightsHiGain[cnt]))
829 {
830 *fLog << err << "ERROR - Wrong number of arguments in line #" << line << ":" << endl;
831 *fLog << str << endl;
832 return kERROR;
833 }
834
835 if (++cnt == len)
836 {
837 len = 0;
838 cnt = 0;
839 }
840 }
841
842 if (cnt != len)
843 {
844 *fLog << err << "ERROR - Size mismatch in weights file " << filename << endl;
845 return kERROR;
846 }
847
848 if (!hi)
849 {
850 *fLog << err << "ERROR - No correct header found in weights file " << filename << endl;
851 return kERROR;
852 }
853
854 *fLog << "done." << endl;
855
856 *fLog << inf << " File contains " << fWindowSizeHiGain << " hi-gain slices ";
857 *fLog << "with a resolution of " << fBinningResolutionHiGain << endl;
858
859 *fLog << inf << " File contains " << fWindowSizeLoGain << " lo-gain slices ";
860 *fLog << "with a resolution of " << fBinningResolutionLoGain << endl;
861
862 CalcBinningResArrays();
863
864 switch (fWindowSizeHiGain)
865 {
866 case 4:
867 SetResolutionPerPheHiGain(0.036);
868 break;
869 case 6:
870 SetResolutionPerPheHiGain(0.021);
871 break;
872 default:
873 *fLog << warn << "Could not set the high-gain extractor resolution per phe for window size "
874 << fWindowSizeHiGain << endl;
875 }
876
877 switch (fWindowSizeLoGain)
878 {
879 case 4:
880 SetResolutionPerPheLoGain(0.005);
881 break;
882 case 6:
883 SetResolutionPerPheLoGain(0.004);
884 break;
885 default:
886 *fLog << warn << "Could not set the low-gain extractor resolution per phe for window size "
887 << fWindowSizeLoGain << endl;
888 }
889
890 fNameWeightsFileSet = filename;
891
892 return kTRUE;
893}
894
895
896//----------------------------------------------------------------------------
897//
898// The default (+ prepending possible "MC_") is read for:
899//
900// - RunType: Pedestal (independant of fAutomaticWeights)
901// - fAutomaticWeights disabled
902//
903// if fAutomaticWeights enabled:
904// - fNameWeightsFile.IsNull()
905// - !fCalibPattern
906// - fCalibPattern->GetPulserColor()==MCalibrationCam::kNONE
907//
908// If automatic weights are enabled, the case above didn't take place and
909// fNameWeightsFile starts with "calibration_weights_"
910// - the color (blue, UV) is replaced by the appropriate one
911// taken from the calibration pattern
912//
913// In most cases a debug output is printed. Further output about the color
914// determination can be switched on with debug level > 5;
915//
916// Returns:
917// kFALSE: No new weights set
918// kTRUE: New weights set
919// kERROR: Error
920//
921Int_t MExtractTimeAndChargeDigitalFilter::GetAutomaticWeights()
922{
923 const Ssiz_t pos = fNameWeightsFile.Last('/')+1;
924 const Ssiz_t len = fNameWeightsFile.Length();
925
926 // Split file name in path and name
927 TString path = fNameWeightsFile(0, pos>=0?pos:len);
928 TString name = fNameWeightsFile(pos>=0?pos:0, len);
929
930 // Remove trailing "MC_" for automatic weights
931 if (fAutomaticWeights && name.BeginsWith("MC_"))
932 name.Remove(0, 3);
933
934 // In case of a pedetsal run no calibration pattern can be available
935 // the default weights are always used.
936 if (fRunHeader->GetRunType()==MRawRunHeader::kRTPedestal)
937 {
938 *fLog << dbg << "Pedestal file... using default weights: " << fNameWeightsFile << endl;
939 return ReadWeightsFile(name, path);
940 }
941
942 // If automatic weights are switched off use default weights
943 if (!fAutomaticWeights)
944 {
945 *fLog << dbg << "Automatic weights switched off... using default weights: " << fNameWeightsFile << endl;
946 return ReadWeightsFile(name, path);
947 }
948
949 // If automatic weights are switched on but no filename is given raise error
950 if (fNameWeightsFile.IsNull())
951 {
952 *fLog << err << "ERROR - Cannot get automatic weights without default filename." << endl;
953 return kERROR;
954 }
955
956 // If this is no pedestal run, automatic weights are requested and a
957 // filename for the weights file is given pedestal-extraction from
958 // cosmics data is assumed.
959 if (!fCalibPattern)
960 {
961 *fLog << dbg << "No decoded calibration pattern available... using default weights: " << fNameWeightsFile << endl;
962 return ReadWeightsFile(name, path);
963 }
964
965 const Bool_t debug = gLog.GetDebugLevel()>5;
966
967 // If no calibration pattern is available do not change the
968 // current weighs or current weights file name.
969 if (fCalibPattern->GetPulserColor()==MCalibrationCam::kNONE)
970 {
971 // If we are extracting data and the calibration pattern is kNONE
972 // we assume that it is a data file without interleaved events
973 // and calibration pattern information available.
974 if ((fRunHeader->GetRunType()!=MRawRunHeader::kRTData && !fRunHeader->IsMonteCarloRun()) || debug)
975 *fLog << dbg << "No calibration color set so far... guessing default: " << fNameWeightsFile << endl;
976
977 return ReadWeightsFile(name, path);
978 }
979
980 if (debug)
981 {
982 *fLog << dbg << endl;
983 *fLog << underline << GetDescriptor() << endl;
984 *fLog << " Trying to get automatic weight for " << fNameWeightsFile << endl;
985 *fLog << " Run type: ";
986 }
987
988 if (name.BeginsWith("calibration_weights_") && fCalibPattern)
989 {
990 if (debug)
991 *fLog << " Calibration with color " << fCalibPattern->GetPulserColorStr() << ", setting ";
992 switch (fCalibPattern->GetPulserColor())
993 {
994 case MCalibrationCam::kBLUE: // 2
995 case MCalibrationCam::kGREEN: // 1
996 if (debug)
997 *fLog << "blue/green, ";
998 name.ReplaceAll("UV", "blue");
999 break;
1000
1001 case MCalibrationCam::kUV: // 3
1002 case MCalibrationCam::kCT1: // 0
1003 if (debug)
1004 *fLog << "UV/CT1, ";
1005 name.ReplaceAll("blue", "UV");
1006 break;
1007 case MCalibrationCam::kNONE:
1008 break;
1009 default: // kNone + etc
1010 *fLog << err << "ERROR - Cannot get automatic weights for " << fCalibPattern->GetPulserColorStr() << endl;
1011 return kERROR;
1012 }
1013 }
1014
1015 return ReadWeightsFile(name, path);
1016}
1017
1018//----------------------------------------------------------------------------
1019//
1020// Create the weights file
1021// Beware that the shape-histogram has to contain the pulse starting at bin 1
1022//
1023Bool_t MExtractTimeAndChargeDigitalFilter::WriteWeightsFile(TString filename, TH1F *shapehi, TH2F *autocorrhi,
1024 TH1F *shapelo, TH2F *autocorrlo )
1025{
1026
1027 const Int_t nbinshi = shapehi->GetNbinsX();
1028 Float_t binwidth = shapehi->GetBinWidth(1);
1029
1030 TH1F *derivativehi = new TH1F(Form("%s%s",shapehi->GetName(),"_der"),
1031 Form("%s%s",shapehi->GetTitle()," derivative"),
1032 nbinshi,
1033 shapehi->GetBinLowEdge(1),
1034 shapehi->GetBinLowEdge(nbinshi)+binwidth);
1035
1036 //
1037 // Calculate the derivative of shapehi
1038 //
1039 for (Int_t i = 1; i<nbinshi+1;i++)
1040 {
1041 derivativehi->SetBinContent(i,
1042 ((shapehi->GetBinContent(i+1)-shapehi->GetBinContent(i-1))/2./binwidth));
1043 derivativehi->SetBinError(i,
1044 (sqrt(shapehi->GetBinError(i+1)*shapehi->GetBinError(i+1)
1045 +shapehi->GetBinError(i-1)*shapehi->GetBinError(i-1))/2./binwidth));
1046 }
1047
1048 //
1049 // normalize the shapehi, such that the integral for fWindowSize slices is one!
1050 //
1051 Float_t sum = 0;
1052 Int_t lasttemp = fBinningResolutionHiGain * (fSignalStartBinHiGain + fWindowSizeHiGain);
1053 lasttemp = lasttemp > nbinshi ? nbinshi : lasttemp;
1054
1055 for (Int_t i=fBinningResolutionHiGain*fSignalStartBinHiGain; i<lasttemp; i++) {
1056 sum += shapehi->GetBinContent(i);
1057 }
1058 sum /= fBinningResolutionHiGain;
1059
1060 shapehi->Scale(1./sum);
1061 derivativehi->Scale(1./sum);
1062
1063 //
1064 // read in the noise auto-correlation function:
1065 //
1066 TMatrix Bhi(fWindowSizeHiGain,fWindowSizeHiGain);
1067
1068 for (Int_t i=0; i<fWindowSizeHiGain; i++){
1069 for (Int_t j=0; j<fWindowSizeHiGain; j++){
1070 Bhi[i][j]=autocorrhi->GetBinContent(i+1,j+1); //+fSignalStartBinHiGain +fSignalStartBinHiGain
1071 }
1072 }
1073 Bhi.Invert();
1074
1075 const Int_t nsizehi = fWindowSizeHiGain*fBinningResolutionHiGain;
1076 fAmpWeightsHiGain.Set(nsizehi);
1077 fTimeWeightsHiGain.Set(nsizehi);
1078
1079 //
1080 // Loop over relative time in one BinningResolution interval
1081 //
1082 Int_t start = fBinningResolutionHiGain*(fSignalStartBinHiGain + 1);
1083
1084 for (Int_t i = -fBinningResolutionHalfHiGain+1; i<=fBinningResolutionHalfHiGain; i++)
1085 {
1086
1087 TMatrix g(fWindowSizeHiGain,1);
1088 TMatrix gT(1,fWindowSizeHiGain);
1089 TMatrix d(fWindowSizeHiGain,1);
1090 TMatrix dT(1,fWindowSizeHiGain);
1091
1092 for (Int_t count=0; count < fWindowSizeHiGain; count++){
1093
1094 g[count][0]=shapehi->GetBinContent(start
1095 +fBinningResolutionHiGain*count+i);
1096 gT[0][count]=shapehi->GetBinContent(start
1097 +fBinningResolutionHiGain*count+i);
1098 d[count][0]=derivativehi->GetBinContent(start
1099 +fBinningResolutionHiGain*count+i);
1100 dT[0][count]=derivativehi->GetBinContent(start
1101 +fBinningResolutionHiGain*count+i);
1102 }
1103
1104 TMatrix m_denom = (gT*(Bhi*g))*(dT*(Bhi*d)) - (dT*(Bhi*g))*(dT*(Bhi*g));
1105 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
1106
1107 TMatrix m_first = dT*(Bhi*d); // ROOT thinks, m_first is still a matrix
1108 Float_t first = m_first[0][0]/denom;
1109
1110 TMatrix m_last = gT*(Bhi*d); // ROOT thinks, m_last is still a matrix
1111 Float_t last = m_last[0][0]/denom;
1112
1113 TMatrix m1 = gT*Bhi;
1114 m1 *= first;
1115
1116 TMatrix m2 = dT*Bhi;
1117 m2 *=last;
1118
1119 TMatrix w_amp = m1 - m2;
1120
1121 TMatrix m_first1 = gT*(Bhi*g);
1122 Float_t first1 = m_first1[0][0]/denom;
1123
1124 TMatrix m_last1 = gT*(Bhi*d);
1125 Float_t last1 = m_last1 [0][0]/denom;
1126
1127 TMatrix m11 = dT*Bhi;
1128 m11 *=first1;
1129
1130 TMatrix m21 = gT*Bhi;
1131 m21 *=last1;
1132
1133 TMatrix w_time= m11 - m21;
1134
1135 for (Int_t count=0; count < fWindowSizeHiGain; count++)
1136 {
1137 const Int_t idx = i+fBinningResolutionHalfHiGain+fBinningResolutionHiGain*count-1;
1138 fAmpWeightsHiGain [idx] = w_amp [0][count];
1139 fTimeWeightsHiGain[idx] = w_time[0][count];
1140 }
1141
1142 } // end loop over i
1143
1144 //
1145 // Low Gain histograms
1146 //
1147 TH1F *derivativelo = NULL;
1148 if (shapelo)
1149 {
1150 const Int_t nbinslo = shapelo->GetNbinsX();
1151 binwidth = shapelo->GetBinWidth(1);
1152
1153 derivativelo = new TH1F(Form("%s%s",shapelo->GetName(),"_der"),
1154 Form("%s%s",shapelo->GetTitle()," derivative"),
1155 nbinslo,
1156 shapelo->GetBinLowEdge(1),
1157 shapelo->GetBinLowEdge(nbinslo)+binwidth);
1158
1159 //
1160 // Calculate the derivative of shapelo
1161 //
1162 for (Int_t i = 1; i<nbinslo+1;i++)
1163 {
1164 derivativelo->SetBinContent(i,
1165 ((shapelo->GetBinContent(i+1)-shapelo->GetBinContent(i-1))/2./binwidth));
1166 derivativelo->SetBinError(i,
1167 (sqrt(shapelo->GetBinError(i+1)*shapelo->GetBinError(i+1)
1168 +shapelo->GetBinError(i-1)*shapelo->GetBinError(i-1))/2./binwidth));
1169 }
1170
1171 //
1172 // normalize the shapelo, such that the integral for fWindowSize slices is one!
1173 //
1174 sum = 0;
1175 lasttemp = fBinningResolutionLoGain * (fSignalStartBinLoGain + fWindowSizeLoGain);
1176 lasttemp = lasttemp > nbinslo ? nbinslo : lasttemp;
1177
1178 for (Int_t i=fBinningResolutionLoGain*fSignalStartBinLoGain; i<lasttemp; i++)
1179 sum += shapelo->GetBinContent(i);
1180
1181 sum /= fBinningResolutionLoGain;
1182
1183 shapelo->Scale(1./sum);
1184 derivativelo->Scale(1./sum);
1185
1186 //
1187 // read in the noise auto-correlation function:
1188 //
1189 TMatrix Blo(fWindowSizeLoGain,fWindowSizeLoGain);
1190
1191 for (Int_t i=0; i<fWindowSizeLoGain; i++){
1192 for (Int_t j=0; j<fWindowSizeLoGain; j++){
1193 Blo[i][j]=autocorrlo->GetBinContent(i+1+fSignalStartBinLoGain,j+1+fSignalStartBinLoGain);
1194 }
1195 }
1196 Blo.Invert();
1197
1198 const Int_t nsizelo = fWindowSizeLoGain*fBinningResolutionLoGain;
1199 fAmpWeightsLoGain.Set(nsizelo);
1200 fTimeWeightsLoGain.Set(nsizelo);
1201
1202 //
1203 // Loop over relative time in one BinningResolution interval
1204 //
1205 Int_t start = fBinningResolutionLoGain*fSignalStartBinLoGain + fBinningResolutionHalfLoGain;
1206
1207 for (Int_t i = -fBinningResolutionHalfLoGain+1; i<=fBinningResolutionHalfLoGain; i++)
1208 {
1209
1210 TMatrix g(fWindowSizeLoGain,1);
1211 TMatrix gT(1,fWindowSizeLoGain);
1212 TMatrix d(fWindowSizeLoGain,1);
1213 TMatrix dT(1,fWindowSizeLoGain);
1214
1215 for (Int_t count=0; count < fWindowSizeLoGain; count++){
1216
1217 g[count][0] = shapelo->GetBinContent(start
1218 +fBinningResolutionLoGain*count+i);
1219 gT[0][count]= shapelo->GetBinContent(start
1220 +fBinningResolutionLoGain*count+i);
1221 d[count][0] = derivativelo->GetBinContent(start
1222 +fBinningResolutionLoGain*count+i);
1223 dT[0][count]= derivativelo->GetBinContent(start
1224 +fBinningResolutionLoGain*count+i);
1225 }
1226
1227 TMatrix m_denom = (gT*(Blo*g))*(dT*(Blo*d)) - (dT*(Blo*g))*(dT*(Blo*g));
1228 Float_t denom = m_denom[0][0]; // ROOT thinks, m_denom is still a matrix
1229
1230 TMatrix m_first = dT*(Blo*d); // ROOT thinks, m_first is still a matrix
1231 Float_t first = m_first[0][0]/denom;
1232
1233 TMatrix m_last = gT*(Blo*d); // ROOT thinks, m_last is still a matrix
1234 Float_t last = m_last[0][0]/denom;
1235
1236 TMatrix m1 = gT*Blo;
1237 m1 *= first;
1238
1239 TMatrix m2 = dT*Blo;
1240 m2 *=last;
1241
1242 TMatrix w_amp = m1 - m2;
1243
1244 TMatrix m_first1 = gT*(Blo*g);
1245 Float_t first1 = m_first1[0][0]/denom;
1246
1247 TMatrix m_last1 = gT*(Blo*d);
1248 Float_t last1 = m_last1 [0][0]/denom;
1249
1250 TMatrix m11 = dT*Blo;
1251 m11 *=first1;
1252
1253 TMatrix m21 = gT*Blo;
1254 m21 *=last1;
1255
1256 TMatrix w_time= m11 - m21;
1257
1258 for (Int_t count=0; count < fWindowSizeLoGain; count++)
1259 {
1260 const Int_t idx = i+fBinningResolutionHalfLoGain+fBinningResolutionLoGain*count-1;
1261 fAmpWeightsLoGain [idx] = w_amp [0][count];
1262 fTimeWeightsLoGain[idx] = w_time[0][count];
1263 }
1264
1265 } // end loop over i
1266 }
1267
1268 ofstream fn(filename.Data());
1269
1270 fn << "# High Gain Weights: " << fWindowSizeHiGain << " " << fBinningResolutionHiGain << endl;
1271 fn << "# (Amplitude) (Time) " << endl;
1272
1273 for (Int_t i=0; i<nsizehi; i++)
1274 fn << "\t" << fAmpWeightsHiGain[i] << "\t" << fTimeWeightsHiGain[i] << endl;
1275
1276 fn << "# Low Gain Weights: " << fWindowSizeLoGain << " " << fBinningResolutionLoGain << endl;
1277 fn << "# (Amplitude) (Time) " << endl;
1278
1279 for (Int_t i=0; i<nsizehi; i++)
1280 fn << "\t" << fAmpWeightsLoGain[i] << "\t" << fTimeWeightsLoGain[i] << endl;
1281
1282 delete derivativehi;
1283 if (derivativelo)
1284 delete derivativelo;
1285
1286 return kTRUE;
1287}
1288
1289void MExtractTimeAndChargeDigitalFilter::Print(Option_t *o) const
1290{
1291 if (IsA()==Class())
1292 *fLog << GetDescriptor() << ":" << endl;
1293
1294 MExtractTimeAndCharge::Print(o);
1295 *fLog << " Time Shift HiGain: " << fTimeShiftHiGain << " LoGain: " << fTimeShiftLoGain << endl;
1296 *fLog << " Window Size HiGain: " << fWindowSizeHiGain << " LoGain: " << fWindowSizeLoGain << endl;
1297 *fLog << " Binning Res HiGain: " << fBinningResolutionHiGain << " LoGain: " << fBinningResolutionHiGain << endl;
1298 *fLog << " Weights File desired: " << (fNameWeightsFile.IsNull()?"-":fNameWeightsFile) << endl;
1299 if (!fNameWeightsFileSet.IsNull())
1300 *fLog << " Weights File set: " << fNameWeightsFileSet << endl;
1301
1302 TString opt(o);
1303 if (!opt.Contains("weights"))
1304 return;
1305
1306 *fLog << endl;
1307 *fLog << inf << "Using the following weights: " << endl;
1308 *fLog << "Hi-Gain:" << endl;
1309 for (Int_t i=0; i<fBinningResolutionHiGain*fWindowSizeHiGain; i++)
1310 *fLog << " " << fAmpWeightsHiGain[i] << " \t " << fTimeWeightsHiGain[i] << endl;
1311
1312 *fLog << "Lo-Gain:" << endl;
1313 for (Int_t i=0; i<fBinningResolutionLoGain*fWindowSizeLoGain; i++)
1314 *fLog << " " << fAmpWeightsLoGain[i] << " \t " << fTimeWeightsLoGain[i] << endl;
1315}
Note: See TracBrowser for help on using the repository browser.