Skip to content
Snippets Groups Projects
Commit 776e1dcc authored by Joshua Falco Beirer's avatar Joshua Falco Beirer Committed by Rachid Mazini
Browse files

Added custom streamer to TFCSEnergyInterpolationPiecewiseLinear

parent acb86c81
No related branches found
No related tags found
No related merge requests found
......@@ -296,7 +296,7 @@
#pragma link C++ class TFCSInvisibleParametrization+;
#pragma link C++ class TFCSInitWithEkin+;
#pragma link C++ class TFCSEnergyInterpolationLinear+;
#pragma link C++ class TFCSEnergyInterpolationPiecewiseLinear+;
#pragma link C++ class TFCSEnergyInterpolationPiecewiseLinear-;
#pragma link C++ class TFCSEnergyInterpolationSpline+;
#pragma link C++ class TFCSParametrizationChain-;
#pragma link C++ class TFCSParametrizationBinnedChain+;
......
......@@ -13,6 +13,8 @@
#include "TCanvas.h"
#include "TGraph.h"
#include "TAxis.h"
//TBuffer include required for custom class streamer
#include "TBuffer.h"
class TFCSEnergyInterpolationPiecewiseLinear: public TFCSParametrization {
......@@ -34,20 +36,21 @@ public:
void InitFromArrayInLogEkin(Int_t np, const Double_t logEkin[], const Double_t response[]);
void InitFromArrayInEkin(Int_t np, const Double_t Ekin[], const Double_t response[]);
virtual FCSReturnCode simulate(TFCSSimulationState& simulstate,const TFCSTruthState* truth, const TFCSExtrapolationState* extrapol) const override;
virtual FCSReturnCode simulate(TFCSSimulationState& simulstate, const TFCSTruthState* truth, const TFCSExtrapolationState* extrapol) const override;
void Print(Option_t *option="") const override;
static void unit_test(TFCSSimulationState* simulstate = nullptr,TFCSTruthState* truth = nullptr, const TFCSExtrapolationState* extrapol = nullptr, TGraph* grlinear = nullptr);
static void unit_test(TFCSSimulationState* simulstate = nullptr, TFCSTruthState* truth = nullptr, const TFCSExtrapolationState* extrapol = nullptr, TGraph* grlinear = nullptr);
private:
ROOT::Math::Interpolator m_linInterpol;
double m_minXValue;
double m_maxXValue;
int m_nPoints;
ROOT::Math::Interpolator m_linInterpol; //! Do not persistify
ClassDefOverride(TFCSEnergyInterpolationPiecewiseLinear, 1) //TFCSEnergyInterpolationPiecewiseLinear
std::vector<Double_t> m_logEkin;
std::vector<Double_t> m_response;
std::pair<float, float> m_MinMaxlogEkin;
ClassDefOverride(TFCSEnergyInterpolationPiecewiseLinear, 2) //TFCSEnergyInterpolationPiecewiseLinear
};
......
......@@ -29,10 +29,15 @@ TFCSEnergyInterpolationPiecewiseLinear::TFCSEnergyInterpolationPiecewiseLinear(c
void TFCSEnergyInterpolationPiecewiseLinear::InitFromArrayInLogEkin(Int_t np, const Double_t logEkin[], const Double_t response[])
{
m_linInterpol.SetData(np, logEkin, response);
m_minXValue = *std::min_element(logEkin, logEkin + np);
m_maxXValue = *std::max_element(logEkin, logEkin + np);
m_nPoints = np;
//save logEkin and response as std::vector class members
//this is required for using the custom streamer
m_logEkin.assign (logEkin, logEkin + np);
m_response.assign(response, response + np);
m_linInterpol.SetData(m_logEkin, m_response);
auto min_max = std::minmax_element(m_logEkin.begin(), m_logEkin.end());
m_MinMaxlogEkin = std::make_pair(*min_max.first, *min_max.second);
}
void TFCSEnergyInterpolationPiecewiseLinear::InitFromArrayInEkin(Int_t np, const Double_t Ekin[], const Double_t response[])
......@@ -52,11 +57,11 @@ FCSReturnCode TFCSEnergyInterpolationPiecewiseLinear::simulate(TFCSSimulationSta
const float logEkin = Ekin > Gaudi::Units::keV ? TMath::Log(Ekin) : TMath::Log(Gaudi::Units::keV);
float Emean;
if(logEkin < m_minXValue){
Emean = m_linInterpol.Eval(m_minXValue)*Einit;
if(logEkin < m_MinMaxlogEkin.first){
Emean = m_linInterpol.Eval(m_MinMaxlogEkin.first)*Einit;
}
else if(logEkin > m_maxXValue){
Emean = (m_linInterpol.Eval(m_maxXValue) + m_linInterpol.Deriv(m_maxXValue) * (logEkin - m_maxXValue))*Einit;
else if(logEkin > m_MinMaxlogEkin.second){
Emean = (m_linInterpol.Eval(m_MinMaxlogEkin.second) + m_linInterpol.Deriv(m_MinMaxlogEkin.second) * (logEkin - m_MinMaxlogEkin.second))*Einit;
}
else{
Emean = m_linInterpol.Eval(logEkin)*Einit;
......@@ -73,6 +78,7 @@ FCSReturnCode TFCSEnergyInterpolationPiecewiseLinear::simulate(TFCSSimulationSta
void TFCSEnergyInterpolationPiecewiseLinear::Print(Option_t *option) const
{
TString opt(option);
bool shortprint = opt.Index("short")>=0;
bool longprint = msgLvl(MSG::DEBUG) || (msgLvl(MSG::INFO) && !shortprint);
......@@ -80,11 +86,25 @@ void TFCSEnergyInterpolationPiecewiseLinear::Print(Option_t *option) const
optprint.ReplaceAll("short","");
TFCSParametrization::Print(option);
if(longprint) ATH_MSG_INFO(optprint <<(OnlyScaleEnergy()?" E()*":" Ekin()*")<<"linInterpol N="<<m_nPoints
<<" "<<m_minXValue<<"<=log(Ekin)<="<<m_maxXValue
<<" "<<TMath::Exp(m_minXValue)<<"<=Ekin<="<<TMath::Exp(m_maxXValue));
if(longprint) ATH_MSG_INFO(optprint <<(OnlyScaleEnergy()?" E()*":" Ekin()*")<<"linInterpol N="<<m_logEkin.size()
<<" "<<m_MinMaxlogEkin.first<<"<=log(Ekin)<="<<m_MinMaxlogEkin.second
<<" "<<TMath::Exp(m_MinMaxlogEkin.first)<<"<=Ekin<="<<TMath::Exp(m_MinMaxlogEkin.second));
}
void TFCSEnergyInterpolationPiecewiseLinear::Streamer(TBuffer &R__b)
{
//Stream an object of class TFCSEnergyInterpolationPiecewiseLinear
if (R__b.IsReading()) {
//read the class buffer
R__b.ReadClassBuffer(TFCSEnergyInterpolationPiecewiseLinear::Class(), this);
//initialize interpolation from saved class members
InitFromArrayInLogEkin(m_logEkin.size(), m_logEkin.data(), m_response.data());
}
else {
R__b.WriteClassBuffer(TFCSEnergyInterpolationPiecewiseLinear::Class(), this);
}
}
void TFCSEnergyInterpolationPiecewiseLinear::unit_test(TFCSSimulationState* simulstate, TFCSTruthState* truth, const TFCSExtrapolationState* extrapol, TGraph* grlinear)
{
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment