Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:23:12

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #include "Acts/Material/AccumulatedSurfaceMaterial.hpp"
0010 
0011 #include "Acts/Material/BinnedSurfaceMaterial.hpp"
0012 #include "Acts/Material/HomogeneousSurfaceMaterial.hpp"
0013 #include "Acts/Utilities/ProtoAxisHelpers.hpp"
0014 
0015 #include <utility>
0016 
0017 // Default Constructor - for homogeneous material
0018 Acts::AccumulatedSurfaceMaterial::AccumulatedSurfaceMaterial(double splitFactor)
0019     : m_splitFactor(splitFactor) {
0020   AccumulatedVector accMat = {{AccumulatedMaterialSlab()}};
0021   m_accumulatedMaterial = {{accMat}};
0022 }
0023 
0024 // Binned Material constructor with split factor
0025 Acts::AccumulatedSurfaceMaterial::AccumulatedSurfaceMaterial(
0026     const BinUtility& binUtility, double splitFactor)
0027     : m_binUtility(binUtility), m_splitFactor(splitFactor) {
0028   std::size_t bins0 = m_binUtility.bins(0);
0029   std::size_t bins1 = m_binUtility.bins(1);
0030   AccumulatedVector accVec(bins0, AccumulatedMaterialSlab());
0031   m_accumulatedMaterial = AccumulatedMatrix(bins1, accVec);
0032 }
0033 
0034 // Assign a material properties object
0035 std::array<std::size_t, 3> Acts::AccumulatedSurfaceMaterial::accumulate(
0036     const Vector2& lp, const MaterialSlab& mp, double pathCorrection) {
0037   if (m_binUtility.dimensions() == 0) {
0038     m_accumulatedMaterial[0][0].accumulate(mp, pathCorrection);
0039     return {0, 0, 0};
0040   }
0041   std::size_t bin0 = m_binUtility.bin(lp, 0);
0042   std::size_t bin1 = m_binUtility.bin(lp, 1);
0043   m_accumulatedMaterial[bin1][bin0].accumulate(mp, pathCorrection);
0044   return {bin0, bin1, 0};
0045 }
0046 
0047 // Assign a material properties object
0048 std::array<std::size_t, 3> Acts::AccumulatedSurfaceMaterial::accumulate(
0049     const Vector3& gp, const MaterialSlab& mp, double pathCorrection) {
0050   if (m_binUtility.dimensions() == 0) {
0051     m_accumulatedMaterial[0][0].accumulate(mp, pathCorrection);
0052     return {0, 0, 0};
0053   }
0054   std::array<std::size_t, 3> bTriple = m_binUtility.binTriple(gp);
0055   m_accumulatedMaterial[bTriple[1]][bTriple[0]].accumulate(mp, pathCorrection);
0056   return bTriple;
0057 }
0058 
0059 // Void average for vacuum assignment
0060 void Acts::AccumulatedSurfaceMaterial::trackVariance(const Vector3& gp,
0061                                                      MaterialSlab slabReference,
0062                                                      bool emptyHit) {
0063   if (m_binUtility.dimensions() == 0) {
0064     m_accumulatedMaterial[0][0].trackVariance(slabReference, emptyHit);
0065     return;
0066   }
0067   std::array<std::size_t, 3> bTriple = m_binUtility.binTriple(gp);
0068   std::vector<std::array<std::size_t, 3>> trackBins = {bTriple};
0069   trackVariance(trackBins, slabReference);
0070 }
0071 
0072 // Average the information accumulated during one event
0073 void Acts::AccumulatedSurfaceMaterial::trackVariance(
0074     const std::vector<std::array<std::size_t, 3>>& trackBins,
0075     MaterialSlab slabReference, bool emptyHit) {
0076   // the homogeneous material case
0077   if (m_binUtility.dimensions() == 0) {
0078     m_accumulatedMaterial[0][0].trackVariance(slabReference, emptyHit);
0079     return;
0080   }
0081   // The touched bins are known, so you can access them directly
0082   if (!trackBins.empty()) {
0083     for (auto bin : trackBins) {
0084       m_accumulatedMaterial[bin[1]][bin[0]].trackVariance(slabReference);
0085     }
0086   } else {
0087     // Touched bins are not known: Run over all bins
0088     for (auto& matVec : m_accumulatedMaterial) {
0089       for (auto& mat : matVec) {
0090         mat.trackVariance(slabReference);
0091       }
0092     }
0093   }
0094 }
0095 
0096 // Void average for vacuum assignment
0097 void Acts::AccumulatedSurfaceMaterial::trackAverage(const Vector3& gp,
0098                                                     bool emptyHit) {
0099   if (m_binUtility.dimensions() == 0) {
0100     m_accumulatedMaterial[0][0].trackAverage(emptyHit);
0101     return;
0102   }
0103 
0104   std::array<std::size_t, 3> bTriple = m_binUtility.binTriple(gp);
0105   std::vector<std::array<std::size_t, 3>> trackBins = {bTriple};
0106   trackAverage(trackBins, emptyHit);
0107 }
0108 
0109 // Average the information accumulated during one event
0110 void Acts::AccumulatedSurfaceMaterial::trackAverage(
0111     const std::vector<std::array<std::size_t, 3>>& trackBins, bool emptyHit) {
0112   // the homogeneous material case
0113   if (m_binUtility.dimensions() == 0) {
0114     m_accumulatedMaterial[0][0].trackAverage(emptyHit);
0115     return;
0116   }
0117 
0118   // The touched bins are known, so you can access them directly
0119   if (!trackBins.empty()) {
0120     for (auto bin : trackBins) {
0121       m_accumulatedMaterial[bin[1]][bin[0]].trackAverage(emptyHit);
0122     }
0123   } else {
0124     // Touched bins are not known: Run over all bins
0125     for (auto& matVec : m_accumulatedMaterial) {
0126       for (auto& mat : matVec) {
0127         mat.trackAverage(emptyHit);
0128       }
0129     }
0130   }
0131 }
0132 
0133 /// Total average creates SurfaceMaterial
0134 std::unique_ptr<const Acts::ISurfaceMaterial>
0135 Acts::AccumulatedSurfaceMaterial::totalAverage() {
0136   if (m_binUtility.bins() == 1) {
0137     // Return HomogeneousSurfaceMaterial
0138     return std::make_unique<HomogeneousSurfaceMaterial>(
0139         m_accumulatedMaterial[0][0].totalAverage().first, m_splitFactor);
0140   }
0141   // Create the properties matrix
0142   MaterialSlabMatrix mpMatrix(
0143       m_binUtility.bins(1),
0144       MaterialSlabVector(m_binUtility.bins(0), MaterialSlab::Nothing()));
0145   // Loop over and fill
0146   for (std::size_t ib1 = 0; ib1 < m_binUtility.bins(1); ++ib1) {
0147     for (std::size_t ib0 = 0; ib0 < m_binUtility.bins(0); ++ib0) {
0148       mpMatrix[ib1][ib0] = m_accumulatedMaterial[ib1][ib0].totalAverage().first;
0149     }
0150   }
0151   // Now return the BinnedSurfaceMaterial
0152   return std::make_unique<const BinnedSurfaceMaterial>(
0153       m_binUtility, std::move(mpMatrix), m_splitFactor);
0154 }