#ifdef __cplusplus extern "C" { #endif #pragma once #include "HROOT-histType.h" #include "HROOTHistTH1.h" #include "HROOTCoreTNamed.h" #include "HROOTCoreTAttLine.h" #include "HROOTCoreTAttFill.h" #include "HROOTCoreTAttMarker.h" #include "HROOTCoreTObject.h" #include "STDDeletable.h" #include "HROOT-coreType.h" #define TH2_DECL_VIRT(Type) \ int Type##_fill2 ( Type##_p p, double x, double y );\ int Type##_fill2w ( Type##_p p, double x, double y, double w );\ void Type##_fillN2 ( Type##_p p, int ntimes, double* x, double* y, double* w, int stride );\ void Type##_fillRandom2 ( Type##_p p, TH1_p h, int ntimes );\ int Type##_findFirstBinAbove2 ( Type##_p p, double threshold, int axis );\ int Type##_findLastBinAbove2 ( Type##_p p, double threshold, int axis );\ void Type##_FitSlicesX ( Type##_p p, TF1_p f1, int firstybin, int lastybin, int cut, const char* option, TObjArray_p arr );\ void Type##_FitSlicesY ( Type##_p p, TF1_p f1, int firstxbin, int lastxbin, int cut, const char* option, TObjArray_p arr );\ double Type##_getCorrelationFactor2 ( Type##_p p, int axis1, int axis2 );\ double Type##_getCovariance2 ( Type##_p p, int axis1, int axis2 );\ double Type##_integral2 ( Type##_p p, int binx1, int binx2, int biny1, int biny2, const char* option );\ TH2_p Type##_rebinX2 ( Type##_p p, int ngroup, const char* newname );\ TH2_p Type##_rebinY2 ( Type##_p p, int ngroup, const char* newname );\ TH2_p Type##_Rebin2D ( Type##_p p, int nxgroup, int nygroup, const char* newname );\ void Type##_SetShowProjectionX ( Type##_p p, int nbins );\ void Type##_SetShowProjectionY ( Type##_p p, int nbins ); #define TH2_DECL_NONVIRT(Type) \ TH1D_p Type##_tH2_ProjectionX ( Type##_p p, const char* name, int firstybin, int lastybin, const char* option );\ TH1D_p Type##_tH2_ProjectionY ( Type##_p p, const char* name, int firstxbin, int lastxbin, const char* option ); #define TH2_DECL_ACCESSOR(Type) \ #define TH2_DEF_VIRT(Type) \ int Type##_fill2 ( Type##_p p, double x, double y ) {\ return ((TYPECASTMETHOD(Type, fill2, TH2))(p))->Fill(x, y);\ }\ \ int Type##_fill2w ( Type##_p p, double x, double y, double w ) {\ return ((TYPECASTMETHOD(Type, fill2w, TH2))(p))->Fill(x, y, w);\ }\ \ void Type##_fillN2 ( Type##_p p, int ntimes, double* x, double* y, double* w, int stride ) {\ ((TYPECASTMETHOD(Type, fillN2, TH2))(p))->FillN(ntimes, x, y, w, stride);\ }\ \ void Type##_fillRandom2 ( Type##_p p, TH1_p h, int ntimes ) {\ ((TYPECASTMETHOD(Type, fillRandom2, TH2))(p))->FillRandom(from_nonconst_to_nonconst(h), ntimes);\ }\ \ int Type##_findFirstBinAbove2 ( Type##_p p, double threshold, int axis ) {\ return ((TYPECASTMETHOD(Type, findFirstBinAbove2, TH2))(p))->FindFirstBinAbove(threshold, axis);\ }\ \ int Type##_findLastBinAbove2 ( Type##_p p, double threshold, int axis ) {\ return ((TYPECASTMETHOD(Type, findLastBinAbove2, TH2))(p))->FindLastBinAbove(threshold, axis);\ }\ \ void Type##_FitSlicesX ( Type##_p p, TF1_p f1, int firstybin, int lastybin, int cut, const char* option, TObjArray_p arr ) {\ ((TYPECASTMETHOD(Type, FitSlicesX, TH2))(p))->FitSlicesX(from_nonconst_to_nonconst(f1), firstybin, lastybin, cut, option, from_nonconst_to_nonconst(arr));\ }\ \ void Type##_FitSlicesY ( Type##_p p, TF1_p f1, int firstxbin, int lastxbin, int cut, const char* option, TObjArray_p arr ) {\ ((TYPECASTMETHOD(Type, FitSlicesY, TH2))(p))->FitSlicesY(from_nonconst_to_nonconst(f1), firstxbin, lastxbin, cut, option, from_nonconst_to_nonconst(arr));\ }\ \ double Type##_getCorrelationFactor2 ( Type##_p p, int axis1, int axis2 ) {\ return ((TYPECASTMETHOD(Type, getCorrelationFactor2, TH2))(p))->GetCorrelationFactor(axis1, axis2);\ }\ \ double Type##_getCovariance2 ( Type##_p p, int axis1, int axis2 ) {\ return ((TYPECASTMETHOD(Type, getCovariance2, TH2))(p))->GetCovariance(axis1, axis2);\ }\ \ double Type##_integral2 ( Type##_p p, int binx1, int binx2, int biny1, int biny2, const char* option ) {\ return ((TYPECASTMETHOD(Type, integral2, TH2))(p))->Integral(binx1, binx2, biny1, biny2, option);\ }\ \ TH2_p Type##_rebinX2 ( Type##_p p, int ngroup, const char* newname ) {\ return from_nonconst_to_nonconst((TH2*)((TYPECASTMETHOD(Type, rebinX2, TH2))(p))->RebinX(ngroup, newname));\ }\ \ TH2_p Type##_rebinY2 ( Type##_p p, int ngroup, const char* newname ) {\ return from_nonconst_to_nonconst((TH2*)((TYPECASTMETHOD(Type, rebinY2, TH2))(p))->RebinY(ngroup, newname));\ }\ \ TH2_p Type##_Rebin2D ( Type##_p p, int nxgroup, int nygroup, const char* newname ) {\ return from_nonconst_to_nonconst((TH2*)((TYPECASTMETHOD(Type, Rebin2D, TH2))(p))->Rebin2D(nxgroup, nygroup, newname));\ }\ \ void Type##_SetShowProjectionX ( Type##_p p, int nbins ) {\ ((TYPECASTMETHOD(Type, SetShowProjectionX, TH2))(p))->SetShowProjectionX(nbins);\ }\ \ void Type##_SetShowProjectionY ( Type##_p p, int nbins ) {\ ((TYPECASTMETHOD(Type, SetShowProjectionY, TH2))(p))->SetShowProjectionY(nbins);\ } #define TH2_DEF_NONVIRT(Type) \ TH1D_p Type##_tH2_ProjectionX ( Type##_p p, const char* name, int firstybin, int lastybin, const char* option ) {\ return from_nonconst_to_nonconst((TH1D*)((TYPECASTMETHOD(Type, tH2_ProjectionX, TH2))(p))->ProjectionX(name, firstybin, lastybin, option));\ }\ \ TH1D_p Type##_tH2_ProjectionY ( Type##_p p, const char* name, int firstxbin, int lastxbin, const char* option ) {\ return from_nonconst_to_nonconst((TH1D*)((TYPECASTMETHOD(Type, tH2_ProjectionY, TH2))(p))->ProjectionY(name, firstxbin, lastxbin, option));\ } #define TH2_DEF_ACCESSOR(Type) \ TH1_DECL_VIRT(TH2) TNAMED_DECL_VIRT(TH2) TATTLINE_DECL_VIRT(TH2) TATTFILL_DECL_VIRT(TH2) TATTMARKER_DECL_VIRT(TH2) TOBJECT_DECL_VIRT(TH2) DELETABLE_DECL_VIRT(TH2) TH2_DECL_VIRT(TH2) TH2_DECL_NONVIRT(TH2) TH2_DECL_ACCESSOR(TH2) #ifdef __cplusplus } #endif