Skip to content

Commit 705fea9

Browse files
author
Petr Jacka
committed
[hist] Using explicit implementation of the THnSparse copy constructor
1 parent 6babed6 commit 705fea9

File tree

3 files changed

+25
-1
lines changed

3 files changed

+25
-1
lines changed

hist/hist/inc/THnSparse.h

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -79,6 +79,8 @@ class THnSparse: public THnBase {
7979
THnSparse(const char *name, const char *title, Int_t dim, const Int_t *nbins,
8080
const std::vector<std::vector<double>> &xbins, Int_t chunksize = 1024 * 16);
8181

82+
THnSparse(const THnSparse& other);
83+
8284
~THnSparse() override;
8385

8486
static THnSparse* CreateSparse(const char* name, const char* title,

hist/hist/src/THnSparse.cxx

Lines changed: 20 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -644,6 +644,26 @@ THnBase(name, title, dim, nbins, xbins),
644644
fBinContent.SetOwner();
645645
}
646646

647+
////////////////////////////////////////////////////////////////////////////////
648+
/// Construct a THnSparse as a copy of "other"
649+
650+
THnSparse::THnSparse(const THnSparse& other):
651+
THnBase(other), fChunkSize(other.fChunkSize), fFilledBins(other.fFilledBins),
652+
fBinContent(other.fBinContent), fBins(other.fBins), fBinsContinued(other.fBinsContinued),
653+
fCompactCoord(nullptr) {
654+
655+
fBinContent.SetOwner();
656+
657+
Int_t dim = other.GetNdimensions();
658+
auto nbins=new Int_t[dim];
659+
for (Int_t i=0; i<dim; i++)
660+
nbins[i]=other.GetAxis(i)->GetNbins();
661+
662+
fCompactCoord = new THnSparseCompactBinCoord(dim, nbins);
663+
delete[] nbins;
664+
}
665+
666+
647667
////////////////////////////////////////////////////////////////////////////////
648668
/// Destruct a THnSparse
649669

hist/hist/test/THn.cxx

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -29,13 +29,15 @@ TEST(THn, Constructors) {
2929
THnD hn_v2("hn_v2", "hn_v2", 3, nbins.data(), edges);
3030
THnD hn_v3("hn_v3", "hn_v3", axes);
3131
THnI hn_v4("hn_v4", "hn_v4", axes);
32+
THnD hn_v5(hn_v1);
3233

3334
THnSparseD hs_v1("hs_v1", "hs_v1", 3, nbins.data(), xmin.data(), xmax.data());
3435
THnSparseD hs_v2("hs_v2", "hs_v2", 3, nbins.data(), edges);
3536
THnSparseD hs_v3("hs_v3", "hs_v3", axes);
3637
THnSparseI hs_v4("hs_v4", "hs_v4", axes);
38+
THnSparseD hs_v5(hs_v1);
3739

38-
std::vector<THnBase*> hns = {&hn_v1, &hn_v2, &hn_v3, &hn_v4, &hs_v1, &hs_v2, &hs_v3, &hs_v4};
40+
std::vector<THnBase*> hns = {&hn_v1, &hn_v2, &hn_v3, &hn_v4, &hn_v5, &hs_v1, &hs_v2, &hs_v3, &hs_v4, &hs_v5};
3941
for (THnBase* hn : hns) {
4042
EXPECT_EQ(hn->GetNdimensions(), 3);
4143
for (int dim = 0; dim < 3; ++dim) {

0 commit comments

Comments
 (0)