Skip to content

Commit

Permalink
Merge pull request root-project#4 from zzxuanyuan/rrootbulkds
Browse files Browse the repository at this point in the history
Fix an error that an user_buffer never gets assigned to a basket.
  • Loading branch information
bbockelm authored Feb 6, 2019
2 parents b1003d5 + df8665a commit c1f9b6b
Show file tree
Hide file tree
Showing 3 changed files with 151 additions and 8 deletions.
153 changes: 147 additions & 6 deletions tree/dataframe/test/datasource_rootbulk.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -27,12 +27,25 @@ const auto kSlots = 4U;

TEST(RRootBulkDS, GenerateData)
{
int i = 0;
Char_t c = 'a';
Int_t i = 0;
UInt_t ui = 0;
Long64_t l = 0;
ULong64_t ul = 0;
Float_t f = 0.0;
Double_t d = 0.0;

// TODO: re-introduce multiple files so we can test globbing.
for (auto &&fileName : {fileName0}) {
RDataFrame tdf(eventCount);
tdf.Define("i", [&i]() { return i++; })
.Snapshot<int32_t>(treeName, fileName, {"i"});
.Define("ui", [&ui]() { return ui++; })
.Define("l", [&l]() { return l++; })
.Define("ul", [&ul]() { return ul++; })
.Define("f", [&f]() { return f++; })
.Define("d", [&d]() { return d++; })
.Define("c", [&c]() { return c++; })
.Snapshot(treeName, fileName);
}

// Manually create a larger file, allowing us to closely control its formatting.
Expand All @@ -43,11 +56,29 @@ TEST(RRootBulkDS, GenerateData)
auto tree = new TTree("t", "A ROOT tree of floats.");
tree->SetAutoFlush(bigFlushSize);
tree->SetBit(TTree::kOnlyFlushAtCluster);
int f = 2;
tree->Branch("myInt", &f, 320000, 1);
Int_t i2 = 2;
UInt_t ui2 = 2;
Long64_t l2 = 2;
ULong64_t ul2 = 2;
Float_t f2 = 2.0;
Double_t d2 = 2.0;
Char_t c2 = 'a';
tree->Branch("myInt", &i2, 320000, 1);
tree->Branch("myUInt", &ui2, 320000, 1);
tree->Branch("myLLong", &l2, 320000, 1);
tree->Branch("myULLong", &ul2, 320000, 1);
tree->Branch("myFloat", &f2, 320000, 1);
tree->Branch("myDouble", &d2, 320000, 1);
tree->Branch("myChar", &c2, 320000, 1);
for (Long64_t ev = 0; ev < bigEventCount; ev++) {
tree->Fill();
f ++;
i2++;
ui2++;
l2++;
ul2++;
f2++;
d2++;
c2++;
}
hfile = tree->GetCurrentFile();
hfile->Write();
Expand All @@ -57,7 +88,6 @@ TEST(RRootBulkDS, GenerateData)
delete hfile;
}


// Test that the bulk APIs can actually read the data generated by RDF.
TEST(RRootBulkDS, DirectBulkRead)
{
Expand Down Expand Up @@ -368,6 +398,117 @@ TEST(RRootBulkDS, BenchmarkBulkDS)
printf("Using bulk RDS directly: Total elapsed time (seconds): %.2f\n", sw.RealTime());
}

// Test different primitives.

TEST(RRootBulkDS, Primitives)
{
auto hfile = TFile::Open(fileName0);
TBufferFile branchbufI(TBuffer::kWrite, 32*1024);
TBufferFile branchbufUI(TBuffer::kWrite, 32*1024);
TBufferFile branchbufL(TBuffer::kWrite, 32*1024);
TBufferFile branchbufUL(TBuffer::kWrite, 32*1024);
TBufferFile branchbufF(TBuffer::kWrite, 32*1024);
TBufferFile branchbufD(TBuffer::kWrite, 32*1024);
TBufferFile branchbufC(TBuffer::kWrite, 32*1024);
TTree *tree = dynamic_cast<TTree*>(hfile->Get("t"));
ASSERT_TRUE(tree);

TBranch *branchI = tree->GetBranch("i");
ASSERT_TRUE(branchI);
TBranch *branchUI = tree->GetBranch("ui");
ASSERT_TRUE(branchUI);
TBranch *branchL = tree->GetBranch("l");
ASSERT_TRUE(branchL);
TBranch *branchUL = tree->GetBranch("ul");
ASSERT_TRUE(branchUL);
TBranch *branchF = tree->GetBranch("f");
ASSERT_TRUE(branchF);
TBranch *branchD = tree->GetBranch("d");
ASSERT_TRUE(branchD);
TBranch *branchC = tree->GetBranch("c");
ASSERT_TRUE(branchC);

Int_t events = eventCount;
Long64_t evt_idx = 0;
Int_t expectedI = 0;
UInt_t expectedUI = 0;
Long64_t expectedL = 0;
ULong64_t expectedUL = 0;
Float_t expectedF = 0.0;
Double_t expectedD = 0.0;
Char_t expectedC = 'a';
while (events) {
auto countI = branchI->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufI);
ASSERT_EQ(countI, eventCount);
auto countUI = branchUI->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufUI);
ASSERT_EQ(countUI, eventCount);
auto countL = branchL->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufL);
ASSERT_EQ(countL, eventCount);
auto countUL = branchUL->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufUL);
ASSERT_EQ(countUL, eventCount);
auto countF = branchF->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufF);
ASSERT_EQ(countF, eventCount);
auto countD = branchD->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufD);
ASSERT_EQ(countD, eventCount);
auto countC = branchC->GetBulkRead().GetEntriesSerialized(evt_idx, branchbufC);
ASSERT_EQ(countC, eventCount);
events = events > eventCount ? (events - eventCount) : 0;

int *entryI = reinterpret_cast<int *>(branchbufI.GetCurrent());
unsigned int *entryUI = reinterpret_cast<unsigned int *>(branchbufUI.GetCurrent());
long long int *entryL = reinterpret_cast<long long int *>(branchbufL.GetCurrent());
unsigned long long int *entryUL = reinterpret_cast<unsigned long long int *>(branchbufUL.GetCurrent());
float *entryF = reinterpret_cast<float *>(branchbufF.GetCurrent());
double *entryD = reinterpret_cast<double *>(branchbufD.GetCurrent());
char *entryC = reinterpret_cast<char *>(branchbufC.GetCurrent());
for (Int_t idx=0; idx < eventCount; idx++) {
Int_t tmpI = *reinterpret_cast<Int_t*>(&entryI[idx]);
UInt_t tmpUI = *reinterpret_cast<UInt_t*>(&entryUI[idx]);
Long64_t tmpL = *reinterpret_cast<Long64_t*>(&entryL[idx]);
ULong64_t tmpUL = *reinterpret_cast<ULong64_t*>(&entryUL[idx]);
float_t tmpF = *reinterpret_cast<Float_t*>(&entryF[idx]);
Double_t tmpD = *reinterpret_cast<Double_t*>(&entryD[idx]);
Char_t tmpC = *reinterpret_cast<Char_t*>(&entryC[idx]);
char *tmp_ptrI = reinterpret_cast<char *>(&tmpI);
char *tmp_ptrUI = reinterpret_cast<char *>(&tmpUI);
char *tmp_ptrL = reinterpret_cast<char *>(&tmpL);
char *tmp_ptrUL = reinterpret_cast<char *>(&tmpUL);
char *tmp_ptrF = reinterpret_cast<char *>(&tmpF);
char *tmp_ptrD = reinterpret_cast<char *>(&tmpD);
char *tmp_ptrC = reinterpret_cast<char *>(&tmpC);
int valI;
unsigned int valUI;
long long int valL;
unsigned long long int valUL;
float valF;
double valD;
char valC;
frombuf(tmp_ptrI, &valI);
EXPECT_EQ(valI, expectedI);
frombuf(tmp_ptrUI, &valUI);
EXPECT_EQ(valUI, expectedUI);
frombuf(tmp_ptrL, &valL);
EXPECT_EQ(valL, expectedL);
frombuf(tmp_ptrUL, &valUL);
EXPECT_EQ(valUL, expectedUL);
frombuf(tmp_ptrF, &valF);
EXPECT_EQ(valF, expectedF);
frombuf(tmp_ptrD, &valD);
EXPECT_EQ(valD, expectedD);
frombuf(tmp_ptrC, &valC);
EXPECT_EQ(valC, expectedC);
expectedI++;
expectedUI++;
expectedL++;
expectedUL++;
expectedF++;
expectedD++;
expectedC++;
evt_idx++;
}
}
}

#endif // R__USE_IMT

#endif // R__B64
1 change: 1 addition & 0 deletions tree/tree/inc/TLeafL.h
Original file line number Diff line number Diff line change
Expand Up @@ -39,6 +39,7 @@ class TLeafL : public TLeaf {

virtual void Export(TClonesArray *list, Int_t n);
virtual void FillBasket(TBuffer &b);
virtual DeserializeType GetDeserializeType() const { return kInPlace; }
const char *GetTypeName() const;
virtual Int_t GetMaximum() const { return (Int_t)fMaximum; }
virtual Int_t GetMinimum() const { return (Int_t)fMinimum; }
Expand Down
5 changes: 3 additions & 2 deletions tree/tree/src/TBranch.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -1203,6 +1203,7 @@ TBasket* TBranch::GetBasketImpl(Int_t basketnumber, TBuffer* user_buffer)
// reference to an existing basket in memory ?
if (basketnumber <0 || basketnumber > fWriteBasket) return 0;
TBasket *basket = (TBasket*)fBaskets.UncheckedAt(basketnumber);
if (user_buffer) basket = nullptr;
if (basket) return basket;
if (basketnumber == fWriteBasket) return 0;

Expand Down Expand Up @@ -1320,7 +1321,7 @@ Int_t TBranch::GetBasketAndFirst(TBasket*&basket, Long64_t &first,
{
Long64_t updatedNext = fNextBasketEntry;
Long64_t entry = fReadEntry;
if (R__likely(fCurrentBasket && fFirstBasketEntry <= entry && entry < fNextBasketEntry)) {
if (R__likely(fCurrentBasket && fFirstBasketEntry <= entry && entry < fNextBasketEntry && !user_buffer)) {
// We have found the basket containing this entry.
// make sure basket buffers are in memory.
basket = fCurrentBasket;
Expand Down Expand Up @@ -1350,7 +1351,7 @@ Int_t TBranch::GetBasketAndFirst(TBasket*&basket, Long64_t &first,
// We have found the basket containing this entry.
// make sure basket buffers are in memory.
basket = (TBasket*) fBaskets.UncheckedAt(fReadBasket);
if (!basket) {
if (!basket || user_buffer) {
basket = GetBasketImpl(fReadBasket, user_buffer);
if (!basket) {
fCurrentBasket = 0;
Expand Down

0 comments on commit c1f9b6b

Please sign in to comment.