ATLAS Offline Software
TFCSEnergyInterpolationLinear.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3 */
4 
9 #include "TCanvas.h"
10 #include "TGraph.h"
11 #include "TAxis.h"
12 #include <iostream>
13 
14 //=============================================
15 //======= TFCSEnergyInterpolationLinear =========
16 //=============================================
17 
19  const char *title)
20  : TFCSParametrization(name, title), m_slope(1), m_offset(0) {}
21 
24  const TFCSTruthState *truth,
25  const TFCSExtrapolationState *) const {
26  const float Emean = m_slope * truth->Ekin() + m_offset;
27 
28  ATH_MSG_DEBUG("set E=" << Emean << " for true Ekin=" << truth->Ekin());
29  simulstate.set_E(Emean);
30 
31  return FCSSuccess;
32 }
33 
34 void TFCSEnergyInterpolationLinear::Print(Option_t *option) const {
35  TString opt(option);
36  bool shortprint = opt.Index("short") >= 0;
37  bool longprint = msgLvl(MSG::DEBUG) || (msgLvl(MSG::INFO) && !shortprint);
38  TString optprint = opt;
39  optprint.ReplaceAll("short", "");
41 
42  if (longprint)
43  ATH_MSG_INFO(optprint << " Emean=" << m_slope << "*Ekin(true) + "
44  << m_offset);
45 }
46 
48  TFCSSimulationState *simulstate, TFCSTruthState *truth,
50  if (!simulstate)
51  simulstate = new TFCSSimulationState();
52  if (!truth)
53  truth = new TFCSTruthState();
54  if (!extrapol)
56 
57  TFCSEnergyInterpolationLinear test("testTFCSEnergyInterpolation",
58  "test linear TFCSEnergyInterpolation");
59  test.set_pdgid(22);
60  test.set_Ekin_nominal(1000);
61  test.set_Ekin_min(1000);
62  test.set_Ekin_max(100000);
63  test.set_eta_nominal(0.225);
64  test.set_eta_min(0.2);
65  test.set_eta_max(0.25);
66  test.set_slope(0.95);
67  test.set_offset(-50);
68  test.Print();
69 
70  truth->set_pdgid(22);
71 
72  TGraph *gr = new TGraph();
73  gr->SetNameTitle("testTFCSEnergyInterpolation",
74  "test linear TFCSEnergyInterpolation");
75  gr->GetXaxis()->SetTitle("Ekin [MeV]");
76  gr->GetYaxis()->SetTitle("<E(reco)>/Ekin(true)");
77 
78  int ip = 0;
79  for (float Ekin = 1000; Ekin <= 100000; Ekin *= 1.1) {
80  // Init LorentzVector for truth. For photon Ekin=E
81  truth->SetPxPyPzE(Ekin, 0, 0, Ekin);
82  if (test.simulate(*simulstate, truth, extrapol) != FCSSuccess) {
83  return;
84  }
85  gr->SetPoint(ip, Ekin, simulstate->E() / Ekin);
86  ++ip;
87  }
88 
89 // Drawing doesn't make sense inside athena and necessary libraries not linked
90 // by default
91 #if defined(__FastCaloSimStandAlone__)
92  TCanvas *c = new TCanvas("testTFCSEnergyInterpolation",
93  "test linear TFCSEnergyInterpolation");
94  gr->Draw("APL");
95  c->SetLogx();
96 #endif
97 }
FCSReturnCode
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
Definition: TFCSParametrizationBase.h:41
TFCSEnergyInterpolationLinear::m_offset
float m_offset
Definition: TFCSEnergyInterpolationLinear.h:39
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
gr
#define gr
TrigInDetValidation_Base.test
test
Definition: TrigInDetValidation_Base.py:146
TFCSExtrapolationState
Definition: TFCSExtrapolationState.h:13
RunActsMaterialValidation.extrapol
extrapol
Definition: RunActsMaterialValidation.py:90
TFCSEnergyInterpolationLinear.h
TFCSEnergyInterpolationLinear::Print
void Print(Option_t *option="") const override
Definition: TFCSEnergyInterpolationLinear.cxx:34
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
TFCSTruthState::Ekin
double Ekin() const
Definition: TFCSTruthState.h:26
covarianceTool.title
title
Definition: covarianceTool.py:542
TFCSEnergyInterpolationLinear::unit_test
static void unit_test(TFCSSimulationState *simulstate=nullptr, TFCSTruthState *truth=nullptr, const TFCSExtrapolationState *extrapol=nullptr)
Definition: TFCSEnergyInterpolationLinear.cxx:47
TFCSEnergyInterpolationLinear
Definition: TFCSEnergyInterpolationLinear.h:10
find_tgc_unfilled_channelids.ip
ip
Definition: find_tgc_unfilled_channelids.py:3
TFCSParametrizationBase::Print
void Print(Option_t *option="") const
Print object information.
Definition: TFCSParametrizationBase.cxx:52
FCSSuccess
@ FCSSuccess
Definition: TFCSParametrizationBase.h:41
TFCSParametrization
Definition: TFCSParametrization.h:10
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
pmontree.opt
opt
Definition: pmontree.py:16
TFCSSimulationState::set_E
void set_E(int sample, double Esample)
Definition: TFCSSimulationState.h:48
TFCSTruthState.h
TFCSExtrapolationState.h
DEBUG
#define DEBUG
Definition: page_access.h:11
TFCSSimulationState::E
double E() const
Definition: TFCSSimulationState.h:42
TFCSSimulationState.h
TFCSEnergyInterpolationLinear::m_slope
float m_slope
Definition: TFCSEnergyInterpolationLinear.h:38
TFCSEnergyInterpolationLinear::simulate
virtual FCSReturnCode simulate(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
Method in all derived classes to do some simulation.
Definition: TFCSEnergyInterpolationLinear.cxx:23
TFCSTruthState
Definition: TFCSTruthState.h:13
python.compressB64.c
def c
Definition: compressB64.py:93
TFCSSimulationState
Definition: TFCSSimulationState.h:32
TFCSEnergyInterpolationLinear::TFCSEnergyInterpolationLinear
TFCSEnergyInterpolationLinear(const char *name=nullptr, const char *title=nullptr)
Definition: TFCSEnergyInterpolationLinear.cxx:18
ISF_FCS::MLogging::msgLvl
bool msgLvl(const MSG::Level lvl) const
Check whether the logging system is active at the provided verbosity level.
Definition: MLogging.h:222
TFCSTruthState::set_pdgid
void set_pdgid(int val)
Definition: TFCSTruthState.h:18