Skip to content

[io] don't miss writing a histogram that is only in a last file with option -n 2 #18679

New issue

Have a question about this project? # for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “#”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? # to your account

Open
wants to merge 5 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
6 changes: 4 additions & 2 deletions io/io/src/TFileMerger.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -542,7 +542,8 @@ Bool_t TFileMerger::MergeOne(TDirectory *target, TList *sourcelist, Int_t type,
keyname, keytitle);
return kTRUE;
}
Bool_t canBeFound = (type & kIncremental) && (current_sourcedir->GetList()->FindObject(keyname) != nullptr);
Bool_t canBeFound = (type & kIncremental) && (current_sourcedir->GetList()->FindObject(keyname) != nullptr) &&
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why is it && and not || ?

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Or more exactly, I don't understand yet why the fact the histogram can be found means that it is not written in the end ...

Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Check out:

https://github.com/ferdymercury/root/blob/1531153ea11a7b54f4eb3c170bbd28e9bc46447f/io/io/src/TFileMerger.cxx#L808-L817

so if canBeFound is true, there is an optimization that spare some write cycles.

We use && to avoid it being true, ie to force writing to file. Using || would go in the different direction.

Do you want me to rename canBeFound to skipPartialWriting ?

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Do you want me to rename canBeFound to skipPartialWriting ?

Not yet. I am still confused.

The original canBeFound meant (in the context of incremental merge) 'histogram can be found in the source' while the new version is histogram can be found in the source and in the target.

The optimization is indeed 'skip partial writing if we can found the histogram again'.

So I don't understand (yet) the semantic of the change. ie. Why is the new criteria the right choice? Is the new criteria instead 'just' making canBeFound always false?

Another avenue of inquiry is 'the original code assume that if canBeFound is true then there will be another change to write the histogram. Why is it no true anymore (it does not seem to be realted to 'can not be found in target')? Is there other variation of the example that also fails (how does the -n X value relates to the number of files in the input list and how many are 'missing' the histograms).

Related: is it possible that the alternative if that at the refresh boundary there needs to be a flush/write as if we were at the end?

Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks for the review. I am not sure about these questions; I just followed jblomer's suggestion. My (limited) understanding is that this change just forces an extra partial write the first time that a new histogram appear in any file. So it does not really harm, but is suboptimal since, if all files have exactly the same histograms, then we could have waited until the end. But it makes it work if there are some files with and some without, independently of the chosen N.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

this change just forces an extra partial write the first time that a new histogram appear in any file.

That would be fair enough. Can we verify that it just one and not number_of_input files? (Related but probably unavoidable is that for the case of just 2 files with all the same histograms the number of write is doubled .... actually the number of double for each histogram that is in more than one file (i.e. usual case)).

Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

See here the test with many more files mtest.cpp.txt.

root -l mtest.cpp.txt+ -b -q 2>&1 | grep MergeOne
Info in TFileMerger::MergeOne: Writing partial result of h1 into target
Info in TFileMerger::MergeOne: Writing partial result of h2 into target

(target->GetList()->FindObject(keyname) != nullptr);

// if (cl->IsTObject())
// obj->ResetBit(kMustCleanup);
Expand Down Expand Up @@ -806,7 +807,8 @@ Bool_t TFileMerger::MergeOne(TDirectory *target, TList *sourcelist, Int_t type,
delete ndir;
}
} else if (!canBeFound) { // Don't write the partial result for TTree and TH1

if (gDebug > 0)
Info("MergeOne", "Writing partial result of %s into target", oldkeyname.Data());
if (!canBeMerged) {
TIter peeknextkey(nextkey);
status = WriteCycleInOrder(oldkeyname, nextkey, peeknextkey, target) && status;
Expand Down
36 changes: 36 additions & 0 deletions io/io/test/TFileMergerTests.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -219,3 +219,39 @@ TEST(TFileMerger, ChangeFile)
gSystem->Unlink("file6640mergerinput_2.root");
gSystem->Unlink("file6640mergeroutput.root");
}

// https://github.com/root-project/root/issues/9022
TEST(TFileMerger, SingleHistFile)
{
auto filename1 = "f1_9022.root";
auto filename2 = "f2_9022.root";
auto outname = "file9022mergeroutput.root";
{
TFile f1(filename1, "RECREATE");
TH1F h("h1", "h1", 1, 0, 1);
h.Write();
f1.Close();
TFile f2(filename2, "RECREATE");
TH1F h2("h2", "h2", 1, 0, 1);
h2.Write();
f2.Close();
}
{
TFileMerger filemerger{false, false};
filemerger.SetMaxOpenedFiles(2);
filemerger.OutputFile(std::unique_ptr<TFile>{TFile::Open(outname, "RECREATE")});

filemerger.AddFile(filename1);
filemerger.AddFile(filename2);

filemerger.Merge();
}
{
TFile file(outname, "READ");
EXPECT_NE(file.Get<TH1>("h1"), nullptr);
EXPECT_NE(file.Get<TH1>("h2"), nullptr);
}
gSystem->Unlink(filename1);
gSystem->Unlink(filename2);
gSystem->Unlink(outname);
}
Loading