Belle II Software  release-05-02-19
SVDOccupancyCalibrationsAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributor: Laura Zani (2019) *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #include <svd/calibration/SVDOccupancyCalibrationsAlgorithm.h>
11 
12 #include <svd/calibration/SVDOccupancyCalibrations.h>
13 
14 #include <framework/logging/Logger.h>
15 
16 #include <TH1F.h>
17 #include <TString.h>
18 
19 #include <iostream>
20 
21 using namespace std;
22 using namespace Belle2;
23 
24 SVDOccupancyCalibrationsAlgorithm::SVDOccupancyCalibrationsAlgorithm(const std::string& str) :
25  CalibrationAlgorithm("SVDOccupancyCalibrationsCollector")
26 {
27  setDescription("Calibration algorithm for SVDOccupancyCalibrations payloads");
28  m_id = str;
29 }
30 
32 {
33 
34  float occCal = 1.;
35  auto payload = new Belle2::SVDOccupancyCalibrations::t_payload(occCal , m_id);
36 
37  auto tree = getObjectPtr<TTree>("HTreeOccupancyCalib");
38 
39  TH1F* hocc = new TH1F("", "", 768, 0, 768);
40 
41  if (!tree) {
42  B2WARNING("No tree object.");
43  } else if (!tree->GetEntries()) {
44  B2WARNING("No data in the tree.");
45  }
46 
47  int layer = 0;
48  int ladder = 0;
49  int sensor = 0;
50  int side = 0;
51 
52  tree->SetBranchAddress("hist", &hocc);
53  tree->SetBranchAddress("layer", &layer);
54  tree->SetBranchAddress("ladder", &ladder);
55  tree->SetBranchAddress("sensor", &sensor);
56  tree->SetBranchAddress("view", &side);
57 
58  for (int i = 0; i < tree->GetEntries(); i++) {
59  tree->GetEntry(i);
60 
61  int nstrips = 768;
62  if (!side && layer != 3) nstrips = 512;
63 
64  for (int iterStrip = 0; iterStrip < nstrips; iterStrip++) {
65  occCal = hocc->GetBinContent(iterStrip + 1);
66 
67  payload->set(layer, ladder, sensor, bool(side), iterStrip, occCal);
68  }
69  }
70 
71  saveCalibration(payload, "SVDOccupancyCalibrations");
72 
73  // probably not needed - would trigger re-doing the collection
74  //if ( ... too large corrections ... ) return c_Iterate;
75  return c_OK;
76 
77 }
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::SVDOccupancyCalibrationsAlgorithm::m_id
std::string m_id
identifier
Definition: SVDOccupancyCalibrationsAlgorithm.h:47
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::SVDOccupancyCalibrations::t_payload
SVDCalibrationsBase< SVDCalibrationsVector< float > > t_payload
typedef of the Occupancy payload of all SVD strips
Definition: SVDOccupancyCalibrations.h:45
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::SVDOccupancyCalibrationsAlgorithm::calibrate
virtual EResult calibrate() override
Run algo on data.
Definition: SVDOccupancyCalibrationsAlgorithm.cc:31