Re: [ROOT] C structure as a branch?

From: Rene Brun (Rene.Brun@cern.ch)
Date: Fri Nov 22 2002 - 08:34:31 MET


Hi,

The file I sent to you yesterday was OK.
I repost it (in the attachement).
With this file, just do:
 root > .L test.cc++
 root > test()

Rene Brun

On Thu, 21 Nov 2002, Yiqun Wang wrote:

> Hi,
> 
> So I try to compile the attached "test.cc". If I do not comment out those
> "ClassDef(FitEvent,1)", I got an error:
> 
> /tmp/cc3seyVF.o: In function `SmdPeak::SmdPeak(void)':
> /home/yqwang/working/K0-short/test.cc:14: undefined reference to `FitEvent
> virtual table'
> collect2: ld returned 1 exit status
> gmake: *** [test] Error 1
> 
> If I comment out "ClassDef(FitEvent,1)", compilation is successful.
> However, when I run the program, I will get this error:
> 
> Error in <TTree::Bronch>: Cannot find class:FitEvent
> 
> and I will not get the branch in the root file.
> 
> Yiqun Wang
> 
> > Don't use typedefs, like in your original example, plain structs like
> > this should work too (struct is an all public class in c++):
> >
> > struct SmdPeak {
> >   Float_t peakPos;
> >   Float_t peakEnerg;
> >   Float_t sigmaPeak;
> >   Int_t   corSimuPeak;
> >   bool    used;
> >
> >   SmdPeak() {}
> >
> >   ClassDef (SmdPeak,1)
> > };
> >
> >
> > -- Fons
> >
> >
> > On Thu, 2002-11-21 at 16:25, Rene Brun wrote:
> >> Hi,
> >>
> >> Create classes instead of struct as shown below
> >>
> >> Rene Brun
> >>
> >>
> >> #include "TTree.h"
> >> #include "TFile.h"
> >> #include "TRandom.h"
> >>
> >> class SmdPeak {
> >>
> >> public:
> >>   Float_t peakPos;
> >>   Float_t peakEnerg;
> >>   Float_t sigmaPeak;
> >>   Int_t   corSimuPeak;
> >>   bool    used;
> >>
> >>   SmdPeak() {}
> >>
> >>   ClassDef (SmdPeak,1)
> >> };
> >>
> >> class FitEvent {
> >>
> >> public:
> >>         UInt_t  runnumber;
> >>         Int_t   event;
> >>         Int_t   numbXPeaks;
> >>         Int_t   numbYPeaks;
> >>         Float_t tower[4][3];
> >>         Float_t smdX[61];
> >>         Float_t smdY[101];
> >>         SmdPeak xPeak[4];
> >>
> >>         ClassDef(FitEvent,1)
> >> };
> >>
> >> void struct2()
> >> {
> >>    TFile f("struct.root","recreate");
> >>    TTree T("T","test struct");
> >>    FitEvent *event = new FitEvent;
> >>
> >>    T.Branch("first","FitEvent",&event);
> >>    for (Int_t i=0;i<1000;i++) {
> >>       event->event = i;
> >>       event->runnumber = (UInt_t)(100*gRandom->Rndm());
> >>       //for (Int_t j=0;j<event.count;j++) event.meas[j] = i+j;
> >>       T.Fill();
> >>    }
> >>    T.Print();
> >>    T.Show(45);
> >>    T.Write();
> >> }
> >>
> >>
> >> On Wed, 20 Nov 2002, Yiqun Wang wrote:
> >>
> >> > Hi,
> >> >
> >> > Is there a way to create a branch that is an array of C structure?
> >> The online reference on
> >> >
> >> > TBranch(const char *name, void *address, const char *leaflist, Int_t
> >> basketsize, Int_t compress) :TNamed(name,leaflist)
> >> >
> >> > only talked about a C structure branch.
> >> >
> >> > My header file is like this:
> >> >
> >> >
> >> > typedef struct{
> >> >   Float_t peakPos;
> >> >   Float_t peakEnerg;
> >> >   Float_t sigmaPeak;
> >> >   Int_t   corSimuPeak;
> >> >   bool    used;
> >> > } SmdPeak;
> >> >
> >> > typedef struct{
> >> >
> >> > 	UInt_t  runnumber;
> >> > 	Int_t   event;
> >> > 	Int_t   numbXPeaks;
> >> > 	Int_t   numbYPeaks;
> >> > 	Float_t tower[4][3];
> >> > 	Float_t smdX[61];
> >> > 	Float_t smdY[101];
> >> > 	SmdPeak xPeak[4];
> >> > SmdPeak xPeak[4];
> >> > } FitEvent;
> >> >
> >> >
> >> >
> > --
> > Org:    CERN, European Laboratory for Particle Physics.
> > Mail:   1211 Geneve 23, Switzerland
> > E-Mail: Fons.Rademakers@cern.ch              Phone: +41 22 7679248 WWW:
> >   http://root.cern.ch/~rdm/            Fax:   +41 22 7679480
> 
> 
> 





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