RE: [ROOT] Scans and Draws

From: Philippe Canal (pcanal@fnal.gov)
Date: Tue Oct 16 2001 - 22:30:47 MEST


Hi Sidoti,

First I point out a few details (I tested with root 3.02.01) and then 
I discussed better alternative to build your trees.

> if I do a scan all the entries are the same:
>   tree->Scan("ptsg3x.c:ptsint.c")

Actually scroll down to entries passed 100 and you will see that the
2 columns seem to be the same for 0 to 100 and are different for 101 to 201.

> if I do a scan requiring a selection --> less entry but again the first
>                                        two columns are the same.

The entries that are the same are printed ....

Try tree->Scan("ptsg3x.c:ptsint.c","abs(ptsg3x.c-ptsint.c)>1","lego")
and see that 10 entries (between 0 and 100) are different and all entries
beyhond 100 are also different.

> If I plot:
> tree->Draw("ptsg3x.c-ptsint.c") -->I see two peaks at -1000 and +1000 (and
                                   also a peak at 0)

I see 4 'spots' being filled: (0,0) with about 50 entries, (-1000,-1000) with
about 50 entries, (0,-1000) with about 10 entries and (-1000,0) with about 100
entries.   

> If I plot with constraints --> I see effectively only one peak in 0
> tree->Draw("ptsg3x.c-ptsint.c","abs(ptsg3x.c-ptsint.c)<1","lego")

As expected, I see 2 'spots' filled (0,0): with about 50 entries, 
(-1000,-1000) with about 50 entries,

> 1) if the entries are NOT the same why the TTree:Scan shows me same
> columns

Most of the first 100 entries are the same but not the last 100.

> 2) are the two branches different?

They are different on about 10 entries.

> 3) I see a peak at 200 when tree->Draw("ptsg3x.obsp") I am expecting only
> one or 0 entry per channel.

This is actually an indication of your core issues.  As you noticed TTree::Draw
has been using 202 values rather than 101 values you would have expected.

This is solely due to the way you updated your tree.
You did:

  TFile hfile("help.root","UPDATE","compare simple int and g3x");
  ....
  tree->Branch("ptsg3x",&comp,,"ev/I:obsp/I:c/F:d0/F:phi0/F:z0/F:ctg/F");
  for (Int_t jentry=0; jentry<nentries;jentry++) {
    ....
    tree->Fill();

In other word, you added a new branch to a already filled tree and then 
proceed to call TTree::Fill.  This had the side effect of filling ALL
existing branch.  This means that an additional 101 entries were added to
the ptsint.  Those entries all have some default/random value inserted.
Now you have a poorly balanced tree with the branch ptsint branch having
202 entries (101 of which are junk) and the ptsg3x branch having only 
101 entries....

You could force the Draw and Scan method to give you the expected result
by limiting the range to the first 101 entries (4th parameter).

A much much better alternative it to never add a branch to an existing TTree
and to use TTree FriendShip instead.  So your new code will look like:

  TFile hfile("help.root","UPDATE","compare simple int and g3x");

  TTree *tree = new TTree("base","Tree containing the base comparaison values");

  tree->Branch("ptsg3x",&comp,,"ev/I:obsp/I:c/F:d0/F:phi0/F:z0/F:ctg/F");
  ....

  tree->AddFriend("cmp");
  hfile->Write();
  hfile.Close();

and then do you analysis as before with:

  TFile *help = TFile::Open("help.root")
  TTree *tree = (TTree*)help->Get("base");
  
Cheers,
Philippe.

-----Original Message-----
From: owner-roottalk@pcroot.cern.ch
[mailto:owner-roottalk@pcroot.cern.ch]On Behalf Of Sidoti Antonio tel.
+39+0461 88 1525
Sent: Tuesday, October 16, 2001 2:06 PM
To: roottalk@pcroot.cern.ch
Subject: [ROOT] Scans and Draws


Hi rooters,
I am having a problem. I want to check if the entries two branches of a
rootuple are the same (see the attached help.root file).

TFile *help = TFile::Open("help.root")
TTree *tree = (TTree*)help->Get("cmp");

if I do a scan all the entries are the same:
    tree->Scan("ptsg3x.c:ptsint.c")

if I do a scan requiring a selection --> less entry but again the first
                                         two columns are the same.
    tree->Scan("ptsg3x.c:ptsint.c","abs(ptsg3x.c-ptsint.c)>1")

If I plot:
tree->Draw("ptsg3x.c-ptsint.c") -->I see two peaks at -1000 and +1000 (and
                                   also a peak at 0)

If I plot with constraints --> I see effectively only one peak in 0
tree->Draw("ptsg3x.c-ptsint.c",abs(ptsg3x.c-ptsint.c)<1")

Some questions:
1) if the entries are NOT the same why the TTree:Scan shows me same
columns
2) are the two branches different?
3) I see a peak at 200 when tree->Draw("ptsg3x.obsp") I am expecting only
one or 0 entry per channel.
...
Any ideas?
Thanks in advance
Antonio Sidoti

P.S I'm using  3.01/06

P.S2
I am filling the root file in the following way...

  /*TFile hfile("help.root","RECREATE","compare simple int and g3x");
  TTree *tree = new TTree("cmp","Tree with different branches");

tree->Branch("ptsint",&compn,"ev/I:obsp/I:c/F:d0/F:phi0/F:z0/F:ctg/F");*/

  TFile hfile("help.root","UPDATE","compare simple int and g3x");
  TTree *tree = (TTree*)hfile->Get("cmp");
  tree->Branch("ptsg3x",&comp,,"ev/I:obsp/I:c/F:d0/F:phi0/F:z0/F:ctg/F");

   for (Int_t jentry=0; jentry<nentries;jentry++) {
    Int_t ientry = LoadTree(jentry); //in case of a TChain, ientry is the
entry number in the current file
    nb = fChain->GetEntry(jentry);   nbytes += nb;
    compn.ev = EV_evnum;
    compn.obsp = MC_obsp;
    compn.c = PTS_c;
    compn.d0 = PTS_d0;
    compn.phi0 = PTS_phi0;
    compn.z0 = PTS_z0;
    compn.ctg = PTS_ctg;
    tree->Fill();
  }
  hfile->Write();
  hfile.Close();
}



_______________________________________________________________________________
	     a n t o n i o   	      s i d o t i

	     e-mail	              sidoti@science.unitn.it
	                              sidoti@fnal.gov
	     www   		      http://higgs.tn.infn.it/~sidoti/
_______________________________________________________________________________
		"Il meglio e` nemico del buono"



This archive was generated by hypermail 2b29 : Tue Jan 01 2002 - 17:51:03 MET