[ROOT] derived classes, TTree, Branches and TClonesArray

From: Benoit Revenu (revenu@cdf.in2p3.fr)
Date: Mon Nov 26 2001 - 16:55:48 MET


Hi,

I would like to use two derived classes and put them in two different
branches of a tree. The two classes (PhysEvent and RandEvent) inherit from
GenEvent.
When I try to link the executable, I have the following answer :
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent::StreamerNVirtual(TBuffer &)'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent::DeclFileName(void)'
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent::IsA(void)const'
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent virtual table'
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent::DeclFileName(void)'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent::Class_Version(void)'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent::IsA(void)const'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent::DeclFileLine(void)'
/opt/Libs//libanalyse.so: undefined reference to`GenEvent::Class_Name(void)'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent virtual table'
/opt/Libs//libanalyse.so: undefined reference to`GenEvent::ShowMembers(TMemberInspector &, char *)'
/opt/Libs//libanalyse.so: undefined reference to`RandEvent::StreamerNVirtual(TBuffer &)'
/opt/Libs//libanalyse.so: undefined reference to`Station::Print(ostream&) const'
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent::Class_Version(void)'
/opt/Libs//libanalyse.so: undefined reference to`GenEvent::Streamer(TBuffer &)
/opt/Libs//libanalyse.so: undefined reference to`PhysEvent::DeclFileLine(void)'

and I don't understand why !
Is there a special way to use ClassDef and/or ClassImp when dealing with
derived classes ?

There are the files in the following order :
1) anaclass.h
2) anaclass.cc
3) Makefile
4) analyse.cc (contains the main routine)
5) AnaLinkDef.h

==============================================================
cat anaclass.h
==============================================================
#ifndef _ANACLASS_H_
#define _ANACLASS_H_

#include <string>
#include <vector>
#include <iostream>

#include "TROOT.h"
#include "TObject.h"
#include "TClonesArray.h"

using namespace std;

class Station : public TObject
{

public:

  Station(){};
  ~Station(){};

  double test;

  void Print(ostream& os) const;

  ClassDef(Station,1)

};

class GenEvent : public TObject
{
public:

  GenEvent();
  virtual ~GenEvent();

  string fileName;

  TClonesArray *tca;

  Int_t Nstat;
  Int_t totNstat;

  void AddStation(const Station&);
  void Clear();

  ClassDef(GenEvent,1)

};

class PhysEvent : public GenEvent
{
public:

  PhysEvent();
  ~PhysEvent();

  Double_t a;

  ClassDef(PhysEvent,1)

};

class RandEvent : public GenEvent
{
public:

  RandEvent();
  ~RandEvent();

  Double_t b;

  ClassDef(RandEvent,1)

};
#endif

=============================================================
cat anaclass.cc
==============================================================
#include "anaclass.h"

ClassImp(Station)
ClassImp(GenEvent)
ClassImp(PhysEvent)
ClassImp(RandEvent)

#define MAXPEVTS 20
#define MAXREVTS 20

/**************** GenEvent ***************/

GenEvent::GenEvent()
{
  tca = new TClonesArray("Station",MAXPEVTS);
  Nstat = 0;
  totNstat = 0;
}

GenEvent::~GenEvent()
{
  if( tca != NULL ) delete tca;
}

void GenEvent::AddStation(const Station& st)
{
  TClonesArray &newstation = *tca;
  new(newstation[Nstat++]) Station(st);
  totNstat++;
}

void GenEvent::Clear()
{
  tca->Clear();
  Nstat = 0;
}

/******************** PhysEvent *******************/

PhysEvent::PhysEvent() : GenEvent() {}

/******************** RandEvent *******************/

RandEvent::RandEvent() : GenEvent() {}

==============================================================
cat Makefile
==============================================================
ROOTDIR=/opt/root/
ROOTBIN=${ROOTDIR}bin/
ROOTLIB=${ROOTDIR}lib/

TMP=/opt/
OBJSDIR=${TMP}Objs/
LIBSDIR=${TMP}Libs/
MOBJS=${OBJSDIR}AnaDict.o ${OBJSDIR}analib.o ${OBJSDIR}anaclass.o
CFLAGS=`${ROOTBIN}root-config --cflags` -Wall -DDPA -g
LIBS=`${ROOTBIN}root-config --libs` -Wl,-rpath ${ROOTLIB}
EXE=./analyse

all: ${EXE}

${OBJSDIR}%.o: %.cc
        ${CXX} -c $< -o $@ ${CFLAGS}

AnaDict.cc: analib.h anaclass.h AnaLinkDef.h
        ${ROOTBIN}rootcint -f AnaDict.cc -c analib.h anaclass.h
AnaLinkDef.h

${LIBSDIR}libanalyse.a: ${MOBJS}
        ar r $@ ${MOBJS}

${LIBSDIR}libanalyse.so: ${MOBJS}
        g++ -shared -o $@ ${MOBJS}

clean:
        rm -f ${OBJSDIR}*.o ${EXE}
        rm -f ${LIBSDIR}*.a ${LIBSDIR}*.so
        rm -f EbDict.* AnaDict.* core

analyse: ${LIBSDIR}libanalyse.a ${LIBSDIR}libanalyse.so
${OBJSDIR}analyse.o
        ${CXX} ${CFLAGS} $(OBJS) ${LIBS} ${OBJSDIR}analyse.o \
        -L${LIBSDIR} -lanalyse -o $@

==============================================================
cat analyse.cc (contains the main routine)
==============================================================
#include "TROOT.h"
#include "TApplication.h"
#include "analib.h"
#include "anaclass.h"
#include <TFile.h>
#include <TTree.h>


TROOT root("Example","Example of data reading for the Auger experiment");

TTree * resultP;
TTree * resultR;
TFile * outfile;
PhysEvent * pevt;
RandEvent * revt;

int main()
{
  TApplication theApp ( "App", NULL, 0 );
  gROOT->Reset( );

  outfile = new TFile( "/tmp/outfile.root", "RECREATE", "Demo root file
with ntuple" );

  resultP = new TTree( "ResultP", "EAS parameters" );
  resultR = new TTree( "ResultR", "EAS parameters" );

  pevt = new PhysEvent();
  revt = new RandEvent();

  int bufsize = 10000;

  Int_t split = 1;
  resultP->Branch( "PhysEvents", "PhysEvent", &pevt, bufsize, split );
  resultR->Branch( "RandEvents", "RandEvent", &revt, bufsize, split );

  resultP->AutoSave();
  resultR->AutoSave();
  resultR->GetBranch("RandEvents")->SetAddress(&revt);
  resultP->GetBranch("PhysEvents")->SetAddress(&pevt);

  Station stat;

  for( int i=0;i<100;i++ )
    {
      stat.test = i;
      revt->AddStation(stat);
      resultR->Fill();
      revt->Clear();
    }

  outfile->Write();
  outfile->Close();

  return 0;
}

==============================================================
cat AnaLinkDef.h
==============================================================
#ifdef __CINT__
 
#pragma link off all globals;
#pragma link off all classes;
#pragma link off all functions;
#pragma link off all typedef;

#pragma link C++ class GenEvent;
#pragma link C++ class PhysEvent;
#pragma link C++ class RandEvent;
#pragma link C++ class Station;

#endif



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