Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DD4hepBinningHelpers.hpp
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file DD4hepBinningHelpers.hpp
1 // This file is part of the Acts project.
2 //
3 // Copyright (C) 2023 CERN for the benefit of the Acts project
4 //
5 // This Source Code Form is subject to the terms of the Mozilla Public
6 // License, v. 2.0. If a copy of the MPL was not distributed with this
7 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
8 
9 #pragma once
10 
16 
17 #include <string>
18 #include <tuple>
19 #include <vector>
20 
21 #include <DD4hep/DD4hepUnits.h>
22 #include <DD4hep/DetElement.h>
23 #include <DD4hep/DetFactoryHelper.h>
24 #include <DD4hep/Objects.h>
25 #include <DDRec/DetectorData.h>
26 #include <XML/Utilities.h>
27 
28 namespace Acts {
29 
30 static std::vector<std::tuple<std::string, BinningValue>> allowedBinnings = {
31  {"x", binX}, {"y", binY}, {"z", binZ}, {"phi", binPhi}, {"r", binR}};
32 
57 void decodeBinning(dd4hep::rec::VariantParameters &variantParams,
58  const xml_comp_t &xmlBinning, const std::string &bname,
59  const std::vector<std::string> &bvals) {
60  // Set the surface binninng parameter to true
61  variantParams.set<int>(std::string(bname + "_dim"), bvals.size());
62  for (const auto &bv : bvals) {
63  // Gather the number of bins, 0 indicates variable binning
64  int nBins = Acts::getAttrValueOr<int>(xmlBinning, std::string("n" + bv), 0);
65  // Gather the bin expansion parameter, expansion of 0 is default
66  int nExpansion =
67  Acts::getAttrValueOr<int>(xmlBinning, std::string(bv + "expansion"), 0);
68  variantParams.set<int>(bname + "_" + bv + "_exp", nExpansion);
69  // Equidistant binning detected
70  if (nBins > 0) {
71  // Set the type identificatio
72  variantParams.set<std::string>(bname + "_" + bv + "_type", "equidistant");
73  // Set the number of bins
74  variantParams.set<int>(bname + "_" + bv + "_n", nBins);
75  // Set min/max paraeter
76  variantParams.set<double>(
77  bname + "_" + bv + "_min",
78  xmlBinning.attr<double>(std::string(bv + "min").c_str()));
79  variantParams.set<double>(
80  bname + "_" + bv + "_max",
81  xmlBinning.attr<double>(std::string(bv + "max").c_str()));
82  } else {
83  // Variable binning detected
84  variantParams.set<std::string>(bname + "_" + bv + "_type", "variable");
85  // Get the number of bins explicitly
86  auto boundaries =
87  xmlBinning.attr<std::string>(std::string(bv + "boundaries").c_str());
88  std::string del = ",";
89  int end = boundaries.find(del);
90  int ib = 0;
91  // Unit conversion
92  double unitScalar = 1.;
93  if (bv != "phi") {
94  unitScalar = Acts::UnitConstants::mm / dd4hep::millimeter;
95  }
96  // Split and convert
97  while (end != -1) {
98  double bR = unitScalar * dd4hep::_toFloat(boundaries.substr(0, end));
99  variantParams.set<double>(
100  bname + "_" + bv + "_b" + std::to_string(ib++), bR);
101  boundaries.erase(boundaries.begin(), boundaries.begin() + end + 1);
102  end = boundaries.find(del);
103  }
104  double bR = unitScalar * std::stod(boundaries.substr(0, end));
105  variantParams.set<double>(bname + "_" + bv + "_b" + std::to_string(ib),
106  bR);
107  // The number of bins are needed to unpack the data
108  variantParams.set<int>(bname + "_" + bv + "_n", ib);
109  }
110  }
111 }
112 
119 std::vector<Acts::Experimental::ProtoBinning> convertBinning(
120  const dd4hep::DetElement &dd4hepElement, const std::string &bname) {
121  std::vector<Experimental::ProtoBinning> protoBinnings;
122  for (const auto &[ab, bVal] : allowedBinnings) {
123  auto type =
124  getParamOr<std::string>(bname + "_" + ab + "_type", dd4hepElement, "");
125  if (not type.empty()) {
126  // Default binning is bound
127  auto bType = Acts::detail::AxisBoundaryType::Bound;
128  // Equidistant or variable binning
129  detail::AxisType aType = type == "equidistant"
130  ? detail::AxisType::Equidistant
131  : detail::AxisType::Variable;
132  int nBins = getParamOr<int>(bname + "_" + ab + "_n", dd4hepElement, 0);
133  int nExpansion =
134  getParamOr<int>(bname + "_" + ab + "_exp", dd4hepElement, 0);
135  if (aType == detail::AxisType::Equidistant) {
136  // Equidistant binning
137  auto min = getParamOr<ActsScalar>(bname + "_" + ab + "_min",
138  dd4hepElement, 0.);
139  auto max = getParamOr<ActsScalar>(bname + "_" + ab + "_max",
140  dd4hepElement, 0.);
141  // Check for closed phi binning
142  if (bVal == binPhi and (max - min) > 1.9 * M_PI) {
143  bType = Acts::detail::AxisBoundaryType::Closed;
144  }
145  protoBinnings.push_back(Experimental::ProtoBinning(
146  bVal, bType, min, max, nBins, nExpansion));
147  } else {
148  // Variable binning
149  std::vector<ActsScalar> edges;
150  for (int ib = 0; ib <= nBins; ++ib) {
151  edges.push_back(getParamOr<ActsScalar>(
152  bname + "_" + ab + "_b" + std::to_string(ib), dd4hepElement, 0.));
153  }
154  // Check for closed phi binning
155  if (bVal == binPhi and (edges.back() - edges.front()) > 1.9 * M_PI) {
156  bType = Acts::detail::AxisBoundaryType::Closed;
157  }
158  protoBinnings.push_back(
159  Experimental::ProtoBinning(bVal, bType, edges, nExpansion));
160  }
161  }
162  }
163  return protoBinnings;
164 }
165 
166 } // namespace Acts