[ROOT] TObjArray of TClonesArrays

From: Ruben Shahoian (Ruben.Shahoyan@cern.ch)
Date: Wed Jan 30 2002 - 13:40:58 MET


Hello,
I need to write and read the branc of the tree, which constists of the
TObjArray of TClonesArray (TObjArray being collector of data
(TClonesArray) of different sensor).
So, I do
 
  detectors = new TObjArray(nsens);
  for (int i=0;i<nsens;i++) {
    TClonesArray *arr = new TClonesArray("TObject",3000);
    detectors->AddAt(arr,i);
  }
  file = new TFile("tree.root","RECREATE");
  tree = new TTree("DataTree","Tree");
  branch = tree->Branch("Det1","TObjArray",&detectors,32000,0);

It works fine for writing.
But then, in the reading method, I do 
 TObjArray *detectors=0;
  branch->SetAddress(&detectors);
  ...
  branch->GetEntry(i);

It reads back everything correctly but it appears that for each event
it creates new TClonesArrays and attaches it to TObjArray (without even
deleting the previous TClonesArrays, although, I would expect that they
simply should be reused).
I've tried both splitting the tree and setting SetAutoDelete(kTRUE),
very soon it swallows up all RAM and starts swapping ...
Is this normal behavoiur?
I attach short macros (run .x testrw.C, demonstrating the problem) 

Best regards,
	Ruben



void testwrite(Int_t nev=1000);
void testread(Int_t nev=-1);

void testrw(Int_t nev=1000)
{
  //
  printf("\nWill write %d tree events: TObjArray of TClonesArrays\n",nev);
  testwrite(nev);
  printf("\nNow reading them back\n");
  testread();
  //
}

void testwrite(Int_t nev)
{
  const int kBufferSize=32000;
  const int kSplitLevel=0;
  //
  Int_t nsens = 2;
  TObjArray *detectors;
  TTree *tree;
  TBranch *branch;
  TFile* file;
  //
  // Init detector and sensors
  detectors = new TObjArray(nsens);
  for (int i=0;i<nsens;i++) {
    TClonesArray *arr = new TClonesArray("TObject",3000);
    detectors->AddAt(arr,i);
  }
  //
  // Create the tree with unsplitted branch for this detector
  file = new TFile("tree.root","RECREATE");
  tree = new TTree("DataTree","Tree");
  branch = tree->Branch("Det1","TObjArray",&detectors,kBufferSize,kSplitLevel);
  //
  for (int i=0;i<nev;i++) {
    // Generate Dummy Data for each sensor
    for (int j=0;j<nsens;j++) {
      TClonesArray *arr = (TClonesArray*)detectors->At(j);
      Int_t ndigits = gRandom->Integer(10000);
      for (int k=0;k<ndigits;k++) new ((*arr)[k]) TObject();
    }
    if (!(i%100)) {
      printf("Ev. %5d : %7s %7s %s\n",i,"Sensor","NDigits","&TClonesArray");
      for (int j=0;j<nsens;j++) {
	TClonesArray *arr = (TClonesArray*)detectors->At(j);
	printf("            %7d %7d %p |\n",j,arr->GetLast()+1,arr);
      }
    }
    //
    tree->Fill();
    // clean sensors
    for (int j=0;j<nsens;j++) ((TClonesArray*)detectors->At(j))->Clear(); 
  }
  tree->Write();
  file->Close();
  detectors->Delete();
  delete detectors;
  //
}

//--------------------------------------------------

void testread(int nev)
{
  //
  TObjArray *detectors;
  TTree *tree;
  TBranch *branch;
  TFile* file;
  //
  // Create the tree with unsplitted branch for this detector
  file = new TFile("tree.root");
  tree=(TTree *)file->Get("DataTree");
  branch = tree->GetBranch("Det1");
  branch->SetAddress(&detectors);
  //  branch->SetAutoDelete(kTRUE);
  //
  Int_t nevtot = tree->GetEntries();
  if (nev<1) nev = nevtot;
  else nev = TMath::Min(nev,nevtot);
  //
  for (int i=0;i<nev;i++) {
    // Read the branch for this detector
    branch->GetEntry(i);
    if (!(i%100)) {
      printf("Ev. %5d : %7s %7s %s\n",i,"Sensor","NDigits","&TClonesArray");
      int nsens = detectors->GetLast()+1;
      for (int j=0;j<nsens;j++) {
	TClonesArray *arr = (TClonesArray*)detectors->At(j);
	printf("            %7d  %7d  %p |\n",j,arr->GetLast()+1,arr);
      }
    }

  }
  file->Close();
  //
}



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