Problem with Dictionnary Generation

From: Frederic Ronga (Frederic.Ronga@iphe.unil.ch)
Date: Thu Dec 17 1998 - 17:26:59 MET


Dear happy ROOTers,
I'm trying to integrate a 'set' of two classes in the ROOT
system. I thought I had done everything the corresponding
tutorial (p. 59: Dictionary Gen. for Interactive...)
tells us, but it doesn't work. When typing the following command:

 rootcint -f ReducedEventCint.cpp -c ReducedEvent.h\
             Pedestal.h ReducedEventLinkDef.h

I get:

 Note: operator new() masked 1c
 Class Pedestal: Streamer() not declared
 Class Pedestal: ShowMembers() not declared
 *** Datamember ReducedEvent::datas: pointer to fundamental type (need manual intervention)
 *** Datamember ReducedEvent::strip: pointer to fundamental type (need manual intervention)

The note is not important, I suppose, and I can manage with the last two warnings, 
but it seems to me that the two remaining ones are a bit worrying.
What am I doing wrong ?

Thanks a lot in advance.

Frederic


***** Technical notes...

I am running under:
  [253]ipndec2-fronga>uname -a
  OSF1 ipndec2 V4.0 878 alpha

and here are the class link definition, headers, implemetations:
(note that I am a brand new ROOTer AND C++-er... Please be lenient)

-------------------------- ReducedEventLinkDef.h -----------------

#ifdef __CINT__

#pragma link off all globals;
#pragma link off all classes;
#pragma link off all functions;

#pragma link C++ class Pedestal-;
#pragma link C++ class ReducedEvent;

#endif


-------------------------- ReducedEvent.h ------------------------

#define _REDUCED_EVENT_

///////////////////////////////////////////////
//                                           //
// Detector event class                      //
// This class is used to store a given event //
// corresponding to a given detector         //
//                                           //
///////////////////////////////////////////////

#include "TObject.h"

#include <math.h>

class Pedestal;

struct reduced
{
     int strip;
     int data;
};

class ReducedEvent : public TObject {
 private:
     int numberOfStrips;                                  // Number of strips in
 this detector
     int firedStrips;                                      // Number of fired st
rips in this event
     int *datas;                                          // Event datas
     int *strip;                                          // Strip number

 public:
     ReducedEvent();
     ReducedEvent(int detNum, const int *eventdatas, const Pedestal *ped, const
int thresh);
     ~ReducedEvent();
     int getDatas() const {return *datas;}
     int getStrips() const {return *strip;}
     int getBinContent(int i) const {return datas[i-1];}
     int getStripNumber(int i) const {return strip[i-1];}
     int getNumberOfStrips() const {return numberOfStrips;}
     int getNumberOfFired() const {return firedStrips;}

     ClassDef(ReducedEvent,1)    // Reduced event class
};
#endif


----------------------- Pedestal.h ---------------------------------

#ifndef _PEDESTAL_
#define _PEDESTAL_

//////////////////////////////////////////////////////////
//                                                      //
// Pedestal class                                       //
// This class is used to compute and store pedestals    //
//                                                      //
//////////////////////////////////////////////////////////

#include "TObject.h"
#include "TH1.h"

class Pedestal : public TObject {
 private:
     int nbchan;
     float *ped;                                          // Pedestal datas
     int *weight;                                         // Number of updates a
lready done for each channel ...
                                                          // ... including build
-up
     TH1F *histo;                                         // Pedestal histogram

 public:
     Pedestal();                                          // Default constructor
     Pedestal(const int *array1, const int *array2, const int thresh, const int
nbchannel);
              // Builds pedestal info from 2 given events (thresh: threshold for
 substracting peaks)
     Pedestal(TH1F *pedestal);
     ~Pedestal();
     void hist();                                         // Histograms pedestal
     int getNumberOfChannels() const {return nbchan;}     // Returns the number
of channels
     TH1F* getHist() const {return histo;}                // Returns pedestal hi
stogram
     // Update pedestals with 2 new events and threshold 'thresh'
     void update(const int *array1, const int *array2, const int thresh);
     void setError(const double *array);
     float getBinContent(const int i) const {return ped[i-1];}
     float getBinError(const int i) const {return histo->GetBinError(i);}

     ClassDef(Pedestal,1);        // Pedestal datas class
};
#endif


------------------------------- ReducedEvent.cpp ------------------------

 #include "ReducedEvent.h"
#include "Pedestal.h"

ClassImp(ReducedEvent)

ReducedEvent::ReducedEvent()
{
     numberOfStrips = 0;
     firedStrips = 0;
}

ReducedEvent::ReducedEvent(const int detNum, const int *eventdatas, const Pedest
al *ped, const int thresh)
{
     int ndet, ngassiplex;
     reduced tmp[ped->getNumberOfChannels()];
     switch (ped->getNumberOfChannels())
     {
     case 1192 : ndet = 5; ngassiplex = 12;  break;
     case 1576 : ndet = 7; ngassiplex = 15; break;
     default : cerr << "Detevt.cc: Unexpected number of channels."; exit(1);
     }

     int nbstrip[ndet];
     nbstrip[0] = nbstrip[1] = 384;
     for (int i = 2; i < ndet-3; i++) nbstrip[i] = 192;
     nbstrip[ndet-1] = 32;

     firedStrips = 0;
     numberOfStrips = nbstrip[detNum];

     int binstart = 0;
     for (int j = 1; j<detNum; j++) binstart = binstart + nbstrip[j];
     binstart = binstart + 1;

     for (int i = binstart; i < binstart + nbstrip[detNum]; i++)            // L
oop on strips
          if ((eventdatas[i] - ped->getBinContent(i)) > thresh*(ped->getBinError
(i)))
          {
               tmp[firedStrips].strip = i;
               tmp[firedStrips].data = (int)trunc(eventdatas[i] - ped->getBinCon
tent(i));
               firedStrips++;
          }
     strip = new int[firedStrips];
     datas = new int[firedStrips];
     for (int i = 0; i < firedStrips; i++)
     {
          strip[i] = tmp[i].strip;
          datas[i] = tmp[i].data;
     }
}

ReducedEvent::~ReducedEvent()
{
     delete [] strip;
     delete [] datas;
}


-------------------------------- Pedestal.cpp -----------------------------

#include "Pedestal.h"

ClassImp(Pedestal)

Pedestal::Pedestal()
{
     nbchan = 0;
}

Pedestal::Pedestal(const int *array1, const int *array2, const int thresh, const
 int nbchannel)
{
     nbchan = nbchannel;
     ped = new float[nbchan];
     weight = new int[nbchan];

     for (int i=0; i<nbchan; i++)
     {
          if (array1[i]-array2[i] > thresh)
          {
               ped[i] = array2[i];
               weight[i] = 1;
          }
          else if (array1[i]-array2[i] < -thresh)
          {
               ped[i] = array1[i];
               weight[i] = 1;
          }
          else
          {
               ped[i] = (array1[i] + array2[i])/2;
               weight[i] = 2;
          }
     }
}


Pedestal::Pedestal(TH1F *pedestal)
{
     if (!pedestal)
     {
          cerr << "Pedestal.cc: No pedestal found." << endl;
          exit(1);
     }
     nbchan = pedestal->GetNbinsX();
     ped = new float[nbchan];
     for (int i = 0; i < nbchan; i++)
          ped[i] = (float)pedestal->GetBinContent(i+1);
     histo = pedestal;
}

Pedestal::~Pedestal()
{
}

void Pedestal::update(const int *array1, const int *array2, const int thresh)
{
     for (int i = 0; i < nbchan; i++)
     {
          if ((!&array1[i]) || (!&array2[i]))
          {
               cerr << "\n Pedestal.cc: cannot update pedestals: bad number of c
hannels.\n";
               exit(3);
          }
          if ((array1[i] - array2[i]) > thresh)
          {
               ped[i] = (ped[i]*weight[i] + array2[i])/(weight[i] + 1);
               weight[i]++;
          }
          else if ((array1[i] - array2[i]) < -thresh)
          {
               ped[i] = (ped[i]*weight[i]+array1[i])/(weight[i]+1);
               weight[i]++;
                  }
          else
          {
               ped[i] = (ped[i]*weight[i] + array1[i] + array2[i])/(weight[i] +
2);
               weight[i] += 2;
          }
     }
}

void Pedestal::hist()
{
     histo = new TH1F("Pedestal","Pedestal",nbchan,0,nbchan);
     for(int i = 0; i < nbchan; i++)
     {
          histo->SetBinContent(i+1,ped[i]);
     }
}

void Pedestal::setError(const double *array)
{
     if (array[0] != nbchan)
     {
          cerr << "Pedestal.cc : Errors array doesn't have right number of chann
els";
          exit(1);
     }

     for (int i=0; i < nbchan; i++)
          histo->SetBinError(i+1,array[i]);
}

-- 
=====================================================================
     /  __  /  /    /  ____/ Frederic RONGA
    /  /   /  /    /  /      Institut de Physique des hautes energies
   /  ____/  __   /  ___/    UNIL - Universite de Lausanne
  /  /      /    /  /        E-mail: Frederic.Ronga@iphe.unil.ch
_/ _/     _/  __/  ____/     URL: http://www-iphe.unil.ch/~fronga
                             Tel. ++41 21 692 37 39
=====================================================================



This archive was generated by hypermail 2b29 : Tue Jan 04 2000 - 00:34:41 MET