RE: [ROOT] tree splitting

From: Philippe Canal (pcanal@fnal.gov)
Date: Mon Apr 23 2001 - 15:47:58 MEST


Hi Susan,

We are in the middle of upgrading the TTreeFormula facilities.  All 3 of
your requests should be implemented in the next few weeks.

Philippe.

-----Original Message-----
From: owner-roottalk@pcroot.cern.ch
[mailto:owner-roottalk@pcroot.cern.ch]On Behalf Of Susan Kasahara
Sent: Sunday, April 22, 2001 9:47 PM
To: Rene Brun
Cc: Susan Kasahara; roottalk@pcroot.cern.ch
Subject: Re: [ROOT] tree splitting


Hi Rene,
I have tried out the Simple example now with the latest CVS version 3.01
of ROOT and find that it does work.  Very impressive!
We are hoping to use this new support of storing polymorphic objects
on a single ROOT TTree branch in our data model, and I have a new set
of questions related to how we can apply selection cuts to objects stored on
this type of branch.
Specifically, (using the Simple example to illustrate the types of cuts
I'd like to apply):
1) I would like to apply selection cuts to the SimpleTree fShape
branch by accessor method, e.g.:
"fShape -> GetVisibility() > 0"
In the past, I've found that I can use TTreeFormula to apply selection
cuts that include the use of accessor methods, but I find that with
the new TTree structure (TTree::SetBranchStyle(1);) that I cannot
successfully apply a cut using the above selection string.  (I attach
the write and read driver code that I'm using to test this.  The
Simple class is the same as that used in previous examples.) It does appear
that applying a selection cut without the use of accessor methods, e.g.:
"fID > 0"
works successfully.

2) I would find it even better to be able to apply selection cuts to the
SimpleTree fShape branch using accessor methods from the derived
class of the object stored there.
For example, for an entry on the fShape branch containing a TTUBE,
I would like to be able to apply a selection cut of the form,
"fShape -> GetRmin() > 100"
I have been able to apply cuts such as this outside of TTreeFormula,
using the TMethodCall class, as in:

//...After retrieving simple object from TTree for a given entry
TShape* shape = simple -> GetShape();
TMethodCall *methodcall = new TMethodCall(shape->IsA(),"GetRmin","");
if (methodcall -> GetMethod()) {
 // method "GetRmin" was found for this object
 Double_t rmin;
 methodcall -> Execute(shape,rmin);
 if (rmin > 100) {
  // Entry passes cuts
 }
}
delete methodcall;

but, I am wondering if the TTreeFormula class can also support this type
of dynamic use of accessor methods for the entries on the polymorphic
fShape branch.

3) I'd also find it desirable to be able to string accessor methods
together,
e.g. apply a selection cut of the form:
"fShape -> GetMaterial() -> GetZ() > 8";
Is support for this something that can be added to the TTreeFormula class?
I notice that the set of TMethodCall::Execute methods includes those
which return basic data types, but not one which returns a pointer to
a TObject.  If there were such a version of TMethodCall::Execute, than it
seems like you could string the calls to the TMethodCall::Execute methods
together, e.g.

//...After retrieving simple object from TTree for a given entry
TShape* shape = simple -> GetShape();
TMethodCall *methodone = new TMethodCall(shape->IsA(),"GetMaterial","");
if (methodone -> GetMethod()) {
  // method "GetMaterial" was found for this object
  // TMethodCall::ReturnType() could be used to determine type of return
value
  TObject* object;
  methodone -> Execute(shape,object);
  TMethodCall *methodtwo = new TMethodCall(object->IsA(),"GetZ","");
  if (methodtwo) {
    // method "GetZ" was found for this object
    Double_t z;
    methodtwo -> Execute(object,z);
    if (z > 8) {
       // Entry passes cuts
    }
  }
}
Is there a fundamental reason why a TMethodCall::Execute method which
returns a pointer to a TObject doesn't exist or is it just that no one
has requested this before?

Thanks for all your help.  I am using the ROOT CVS version downloaded
today (4/22/01) and gcc 2.95.2.
-Sue

************* Driver code to write data **************************

//
// Driver code to open file, create tree, and write 3 Simple objects to
//
// tree.
//
#include "TROOT.h"
#include "TFile.h"
#include "TTree.h"
#include "TTUBE.h"
#include "Simple.h"

TROOT test("test","Program to test simple classes");
int main(int argc, char **argv)
{

// Open output file to hold events
  TFile *simplefile = new TFile("simple.root","RECREATE","Simple root
file");

// Create a ROOT tree to hold the simple data
  TTree *simpletree = new TTree("SimpleTree","Simple Tree");

  Simple *simple = 0;
  // Split branches to store fShape objects & fID on separate branches
  TTree::SetBranchStyle(1);
  simpletree->Branch("Simple","Simple",&simple,16000,1);

  for (Int_t ient=0; ient < 3; ient++) {
// Create 3 Simple objects containing TTube's and store each in tree
    TTUBE *tube = new TTUBE("TUBE","TUBE","void",150,200,400);
    // tube is adopted & deleted by Simple
    simple = new Simple(ient,tube);
    simple -> Print();
    simpletree -> Fill();
    delete simple;
  }

// Write the simple tree to file
  simplefile->Write();
// Print final summary of simple tree contents
  simpletree->Print();

  return 0;
}

****** Driver code to read data makes use of TTreeFormula *********
//
// Code to open root file containing tree, Get tree, and attempt to read in
//
// 3 Simple objects from tree making use of TTreeFormula to apply selection
// cuts.
//
#include "TROOT.h"
#include "TFile.h"
#include "TTree.h"
#include "TTUBE.h"
#include "Simple.h"
#include "TTreeFormula.h"

TROOT test("test","Program to test simple classes");

int main(int argc, char **argv)
{

// Open file to read Simple entries
  TFile *file = new TFile("simple.root","READ");

// Retrieve ROOT tree holding the data
  TTree *tree = dynamic_cast<TTree*>(file -> Get("SimpleTree"));

  // const char* select = "fID > 0"; // this works
  // The following selection cut should select all 3 entries, but doesn't
work
  const char* select = "fShape->GetVisibility() > 0";

  TTreeFormula* treeformula = new TTreeFormula("myselection",select,tree);
  Int_t npass = 0;

  Int_t nent = (Int_t)tree -> GetEntries();

  for (Int_t ient=0; ient < nent; ient++) {
    Simple *simple = new Simple();
    tree -> SetBranchAddress("Simple",&simple);

    // change TTree fReadEntry to entry of interest
    tree -> LoadTree(ient);

    // Ask the TTreeFormula to evaluate if this entry passes selection cut
    if (( treeformula -> EvalInstance()) > 0) {
      // entry passes the cuts, retrieve the entire record
      tree -> GetEntry(ient);
      npass++;
      simple -> Print();
    }

    delete simple;
  }

  // Finished all entries
  cout << "Total entries in tree = " << nent << ".  Total passing cuts = "
       << npass << endl;

  return 0;

}

************* Simple class same as before *****************

#ifndef SIMPLE_H
#define SIMPLE_H
//
//
// Simple class
//
#include "TShape.h"
#include <iostream.h>

class Simple : public TObject {

private:
   Int_t   fID;         // id number
   TShape* fShape;      // pointer to base class shape

public:

   Simple() : fID(0), fShape(0) { }
   Simple(Int_t id, TShape* shape): fID(id), fShape(shape) { }
   virtual ~Simple();
   virtual void Print(Option_t *option = "") const;

   ClassDef(Simple,1)  //Simple class
};

#endif                         // SIMPLE_H

//
// Simple.cxx
//
#include <iostream.h>
#include "Simple.h"

ClassImp(Simple)

Simple::~Simple() {
// Destructor
  if (fShape) {
    delete fShape;
    fShape =0;
  }
}

void Simple::Print(Option_t *option) const {
  // Print the contents
  cout << "fID= " << fID << endl;
  fShape -> Print();

}




Rene Brun wrote:

> Hi Sue,
> Thanks for testing the new Branch mechanism.
> I have fixed this problem happening when deleting and recreating
> the object at address &simple in the Fill loop.
> The fix is now in CVS.
>
> Rene Brun



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