RE: [ROOT] TTree::Scan() and TTree::Draw() selection different?

From: Philippe Canal (pcanal@fnal.gov)
Date: Thu Sep 05 2002 - 22:53:07 MEST


Hi Hannes,

Please note that the result you are currently getting from h1.Draw are very
likely to be wrong.
   h1.Draw("Event", "(abs(Etac-Eta2dsk)>0.1) && (Eta2dsk <99) &&
(Ncsc==1)");
currently make the following comparison (due to a bug related to the fact
that there is only one element in Etac)
	(abs(Etac[0]-Eta2dsk[0])>0.1) && (Eta2dsk[0] <99) && (Ncsc==1)
and
	(abs(Etac[1]-Eta2dsk[1])>0.1) && (Eta2dsk[1] <99) && (Ncsc==1)
and since Etac[1] does not exist (sizeof(Etac)==Ncsc), you get a random
result.

If the bug was corrected (not trivial to do, so it will take a little bit
before it is fixed), you would then have gottent ONLY the following test:
	(abs(Etac[0]-Eta2dsk[0])>0.1) && (Eta2dsk[0] <99) && (Ncsc==1)
This is because in TTree::Draw when you use several arrays in the same
expression, all unspecified indexes are varied together.  This means that if
you used:
   h1.Draw("Event", "(abs(Etac-Eta2dsk)>0.1) && (Eta2dsk <99) &&
(Ncsc==2)");
the following test would have been done:
	(abs(Etac[0]-Eta2dsk[0])>0.1) && (Eta2dsk[0] <99) && (Ncsc==2)
	(abs(Etac[1]-Eta2dsk[1])>0.1) && (Eta2dsk[1] <99) && (Ncsc==2)

Anyway, since you specified that Etac has only one elements (Nscs==1), I
assume that you meant to have the comparaisons:
	(abs(Etac[0]-Eta2dsk[0])>0.1) && (Eta2dsk[0] <99) && (Ncsc==1)
and
	(abs(Etac[0]-Eta2dsk[1])>0.1) && (Eta2dsk[1] <99) && (Ncsc==1)

In which case you should use the command:

   h1.Draw("Event", "(abs(Etac[0]-Eta2dsk)>0.1) && (Eta2dsk <99) &&
(Ncsc==1)");

You should find that there is no Event that meets this critiria on your
sample file.

Cheers,
Philippe
-----Original Message-----
From: owner-roottalk@pcroot.cern.ch
[mailto:owner-roottalk@pcroot.cern.ch]On Behalf Of Hannes Sakulin
Sent: Friday, August 16, 2002 5:54 AM
To: roottalk@pcroot.cern.ch
Subject: [ROOT] TTree::Scan() and TTree::Draw() selection different?


Hello,

I wanted to use TTree::Scan() to dump a few events that I had Drawn from a
Tree. It seems that TTree::Draw() and TTree::Scan() use a different
selection. While Draw() finds the events (in the example it is 1 event),
Scan() does not find them. Is there any explanation for this?
root [0] TFile f("gmt_1.root")
root [1] h1.Draw("Event", "(abs(Etac-Eta2dsk)>0.1) && (Eta2dsk <99) &&
(Ncsc==1)")
<TCanvas::MakeDefCanvas>: created default TCanvas with name c1
(Int_t)1
root [2] h1.Scan("Event", "(abs(Etac-Eta2dsk)>0.1) && (Eta2dsk <99) &&
(Ncsc==1)")
************************
*    Row   *     Event *
************************
************************
==> 0 selected entries
(Int_t)0

The behavior is the same in versions 3.02/07, 2.03/02, 3.03/07. You can find
the Tree in ~hsakulin/public/gmt_1.root.

Cheers,
  Hannes.
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
Hannes Sakulin
Institute for High Energy Physics, Vienna, Austria
CERN / EP                           Phone: +41 22 767 7372
CH-1211 Geneva 23, Switzerland      Fax:   +41 22 767 8940
E-mail: Hannes.Sakulin@cern.ch
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~



This archive was generated by hypermail 2b29 : Sat Jan 04 2003 - 23:51:07 MET