Skip to content

[tree] return NaN instead of zero when out of bonds in TreeFormula #18766

New issue

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

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

Already on GitHub? Sign in to your account

Merged
merged 3 commits into from
May 17, 2025
Merged
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
15 changes: 12 additions & 3 deletions roottest/root/treeformula/array/execOutOfBounds.C
Original file line number Diff line number Diff line change
Expand Up @@ -4,11 +4,20 @@ int execOutOfBounds() {
{
#endif
TFile *_file0 = TFile::Open("run133.root");
TTree *t; _file0->GetObject("t",t);
auto t = _file0->Get<TTree>("t");
// ThetaDeg is a dynamic array with indices going from 0 to (ringtouche_DE1 - 1)
TTreeFormula tf("tf", "ThetaDeg", t);
t->GetEntry(0); // ringtouche for entry 0 is 0
auto res = tf.EvalInstance(1); // ThetaDeg[1] goes out of bonds
if (!TMath::IsNaN(res)) {
printf("Error: evaluated instance is %f rather than NaN\n", res);
return 1;
}
TH2F *h2 = new TH2F("h2","h2",5,0,5,100,0,100);
t->Draw("ThetaDeg[ringtouche_DE1]:ringtouche_DE1>>h2","","colz goff");
if ( (int)(h2->GetMean(1)*1000) != 3) {
printf("Error: x mean is %f rather than 0.003\n",h2->GetMean(1));
// Filling NaNs increases entry number but does not change any bin content
if ( (int)(h2->GetMean(1)*1000) != 0) {
printf("Error: x mean is %f rather than 0.000\n",h2->GetMean(1));
return 1;
}
if ( (int)(h2->GetMean(2)*10) != 0) {
Expand Down
77 changes: 69 additions & 8 deletions tree/treeplayer/src/TTreeFormula.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -3759,7 +3759,7 @@ const char* TTreeFormula::EvalStringInstance(Int_t instance)
const Int_t real_instance = GetRealInstance(instance,0); \
\
if (instance==0) fNeedLoading = true; \
if (real_instance>=fNdata[0]) return 0; \
if (real_instance>=fNdata[0]) return TMath::SignalingNaN(); \
\
/* Since the only operation in this formula is reading this branch, \
we are guaranteed that this function is first called with instance==0 and \
Expand Down Expand Up @@ -3793,7 +3793,7 @@ const char* TTreeFormula::EvalStringInstance(Int_t instance)
#define TREE_EVAL_INIT \
const Int_t real_instance = GetRealInstance(instance,0); \
\
if (real_instance>=fNdata[0]) return 0; \
if (real_instance>=fNdata[0]) return TMath::SignalingNaN(); \
\
if (fAxis) { \
char * label; \
Expand Down Expand Up @@ -3845,13 +3845,13 @@ const char* TTreeFormula::EvalStringInstance(Int_t instance)
} \
} \
} \
if (real_instance>=fNdata[code]) return 0;
if (real_instance>=fNdata[code]) return TMath::SignalingNaN();

#define TREE_EVAL_INIT_LOOP \
/* Now let calculate what physical instance we really need. */ \
const Int_t real_instance = GetRealInstance(instance,code); \
\
if (real_instance>=fNdata[code]) return 0;
if (real_instance>=fNdata[code]) return TMath::SignalingNaN();


template<typename T> T Summing(TTreeFormula *sum) {
Expand Down Expand Up @@ -3993,14 +3993,75 @@ template<> inline Long64_t TTreeFormula::GetConstant(Int_t k) { return (Long64_t
/// \tparam T The type used to interpret the numbers then used for the operations
/// \param instance iteration instance
/// \param stringStackArg formula as string
/// \return the result of the evaluation
/// \return the result of the evaluation, or a signaling NaN if out of bounds
///
/// \warning Care has to be taken before calling this function with std::vector
/// or dynamically sized objects, rather than plain fixed-size arrays.
/// For example, this works without problems:
/// ~~~{.cpp}
/// TTree t("t", "t");
/// Float_t x[2]{};
/// t.Branch("xa", &x, "x[2]/F");
/// x[1] = 1;
/// t.Fill();
/// x[1] = 2;
/// t.Fill();
/// t.Scan();
/// TTreeFormula tfx("tfx", "xa[1]", &t);
/// t.GetEntry(0);
/// tfx.EvalInstance()
/// t.GetEntry(1);
/// tfx.EvalInstance()
/// ~~~
/// But the following fails (independently on whether the size changed or not between entries):
/// ~~~{.cpp}
/// TTree t("t", "t");
/// vector<Short_t> v;
/// t.Branch("vec", &v);
/// v.push_back(2);
/// v.push_back(3);
/// t.Fill();
/// v.clear();
/// v.push_back(4);
/// v.push_back(5);
/// t.Fill();
/// t.Scan();
/// TTreeFormula tfv1("tfv1", "vec[1]", &t);
/// TTreeFormula tfv("tfv", "vec", &t);
/// t.GetEntry(0);
/// tfv1.EvalInstance()
/// tfv.EvalInstance(1)
/// t.GetEntry(1);
/// tfv1.EvalInstance()
/// tfv.EvalInstance(1)
/// ~~~
/// To prevent this, when working with objects with dynamic size for each entry, one needs
/// to mimick what TTree::Scan does, i.e. to check the value of
/// `GetNdata()` before calling `EvalInstance()`:
/// ~~~{.cpp}
/// t.GetEntry(0);
/// if (tfv1.GetNdata() > 0)
/// tfv1.EvalInstance()
/// if (tfv.GetNdata() > 1)
/// tfv.EvalInstance(1)
/// t.GetEntry(1);
/// if (tfv1.GetNdata() > 0)
/// tfv1.EvalInstance()
/// if (tfv.GetNdata() > 1)
/// tfv.EvalInstance(1)
/// ~~~
/// Note that for `tfv1`, even if the index is fixed in the formula and even if each entry
/// had the same std::vector size, since the formula contains a branch with theoretically variable size,
/// one must check GetNData() as there might 0 or 1 answers. Since even with fixed index,
/// the collection might be too small to fulfill it.
/// TTreeFormula::GetMultiplicity tells you (indirectly) whether you need to call GetNData or not for a given formula.

template<typename T>
T TTreeFormula::EvalInstance(Int_t instance, const char *stringStackArg[])
{
// Note that the redundancy and structure in this code is tailored to improve
// efficiencies.
if (TestBit(kMissingLeaf)) return 0;
if (TestBit(kMissingLeaf)) return TMath::SignalingNaN();
if (fNoper == 1 && fNcodes > 0) {

switch (fLookupType[0]) {
Expand Down Expand Up @@ -4059,7 +4120,7 @@ T TTreeFormula::EvalInstance(Int_t instance, const char *stringStackArg[])
}
return fx->EvalInstance<T>(instance);
}
default: return 0;
default: return TMath::SignalingNaN();
}
}

Expand Down Expand Up @@ -4443,7 +4504,7 @@ T TTreeFormula::EvalInstance(Int_t instance, const char *stringStackArg[])
Long64_t treeEntry = br->GetTree()->GetReadEntry();
R__LoadBranch(br,treeEntry,true);
}
if (real_instance>=fNdata[string_code]) return 0;
if (real_instance>=fNdata[string_code]) return TMath::SignalingNaN();
}
pos2++;
if (fLookupType[string_code]==kDirect) {
Expand Down
Loading