Belle II Software  release-08-01-10
eclAutocovarianceCalibrationC3Collector.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 //This module`
10 #include <ecl/modules/eclAutocovarianceCalibrationC3Collector/eclAutocovarianceCalibrationC3Collector.h>
11 
12 #include <iostream>
13 //Framework
14 #include <framework/dataobjects/EventMetaData.h>
15 
16 //ECL
17 #include <ecl/dataobjects/ECLDsp.h>
18 #include <ecl/dbobjects/ECLCrystalCalib.h>
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 //-----------------------------------------------------------------
24 // Register the Modules
25 //-----------------------------------------------------------------
26 REG_MODULE(eclAutocovarianceCalibrationC3Collector);
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
31 // constructor
32 eclAutocovarianceCalibrationC3CollectorModule::eclAutocovarianceCalibrationC3CollectorModule() : CalibrationCollectorModule(),
33  m_ECLAutocovarianceCalibrationC1Threshold("ECLAutocovarianceCalibrationC1Threshold")
34 {
35  // Set module properties
36  setDescription("Module to compute covarience matrix using waveforms from delayed Bhabha events");
38 }
39 
41 {
42 
45  CovarianceMatrixInfoVsCrysID = new TH2F("CovarianceMatrixInfoVsCrysID", "", ECLElementNumbers::c_NCrystals, 0,
47 
48  registerObject<TH2F>("CovarianceMatrixInfoVsCrysID", CovarianceMatrixInfoVsCrysID);
49 
50  m_eclDsps.isRequired();
51 
52 }
53 
55 {
57  B2INFO("eclAutocovarianceCalibrationC3Collector: Experiment = " << m_evtMetaData->getExperiment() << " run = " <<
58  m_evtMetaData->getRun());
59 
61 
62 }
63 
64 
66 {
67 
68  std::vector<float> m_tempArray;
70  const int NumDsp = m_eclDsps.getEntries();
71 
72  //Random Trigger Events have waveform for each crystal
73  if (NumDsp == ECLElementNumbers::c_NCrystals) {
74 
75  for (auto& aECLDsp : m_eclDsps) {
76 
77  const int id = aECLDsp.getCellId() - 1;
78 
79  //Peak to peak amplitude used to gauge noise level
80  float PeakToPeak = (float) aECLDsp.computePeaktoPeakAmp();
81 
82  if (PeakToPeak < m_PeakToPeakThresholds[id]) {
83 
84  float baseline = 0;
85  for (int i = 0; i < m_BaselineLimit; i++) baseline += aECLDsp.getDspA()[i];
86  baseline = baseline * (1. / ((float) m_BaselineLimit));
87 
88  // Computing baseline subtracted waveform
89  m_tempArray.clear();
90  for (int i = 0; i < m_nADCWaveformPoints; i++) m_tempArray.push_back(aECLDsp.getDspA()[i] - baseline);
91 
92  /*
93  Info on convariance matrix calucation:
94 
95  We store only nADCWaveformPoints numbers because there are several symmetries:
96  - In general the matrix is nADCWaveformPointsxnADCWaveformPoints, however, the matrix symmetric so we only compute the upper part of the matrix.
97  - We also assume the noise is random in time. With this assumption the values of each matrix diagonal will be the same. This allows us to define the matrix with only nADCWaveformPoints numbers: One for each diagonal (starting with the main diagonal)
98  I divide by (nADCWaveformPoints - (i-j)) to average over the number of off-diagonal elements for each row.
99  Using the final nADCWaveformPoints numbers the full nADCWaveformPointsxnADCWaveformPoints matrix in constructed by:
100  Element 0 is the nADCWaveformPoints main diagonals
101  Element 1 is the value of the 30 off diagonals
102  */
103 
104  for (int i = 0; i < m_nADCWaveformPoints; i++) {
105 
106  float value_i = m_tempArray[i];
107 
108  for (int j = i; j < m_nADCWaveformPoints; j++) {
109 
110  int tempIndex = abs(i - j);
111 
112  m_CovarianceMatrixInfoVsCrysIDHistogram[id][tempIndex] += (value_i * m_tempArray[j]);
113 
114  }
115  }
117  }
118  }
119  }
120 
121 }
122 
124 {
125  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++) {
126  for (int j = 0; j < (m_nADCWaveformPoints + 1); j++) {
127  getObjectPtr<TH2>("CovarianceMatrixInfoVsCrysID")->SetBinContent(i + 1, j + 1, m_CovarianceMatrixInfoVsCrysIDHistogram[i][j]);
128  }
129  }
130 }
Calibration collector module base class.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
StoreArray< ECLDsp > m_eclDsps
Required input array of ECLDSPs.
float m_CovarianceMatrixInfoVsCrysIDHistogram[ECLElementNumbers::c_NCrystals][m_nADCWaveformPoints+1]
container for coveriance matrix
TH2F * CovarianceMatrixInfoVsCrysID
result returned by collector that contains the coveriance matrix for each crystal
DBObjPtr< ECLCrystalCalib > m_ECLAutocovarianceCalibrationC1Threshold
thresholds obtained from C1 stage
void collect() override
Select events and crystals and accumulate histograms.
StoreObjPtr< EventMetaData > m_evtMetaData
dataStore EventMetaData
std::vector< float > m_PeakToPeakThresholds
vector of thresholds obtained from DB object
void closeRun() override
Transfer fom array container to ROOT histogram.
void prepare() override
Define histograms and read payloads from DB.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.