From f6497d9d6af01ff415cf6aa40371cfd8dfbf5465 Mon Sep 17 00:00:00 2001 From: hristov Date: Wed, 13 Aug 2014 17:53:45 +0200 Subject: [PATCH] ALIROOT-5559 Fix for memory jump in CPass0 (LHC10 processing). A custom streamer is implemented, where the "big data members" are written directly to the file in parallel with the original object. The custom streamers are arror prone: every time we have a change in the data members, in addition to the incremented version the streamer has to be adapter. --- TPC/Calib/AliTPCcalibAlign.cxx | 150 +++++++++++++++++++++++++++++++++ TPC/Calib/AliTPCcalibAlign.h | 4 +- TPC/TPCcalibLinkDef.h | 2 +- 3 files changed, 154 insertions(+), 2 deletions(-) diff --git a/TPC/Calib/AliTPCcalibAlign.cxx b/TPC/Calib/AliTPCcalibAlign.cxx index 1bc2fff434e..bbb60b4def3 100644 --- a/TPC/Calib/AliTPCcalibAlign.cxx +++ b/TPC/Calib/AliTPCcalibAlign.cxx @@ -139,6 +139,8 @@ #include "Riostream.h" #include "TRandom.h" #include + +#include "AliSysInfo.h" using namespace std; AliTPCcalibAlign* AliTPCcalibAlign::fgInstance = 0; @@ -3816,3 +3818,151 @@ void AliTPCcalibAlign::MakeReportDyPhi(TFile */*output*/){ } +//______________________________________________________________________________ +void AliTPCcalibAlign::Streamer(TBuffer &R__b) +{ + // Stream an object of class AliTPCcalibAlign. + Bool_t isDebug=AliLog::GetDebugLevel("","AliTPCcalibAlign")>0; + if (isDebug) AliSysInfo::SetVerbose(kTRUE); + + if (R__b.IsReading()) { + // R__b.ReadClassBuffer(AliTPCcalibAlign::Class(),this); + UInt_t R__s, R__c; + Version_t R__v = R__b.ReadVersion(&R__s, &R__c); + + AliTPCcalibBase::Streamer(R__b); // Stream the base class + + // Read the "big data members" from the same Root directory + if (gDirectory){ + for (Int_t i=0; i<2; ++i){ + TString hisName=TString::Format("AliTPCcalibAlign.%s.fClusterDelta_%d",GetName(),i); + if (gDirectory->Get(hisName.Data())){ + fClusterDelta[i] = dynamic_cast(gDirectory->Get((hisName).Data())); + } + } + + for (Int_t i=0; i<4; ++i){ + TString hisName=TString::Format("AliTPCcalibAlign.%s.fTrackletDelta_%d",GetName(),i); + if (gDirectory->Get(hisName.Data())){ + fTrackletDelta[i] = dynamic_cast(gDirectory->Get((hisName).Data())); + } + } + } + + // If the "big data members"were not in the file, try to read them from the object + // This is suppose to restore the backward compatibility + for (Int_t i=0; i<2; ++i) { + if (!fClusterDelta[i]) R__b >> fClusterDelta[i]; + } + + for (Int_t i=0; i<4; ++i) { + if (!fTrackletDelta[i]) R__b >> fTrackletDelta[i]; + } + + // Read all the other data members. This is error prone, + // please make sure this part is updated if you change the data members + fDphiHistArray.Streamer(R__b); + fDthetaHistArray.Streamer(R__b); + fDyHistArray.Streamer(R__b); + fDzHistArray.Streamer(R__b); + + fDyPhiHistArray.Streamer(R__b); + fDzThetaHistArray.Streamer(R__b); + + fDphiZHistArray.Streamer(R__b); + fDthetaZHistArray.Streamer(R__b); + fDyZHistArray.Streamer(R__b); + fDzZHistArray.Streamer(R__b); + + fFitterArray12.Streamer(R__b); + fFitterArray9.Streamer(R__b); + fFitterArray6.Streamer(R__b); + + fMatrixArray12.Streamer(R__b); + fMatrixArray9.Streamer(R__b); + fMatrixArray6.Streamer(R__b); + + fCombinedMatrixArray6.Streamer(R__b); + + R__b.ReadFastArray(fPoints,72*72); + + R__b >> fNoField; + R__b >> fXIO; + R__b >> fXmiddle; + R__b >> fXquadrant; + + fArraySectorIntParam.Streamer(R__b); + fArraySectorIntCovar.Streamer(R__b); + + R__b >> fSectorParamA; + R__b >> fSectorCovarA; + R__b >> fSectorParamC; + R__b >> fSectorCovarC; + + R__b >> fUseInnerOuter; + R__b >> fgkMergeEntriesCut; // Should we write the static member? + } else { + // R__b.WriteClassBuffer(AliTPCcalibAlign::Class(),this); + if (isDebug) AliSysInfo::AddStamp("aliengStreamer::Start"); + + R__b.WriteVersion(AliTPCcalibAlign::IsA()); + AliTPCcalibBase::Streamer(R__b); // Stream the base class + + // Write the "big data members directly in the file in parallel with the object itself + for (Int_t i=0; i<2; ++i){ + if (fClusterDelta[i]) fClusterDelta[i]->Write(TString::Format("AliTPCcalibAlign.%s.fClusterDelta_%d",GetName(),i).Data()); + } + + if (isDebug) AliSysInfo::AddStamp("alignStreamer::fClusterDelta"); + + for (Int_t i=0; i<4; ++i){ + if (fTrackletDelta[i]) fTrackletDelta[i]->Write(TString::Format("AliTPCcalibAlign.%s.fTrackletDelta_%d",GetName(),i).Data()); + } + + if (isDebug) AliSysInfo::AddStamp("alignStreamer::fTrackletDelta"); + + // Write all the other data members. This is error prone, + // please make sure this part is updated if you change the data members + fDphiHistArray.Streamer(R__b); + fDthetaHistArray.Streamer(R__b); + fDyHistArray.Streamer(R__b); + fDzHistArray.Streamer(R__b); + + fDyPhiHistArray.Streamer(R__b); + fDzThetaHistArray.Streamer(R__b); + + fDphiZHistArray.Streamer(R__b); + fDthetaZHistArray.Streamer(R__b); + fDyZHistArray.Streamer(R__b); + fDzZHistArray.Streamer(R__b); + + fFitterArray12.Streamer(R__b); + fFitterArray9.Streamer(R__b); + fFitterArray6.Streamer(R__b); + + fMatrixArray12.Streamer(R__b); + fMatrixArray9.Streamer(R__b); + fMatrixArray6.Streamer(R__b); + + fCombinedMatrixArray6.Streamer(R__b); + + R__b.WriteFastArray(fPoints,72*72); + + R__b << fNoField; + R__b << fXIO; + R__b << fXmiddle; + R__b << fXquadrant; + + fArraySectorIntParam.Streamer(R__b); + fArraySectorIntCovar.Streamer(R__b); + + R__b << fSectorParamA; + R__b << fSectorCovarA; + R__b << fSectorParamC; + R__b << fSectorCovarC; + + R__b << fUseInnerOuter; + R__b << fgkMergeEntriesCut; // Should we write the static member? + if (isDebug) AliSysInfo::AddStamp("alignStreamer::DefaultStreamer"); + } +} diff --git a/TPC/Calib/AliTPCcalibAlign.h b/TPC/Calib/AliTPCcalibAlign.h index 9cf738cad3c..c48d1bf97bf 100644 --- a/TPC/Calib/AliTPCcalibAlign.h +++ b/TPC/Calib/AliTPCcalibAlign.h @@ -206,7 +206,9 @@ protected: private: AliTPCcalibAlign& operator=(const AliTPCcalibAlign&);// not implemented - ClassDef(AliTPCcalibAlign,6) + // IMPORTANT: If you change the data members, + // please do not forget to increment the ClassDef and to update the Streamer in AliTPCcalibAlign.cxx + ClassDef(AliTPCcalibAlign,7) }; diff --git a/TPC/TPCcalibLinkDef.h b/TPC/TPCcalibLinkDef.h index 584e2ab7878..afe28e452b4 100644 --- a/TPC/TPCcalibLinkDef.h +++ b/TPC/TPCcalibLinkDef.h @@ -22,7 +22,7 @@ // Duplicted also (most probably) AliTPCCalPad --- to be cleaned up #pragma link C++ class AliTPCcalibCalib+; // Re-applying calib on cluster level - refitting of the tracks -#pragma link C++ class AliTPCcalibAlign+; // (Sector)-alignment calibration +#pragma link C++ class AliTPCcalibAlign-; // (Sector)-alignment calibration // Histogram cluster to track residuals in order to create later on distortion maps // --- update documentation (current documentation is obsolete) -- 2.43.0