Re: [ROOT] Problems with TSelector

From: Rene Brun (Rene.Brun@cern.ch)
Date: Fri Dec 19 2003 - 10:12:41 MET


Hi Christian,

In the ProcessFill function, you should call GetEntry for
the current Tree in the chain, not for the chain itself.
I have modified your code (attachement) to get the correct behaviour.
I have renamed your variable fTree to fChain to show what was the
confusion.
You can also make a variant with a member called fTree in your class
in addition to fChain and compute ftree in Notify.
I will clarify this point in the documentation

Rene Brun

Christian Vogel wrote:
> 
> Hi,
> 
> On Thu, Dec 18, 2003 at 06:39:39PM +0100, Rene Brun wrote:
> > see : http://root.cern.ch/root/roottalk/RoottalkRules.html
> 
> I also tried to create a small minimal test-file, so please see
> the attachments of this email. I also found out that I only read out the
> first file of the chain:
> 
> I generate 3 files with a tree "t1" in it. The Tree is just one Branch,
> a integer, starting at a base of 100 (200,300) and only 4 events are
> stored in each tree:
> 
> emil:tselector_bug$ ./tselector_bug
> Writing file tree1.root, base is 100
> Writing file tree2.root, base is 200
> Writing file tree3.root, base is 300
> new BugSelector!
> Begin->Init(0xbffff8f4)
> Notify()
> Notify()   <----------- 2nd Notify?
> ProcessFill(0) -> 100
> ProcessFill(1) -> 101
> ProcessFill(2) -> 102
> ProcessFill(3) -> 103
> Notify()
> Notify()   <---------- 2nd Notify?
> ProcessFill(0) -> 100  <---- Still reading tree1.root
> Notify()
> Notify()   <----- two notifies per Event?
> ProcessFill(1) -> 101
> Notify()
> Notify()
> ProcessFill(2) -> 102
> Notify()
> Notify()
> ProcessFill(3) -> 103
> Notify()
> Notify()
> ProcessFill(0) -> 100
> Notify()
> Notify()
> ProcessFill(1) -> 101
> Notify()
> Notify()
> ProcessFill(2) -> 102
> Notify()
> Notify()
> ProcessFill(3) -> 103
> bye BugSelector!
> 
>         Chris
> 
> --
> Rather than a beep
> Or a rude error message,
> These words: "File not found."
> -- Len Dvorkin
> 
>   --------------------------------------------------------------------------------
> 
>    tselector_bug.ccName: tselector_bug.cc
>                    Type: Plain Text (text/plain)
> 
>    MakefileName: Makefile
>            Type: Plain Text (text/plain)

#include <stdio.h>

#include <TTree.h>
#include <TChain.h>
#include <TFile.h>
#include <TSelector.h>
#include <TApplication.h>

void write_tree(char *fn,int base){
	Int_t dat;

	fprintf(stderr,"Writing file %s, base is %d\n",fn,base);

	TFile *of=new TFile(fn,"RECREATE");
	TTree *ot=new TTree("t1","Tree 1");
	ot->Branch("DatBranch",&dat,"dat/I");
	for(int i=0;i<4;i++){
		dat=i+base;
		ot->Fill();
	}
	ot->Write();
	delete ot;
	delete of;
}

class BugSelector : public TSelector {
public:
	Int_t dat;
	TTree *fChain;
	BugSelector(){ fprintf(stderr,"new BugSelector!\n"); };
	~BugSelector(){ fprintf(stderr,"bye BugSelector!\n"); };
	virtual Bool_t Notify(){
		fprintf(stderr,"Notify()\n");
		return kTRUE;
	};
	virtual void Begin(TTree *t){ fprintf(stderr,"Begin->"); Init(t); };
	virtual void Init(TTree *t){
		fprintf(stderr,"Init(%p), Version=%d, Status=%d\n",t,Version(),GetStatus());
		fChain=t;
		t->SetBranchAddress("DatBranch",&dat);
	};
        virtual void ProcessFill(Int_t entry){
                fChain->GetTree()->GetEntry(entry);
                fprintf(stderr,"ProcessFill(%d) -> %d\n",entry,dat);
        };

};

int main(int argc,char **argv){
        TApplication app("app",0,0);
	write_tree("tree1.root",100);
	write_tree("tree2.root",200);
	write_tree("tree3.root",300);

	TChain ch("t1","My Chain");
	ch.Add("tree1.root");
	ch.Add("tree2.root");
	ch.Add("tree3.root");

	BugSelector b;
	ch.Process(&b);
}



This archive was generated by hypermail 2b29 : Thu Jan 01 2004 - 17:50:17 MET