Belle II Software  release-05-02-19
SVDTimeValidationAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Luigi Corona, Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #include <svd/calibration/SVDTimeValidationAlgorithm.h>
11 
12 #include <svd/dbobjects/SVDCoGCalibrationFunction.h>
13 #include <svd/calibration/SVDCoGTimeCalibrations.h>
14 
15 #include <TF1.h>
16 #include <TProfile.h>
17 #include <TH2F.h>
18 #include <framework/logging/Logger.h>
19 #include <iostream>
20 #include <TString.h>
21 #include <TFitResult.h>
22 
23 using namespace std;
24 using namespace Belle2;
25 
26 SVDTimeValidationAlgorithm::SVDTimeValidationAlgorithm() :
27  CalibrationAlgorithm("SVDTimeValidationCollector")
28 {
29  setDescription("SVDTimeValidation calibration algorithm");
30 }
31 
33 {
34 
35  gROOT->SetBatch(true);
36 
37  auto hEventT0 = getObjectPtr<TH1F>("hEventT0");
38  float eventT0_mean = 0;
39  float eventT0_rms = 0;
40  if (hEventT0) {
41  eventT0_mean = hEventT0->GetMean();
42  eventT0_rms = hEventT0->GetRMS();
43  } else
44  B2ERROR("Histogram with Event T0 not found");
45 
46  B2DEBUG(27, "Histogram: " << hEventT0->GetName() <<
47  " Entries (n. clusters): " << hEventT0->GetEntries() <<
48  " Mean: " << eventT0_mean);
49 
51  for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
52  for (auto ladder : geoCache.getLadders(layer)) {
53  for (Belle2::VxdID sensor : geoCache.getSensors(ladder)) {
54  for (int view = SVDHistograms<TH1F>::VIndex ; view < SVDHistograms<TH1F>::UIndex + 1; view++) {
55  char side = 'U';
56  if (view == 0)
57  side = 'V';
58  auto layer_num = sensor.getLayerNumber();
59  auto ladder_num = sensor.getLadderNumber();
60  auto sensor_num = sensor.getSensorNumber();
61  auto hClsTimeOnTracks = getObjectPtr<TH1F>(Form("clsTimeOnTracks__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
62  float clsTimeOnTracks_mean = 0.;
63  if (hClsTimeOnTracks)
64  clsTimeOnTracks_mean = hClsTimeOnTracks->GetMean();
65  else
66  B2ERROR("Histogram " << Form("clsTimeOnTracks__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side) << " not found");
67  auto deviation = (clsTimeOnTracks_mean - eventT0_mean) / eventT0_rms;
68 
69  B2DEBUG(27, "Histogram: " << hClsTimeOnTracks->GetName() <<
70  " Entries (n. clusters): " << hClsTimeOnTracks->GetEntries() <<
71  " Mean: " << clsTimeOnTracks_mean <<
72  " Deviation: " << deviation << " EventT0 RMS");
73  if (abs(deviation) > m_allowedDeviationMean)
74  B2ERROR("Histogram: " << hClsTimeOnTracks->GetName() << " deviates from EventT0 by" << deviation << " times the EventT0 RMS");
75 
76  }
77  }
78  }
79  }
80  return c_OK;
81 }
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::SVDHistograms
template class for SVd histograms
Definition: SVDHistograms.h:36
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::SVDTimeValidationAlgorithm::m_allowedDeviationMean
float m_allowedDeviationMean
Allowed deviation of clsOnTracks histo wrt EventT0 histo in number of EventT0 RMS.
Definition: SVDTimeValidationAlgorithm.h:68
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::SVDTimeValidationAlgorithm::calibrate
virtual EResult calibrate() override
Run algo on data.
Definition: SVDTimeValidationAlgorithm.cc:32
Belle2::VXD::GeoCache::getInstance
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:215
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::VXD::SensorInfoBase::SVD
@ SVD
SVD Sensor.
Definition: SensorInfoBase.h:45
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::VXD::GeoCache
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:41
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47