"manual" schema evolution

From: Chiara Zampolli <Chiara.Zampolli_at_bo.infn.it>
Date: Thu, 12 Nov 2009 11:16:28 +0100


Hello all,

    Following the slides showed By Philippe at last CHEP (http://indicobeta.cern.ch/getFile.py/access?contribId=210&sessionId=59&resId=3&materialId=slides&confId=35523), I'm trying to define my own rules to allow some schema evolution. I'm attching the headers file of the old version of the class (AliESDfriendTrack.h.Orig), and of the new one (AliESDfriendTrack.h), and the LinkDef file (ESDLinkDef.h) where I define the rules. The changes in the class affects some arrays of integeres that now I would like to have as pointes. So:

    old header file:                                  ---->    new
    header file:
    Int_t fITSindex[kMaxITScluster];     ---->    Int_t fnMaxITScluster;
    Int_t* fITSindices; //fnMaxITScluster
    Int_t fTPCindex[kMaxTPCcluster];  ---->    Int_t fnMaxITScluster;
    Int_t* fITSindices; //fnMaxITScluster
    Int_t fTRDindex[kMaxTRDcluster];  ---->    Int_t fnMaxITScluster;
    Int_t* fITSindices; //fnMaxITScluster


In total, I have defined six rules. In the first three:

    #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="" version="[-3]" target="fnMaxITScluster" targetType="Int_t"\     code="{fnMaxITScluster = 12;}"
    #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="" version="[-3]" target="fnMaxTPCcluster" targetType="Int_t"\     code="{fnMaxTPCcluster = 160;}"
    #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="" version="[-3]" target="fnMaxTRDcluster" targetType="Int_t"\     code="{fnMaxTRDcluster = 180;}"

I use hardwired numbers, but I think I could use the kMaxITS/TPC/TRDcluster variables, but I'm not sure... In the last three:

    #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="Int_t fITSindex" version="[-3]" target="fITSindices:     fnMaxITScluster" targetType="Int_t*; Int_t"\     code="{fITSindices = new Int_t[fnMaxITScluster]; for (Int_t i=0; i<     fnMaxITScluster; i++){fITSindices[i] = onfile.fITSindex[i];}}"     //memcpy(fITSindices,
    &(onfile.fITSindex),fnMaxITScluster*sizeof(Int_t); }"     #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="Int_t fTPCindex" version="[-3]" target="fTPCindices;     fnMaxTPCcluster" targetType="Int_t*; Int_t"\     code="{fTPCindices = new Int_t[fnMaxTPCcluster]; memcpy(fTPCindices,     &(onfile.fTPCindex),fnMaxTPCcluster*sizeof(Int_t); }"     #pragma read sourceClass="AliESDfriendTrack"     targetClass="AliESDfriendTrack" \
    source="Int_t fTRDindex" version="[-3]" target="fTRDindices;     fnMaxTRDcluster" targetType="Int_t*; Int_t"\     code="{fTRDindices = new Int_t[fnMaxTRDcluster]; memcpy(fTRDindices,     &(onfile.fTRDindex),fnMaxTRDcluster*sizeof(Int_t); }"

I've put two targets, since I use two variables of the new implementation of the class, but I'm not sure that the definition of their type is then done correctly in the way I do. And, anyway, when I try to read an old object, the old content is not copied in the new variables - just the inizialization values (-2, see .cxx files) are "read". Moreover, I see no track of my rules inthe LinkDef file, so I guess there's something not working as it should...

Do you have any idea what I'm doing wrong?

Thanks for your help, and cheers,

        Chiara

/**************************************************************************

* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
* *
* Author: The ALICE Off-line Project. *
* Contributors are mentioned in the code where appropriate. *
* *
* Permission to use, copy, modify and distribute this software and its *
* documentation strictly for non-commercial purposes is hereby granted *
* without fee, provided that the above copyright notice appears in all *
* copies and that both the copyright notice and this permission notice *
* appear in the supporting documentation. The authors make no claims *
* about the suitability of this software for any purpose. It is *
* provided "as is" without express or implied warranty. *
**************************************************************************/ //------------------------------------------------------------------------- // Implementation of the AliESDfriendTrack class // This class keeps complementary to the AliESDtrack information // Origin: Iouri Belikov, CERN, Jouri.Belikov_at_cern.ch //------------------------------------------------------------------------- #include "AliTrackPointArray.h" #include "AliESDfriendTrack.h"

#include "TObjArray.h"
#include "AliKalmanTrack.h"

ClassImp(AliESDfriendTrack)

AliESDfriendTrack::AliESDfriendTrack(): TObject(),
f1P(0),
fPoints(0),
fCalibContainer(0),
fITStrack(0),
fTRDtrack(0),
fTPCOut(0),
fITSOut(0),
fTRDIn(0)
{

  //
  // Default constructor
  //

  Int_t i;
  for (i=0; i<kMaxITScluster; i++) fITSindex[i]=-2;
  for (i=0; i<kMaxTPCcluster; i++) fTPCindex[i]=-2;
  for (i=0; i<kMaxTRDcluster; i++) fTRDindex[i]=-2;
}

AliESDfriendTrack::AliESDfriendTrack(const AliESDfriendTrack &t):

TObject(t),
f1P(t.f1P),
fPoints(0),

fCalibContainer(0),
fITStrack(0),
fTRDtrack(0),
fTPCOut(0),
fITSOut(0),
fTRDIn(0)
{
  //
  // Copy constructor
  //

  Int_t i;
  for (i=0; i<kMaxITScluster; i++) fITSindex[i]=t.fITSindex[i];
  for (i=0; i<kMaxTPCcluster; i++) fTPCindex[i]=t.fTPCindex[i];
  for (i=0; i<kMaxTRDcluster; i++) fTRDindex[i]=t.fTRDindex[i];
  if (t.fPoints) fPoints=new AliTrackPointArray(*t.fPoints);   if (t.fCalibContainer) {
     fCalibContainer = new TObjArray(5);
     Int_t no=t.fCalibContainer->GetEntriesFast();
     for (i=0; i<no; i++) {
       TObject *o=t.fCalibContainer->At(i);
       fCalibContainer->AddLast(o->Clone());
     }  

  }
  if (t.fTPCOut) fTPCOut = new AliExternalTrackParam(*(t.fTPCOut));
  if (t.fITSOut) fITSOut = new AliExternalTrackParam(*(t.fITSOut));
  if (t.fTRDIn)  fTRDIn = new AliExternalTrackParam(*(t.fTRDIn));
}

AliESDfriendTrack::~AliESDfriendTrack() {

  //
  // Simple destructor
  //

   delete fPoints;
   if (fCalibContainer) fCalibContainer->Delete();    delete fCalibContainer;
   delete fITStrack;
   delete fTRDtrack;
   delete fTPCOut;
   delete fITSOut;
   delete fTRDIn;
}

void AliESDfriendTrack::AddCalibObject(TObject * calibObject){

  //
  // add calibration object to array -
  // track is owner of the objects in the container 
  //

  if (!fCalibContainer) fCalibContainer = new TObjArray(5);   fCalibContainer->AddLast(calibObject); }

TObject * AliESDfriendTrack::GetCalibObject(Int_t index){

  //
  //
  //

  if (!fCalibContainer) return 0;
  if (index>=fCalibContainer->GetEntriesFast()) return 0;   return fCalibContainer->At(index);
}

void AliESDfriendTrack::SetTPCOut(const AliExternalTrackParam &param) {

  // 
  // backup TPC out track
  //

  delete fTPCOut;
  fTPCOut=new AliExternalTrackParam(param); }
void AliESDfriendTrack::SetITSOut(const AliExternalTrackParam &param) {
  //
  // backup ITS out track
  //

  delete fITSOut;
  fITSOut=new AliExternalTrackParam(param); }
void AliESDfriendTrack::SetTRDIn(const AliExternalTrackParam &param) {
  //
  // backup TRD in track
  //

  delete fTRDIn;
  fTRDIn=new AliExternalTrackParam(param); }

#ifndef ALIESDFRIENDTRACK_H
#define ALIESDFRIENDTRACK_H

//-------------------------------------------------------------------------
//                     Class AliESDfriendTrack
//               This class contains ESD track additions
//       Origin: Iouri Belikov, CERN, Jouri.Belikov_at_cern.ch 
//-------------------------------------------------------------------------

#include <TObject.h>
#include <AliExternalTrackParam.h>

class AliTrackPointArray;
class AliKalmanTrack;
class TObjArrray;

//_____________________________________________________________________________
class AliESDfriendTrack : public TObject { public:
  enum {
    kMaxITScluster=12,
    kMaxTPCcluster=160,
    kMaxTRDcluster=180

  };
  AliESDfriendTrack();
  AliESDfriendTrack(const AliESDfriendTrack &t);   virtual ~AliESDfriendTrack();

  void Set1P(Float_t p) {f1P=p;}
  void SetTrackPointArray(AliTrackPointArray *points) {     fPoints=points;
  }

  Float_t Get1P() const  {return f1P;}
  Int_t *GetITSindices() {return fITSindices;}
  Int_t *GetTPCindices() {return fTPCindices;}
  Int_t *GetTRDindices() {return fTRDindices;}
  const AliTrackPointArray *GetTrackPointArray() const {return fPoints;}

  void SetITStrack(AliKalmanTrack *t) {fITStrack=t;}   void SetTRDtrack(AliKalmanTrack *t) {fTRDtrack=t;}   AliKalmanTrack *GetTRDtrack() {return fTRDtrack;}   AliKalmanTrack *GetITStrack() {return fITStrack;}   void AddCalibObject(TObject * calibObject);   TObject * GetCalibObject(Int_t index);   //
  // parameters backup

  void SetTPCOut(const AliExternalTrackParam &param);
  void SetITSOut(const AliExternalTrackParam &param);
  void SetTRDIn(const AliExternalTrackParam  &param);
  //
  const AliExternalTrackParam * GetTPCOut(){return  fTPCOut;} 
  const AliExternalTrackParam * GetITSOut() { return fITSOut;} 
  const AliExternalTrackParam * GetTRDIn() { return fTRDIn;} 

protected:
  Float_t f1P;                     // 1/P (1/(GeV/c))
  Int_t fnMaxITScluster; // max n. of ITS clusters
  Int_t fnMaxTPCcluster; // max n. of TPC clusters   Int_t fnMaxTRDcluster; // max n. of TRD clusters
  //  Int_t fITSindex[kMaxITScluster]; // indices of the ITS clusters 
  //Int_t fTPCindex[kMaxTPCcluster]; // indices of the TPC clusters
  //Int_t fTRDindex[kMaxTRDcluster]; // indices of the TRD clusters
  Int_t* fITSindices; //[fnMaxITScluster] indices of the ITS clusters 
  Int_t* fTPCindices; //[fnMaxTPCcluster] indices of the TPC clusters   Int_t* fTRDindices; //[fnMaxTRDcluster] indices of the TRD clusters   AliTrackPointArray *fPoints;//Array of track space points in the global frame
  TObjArray      *fCalibContainer; //Array of objects for calibration    
  AliKalmanTrack *fITStrack; //! pointer to the ITS track (debug purposes) 
  AliKalmanTrack *fTRDtrack; //! pointer to the TRD track (debug purposes) 
  //
  //
  AliExternalTrackParam * fTPCOut; // tpc outer parameters
  AliExternalTrackParam * fITSOut; // its outer parameters
  AliExternalTrackParam * fTRDIn;  // trd inner parameters

private:
  AliESDfriendTrack &operator=(const AliESDfriendTrack & /* t */) {return *this;}

  ClassDef(AliESDfriendTrack,4) //ESD friend track };

#endif

#ifndef ALIESDFRIENDTRACK_H
#define ALIESDFRIENDTRACK_H

//-------------------------------------------------------------------------
//                     Class AliESDfriendTrack
//               This class contains ESD track additions
//       Origin: Iouri Belikov, CERN, Jouri.Belikov_at_cern.ch 
//-------------------------------------------------------------------------

#include <TObject.h>
#include <AliExternalTrackParam.h>

class AliTrackPointArray;
class AliKalmanTrack;
class TObjArrray;

//_____________________________________________________________________________
class AliESDfriendTrack : public TObject { public:
  enum {
    kMaxITScluster=12,
    kMaxTPCcluster=160,
    kMaxTRDcluster=180

  };
  AliESDfriendTrack();
  AliESDfriendTrack(const AliESDfriendTrack &t);   virtual ~AliESDfriendTrack();

  void Set1P(Float_t p) {f1P=p;}
  void SetTrackPointArray(AliTrackPointArray *points) {     fPoints=points;
  }

  Float_t Get1P() const  {return f1P;}
  Int_t *GetITSindices() {return fITSindex;}
  Int_t *GetTPCindices() {return fTPCindex;}
  Int_t *GetTRDindices() {return fTRDindex;}
  const AliTrackPointArray *GetTrackPointArray() const {return fPoints;}

  void SetITStrack(AliKalmanTrack *t) {fITStrack=t;}   void SetTRDtrack(AliKalmanTrack *t) {fTRDtrack=t;}   AliKalmanTrack *GetTRDtrack() {return fTRDtrack;}   AliKalmanTrack *GetITStrack() {return fITStrack;}   void AddCalibObject(TObject * calibObject);   TObject * GetCalibObject(Int_t index);   //
  // parameters backup

  void SetTPCOut(const AliExternalTrackParam &param);
  void SetITSOut(const AliExternalTrackParam &param);
  void SetTRDIn(const AliExternalTrackParam  &param);
  //
  const AliExternalTrackParam * GetTPCOut(){return  fTPCOut;} 
  const AliExternalTrackParam * GetITSOut() { return fITSOut;} 
  const AliExternalTrackParam * GetTRDIn() { return fTRDIn;} 

protected:
  Float_t f1P;                     // 1/P (1/(GeV/c))
  Int_t fITSindex[kMaxITScluster]; // indices of the ITS clusters 
  Int_t fTPCindex[kMaxTPCcluster]; // indices of the TPC clusters   Int_t fTRDindex[kMaxTRDcluster]; // indices of the TRD clusters   AliTrackPointArray *fPoints;//Array of track space points in the global frame
  TObjArray      *fCalibContainer; //Array of objects for calibration    
  AliKalmanTrack *fITStrack; //! pointer to the ITS track (debug purposes) 
  AliKalmanTrack *fTRDtrack; //! pointer to the TRD track (debug purposes) 
  //
  //
  AliExternalTrackParam * fTPCOut; // tpc outer parameters
  AliExternalTrackParam * fITSOut; // its outer parameters
  AliExternalTrackParam * fTRDIn;  // trd inner parameters

private:
  AliESDfriendTrack &operator=(const AliESDfriendTrack & /* t */) {return *this;}

  ClassDef(AliESDfriendTrack,3) //ESD friend track };

#endif

#ifdef __CINT__
/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *

/* $Id: ESDLinkDef.h 34989 2009-09-25 08:11:29Z kleinb $ */

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

#pragma link C++ enum AliESDEvent::ESDListIndex;

#pragma link C++ class  AliESD+;
#pragma link C++ class  AliESDEvent+;
#pragma link C++ class  AliESDInputHandler+;
#pragma link C++ class  AliESDInputHandlerRP+;
#pragma link C++ class  AliESDRun+;
#pragma link C++ class  AliESDHeader+;
#pragma link C++ class  AliESDZDC+;
#pragma link C++ class  AliESDCaloTrigger+;
#pragma link C++ class  AliESDfriend+;
#pragma read sourceClass="AliESDtrack" targetClass="AliESDtrack" source="UChar_t fTRDpidQuality"  version="[-49]" target="fTRDntracklets" targetType="UChar_t" code="{fTRDntracklets = onfile.fTRDpidQuality;}"
// see http://root.cern.ch/svn/root/trunk/io/doc/DataModelEvolution.txt #pragma link C++ class AliESDtrack+;

#pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="" version="[-3]" target="fnMaxITScluster" targetType="Int_t"\    code="{fnMaxITScluster = 12;}"
#pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="" version="[-3]" target="fnMaxTPCcluster" targetType="Int_t"\    code="{fnMaxTPCcluster = 160;}"
#pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="" version="[-3]" target="fnMaxTRDcluster" targetType="Int_t"\    code="{fnMaxTRDcluster = 180;}"

#pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="Int_t fITSindex" version="[-3]" target="fITSindices: fnMaxITScluster" targetType="Int_t*; Int_t"\    code="{fITSindices = new Int_t[fnMaxITScluster]; for (Int_t i=0; i< fnMaxITScluster; i++){fITSindices[i] = onfile.fITSindex[i];}}" //memcpy(fITSindices, &(onfile.fITSindex), fnMaxITScluster*sizeof(Int_t); }" #pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="Int_t fTPCindex" version="[-3]" target="fTPCindices; fnMaxTPCcluster" targetType="Int_t*; Int_t"\    code="{fTPCindices = new Int_t[fnMaxTPCcluster]; memcpy(fTPCindices, &(onfile.fTPCindex), fnMaxTPCcluster*sizeof(Int_t); }" #pragma read sourceClass="AliESDfriendTrack" targetClass="AliESDfriendTrack" \

   source="Int_t fTRDindex" version="[-3]" target="fTRDindices; fnMaxTRDcluster" targetType="Int_t*; Int_t"\    code="{fTRDindices = new Int_t[fnMaxTRDcluster]; memcpy(fTRDindices, &(onfile.fTRDindex), fnMaxTRDcluster*sizeof(Int_t); }" #pragma link C++ class AliESDfriendTrack+;

#pragma link C++ class  AliESDMuonTrack+;
#pragma link C++ class  AliESDPmdTrack+;
#pragma link C++ class  AliESDTrdTrack+;
#pragma link C++ class  AliESDHLTtrack+;
#pragma link C++ class  AliESDv0+;
#pragma link C++ class  AliESDcascade+;
#pragma link C++ class  AliVertex+;
#pragma link C++ class  AliESDVertex+;
#pragma link C++ class  AliESDpid+;
#pragma link C++ class  AliESDkink+;
#pragma link C++ class  AliESDV0Params+;
#pragma link C++ class  AliESDCaloCluster+;
#pragma link C++ class  AliESDMuonCluster+;
#pragma link C++ class  AliESDMuonPad+;

#pragma link C++ class  AliKFParticleBase+;
#pragma link C++ class AliKFParticle+;
#pragma link C++ class AliKFVertex+;
#pragma link C++ class  AliKalmanTrack+;
#pragma link C++ class  AliNeutralTrackParam+;
#pragma link C++ class  AliVertexerTracks+;
#pragma link C++ class  AliStrLine+;
#pragma link C++ class  AliTrackPointArray+;
#pragma link C++ class  AliTrackPoint+;

#pragma link C++ class AliTrackPointArray+; #pragma link C++ class AliTrackPoint+;

#pragma link C++ class  AliESDFMD+;
#pragma link C++ class  AliFMDMap+;
#pragma link C++ class  AliFMDFloatMap+;

#pragma link C++ class  AliESDVZERO+;

#pragma link C++ class AliESDTZERO+;
#pragma link C++ class AliESDACORDE+;

#pragma link C++ class AliESDMultITS+;
#pragma link C++ class AliMultiplicity+;

#pragma link C++ class AliSelector+;

#pragma link C++ class AliRawDataErrorLog+;

#pragma link C++ class AliMeanVertex+;
#pragma link C++ class AliESDCaloCells+;

#pragma link C++ class AliTriggerIR+;

#pragma link C++ class AliESDVZEROfriend+;

#pragma link C++ class AliTPCpidESD+;
#pragma link C++ class AliTOFpidESD+;

#pragma link C++ class  AliTriggerScalersESD+;
#pragma link C++ class  AliTriggerScalersRecordESD+;
#pragma link C++ class  AliESDZDCScalers+;
#pragma link C++ class  AliESDHandler+;

#pragma link C++ function operator*(const AliFMDMap&,const AliFMDMap&);
#pragma link C++ function operator/(const AliFMDMap&,const AliFMDMap&);
#pragma link C++ function operator+(const AliFMDMap&,const AliFMDMap&); #pragma link C++ function operator-(const AliFMDMap&,const AliFMDMap&);   

#endif

Received on Thu Nov 12 2009 - 11:16:38 CET

This archive was generated by hypermail 2.2.0 : Thu Nov 12 2009 - 17:50:04 CET