ATLAS Offline Software
Loading...
Searching...
No Matches
EfieldInterpolator.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3 */
15
16#ifndef PIXELDIGITIZATION_EFIELDINTERPOLATOR_H
17#define PIXELDIGITIZATION_EFIELDINTERPOLATOR_H
18
20#include <vector>
21#include <string>
22#include <string_view>
23
24class TH1;
25class TH1D;
26class TString;
27
28// Different options to retrieve E field
29// default: interspline
31 TCAD, //use Precomputed TCAD map
32 interspline, //use interpolation based on cubic splines
33 interinvdist, //use inverse distance weighted estimate
34 linearField //linear field according to bias voltage
35};
36
38public:
39 EfieldInterpolator(const std::string& type, const std::string& name, const IInterface* parent);
41 void setLayer(int layer);
42 //Recommended constructor
43 StatusCode loadTCADlist(const std::string& TCADfileListToLoad);
44 //defFct
45
46 virtual StatusCode initialize() override;
47 virtual StatusCode finalize() override;
48
49 // Member Functions
50 const std::string loadTCADfiles(const std::string& targetList = "");
51 const std::string createInterpolationFromTCADtree(const std::string& fTCAD);//TTree* tTCAD);
52 bool initializeFromFile(const std::string& finpath);
53 bool initializeFromDirectory(const std::string& fpath);
54 double estimateEfield(std::vector<double> vvol, const std::vector<double>& vflu, const std::vector<std::vector<double> >& vfluvvol,
55 double aimFlu, double aimVol, std::string_view prepend = "", bool debug = false);
56 double estimateEfieldInvDistance(const std::vector<double> & vvol, const std::vector<double> & vflu,
57 const std::vector<std::vector<double> > & vfluvvol, double aimFlu, double aimVol,
58 double measure = 1.);
59
60 TH1D* createEfieldProfile(double aimFluence, double aimVoltage);
61 TH1D* getEfield(double aimFluence, double aimVoltage);
62 TH1D* loadEfieldFromDat(const std::string& fname, bool fillEdges = true);
63 static void scaleIntegralTo(TH1* hin, double aimInt, int first = 1, int last = -1);
64 void reliabilityCheck(double aimFluence, const std::vector<double>& fluences, double aimVoltage,
65 const std::vector<double>& voltages);
66private:
67 // Member variables
68 Gaudi::Property<bool> m_initialized
69 {
70 this, "initialized", false, "Flag whether already initalized"
71 };
72
73 Gaudi::Property<bool> m_useSpline
74 {
75 this, "useSpline", true, "Flag whether to use cubic splines for interpolation"
76 };
77
78 Gaudi::Property<int> m_sensorDepth
79 {
80 this, "sensorDepth", 200, "Depth of E fields in sensors in um"
81 };
82
84 TH1D* m_efieldProfile{}; //Final efield profile
85 std::string m_fInter; //path to .root file for saving interpolation TTree, i.e. ordered by pixeldepth z
86 std::vector<std::vector<TString> > list_files(const TString& fileList_TCADsamples);
87 static double extrapolateLinear(double x1, double y1, double x2, double y2, double xaim);
88 int fillXYvectors(const std::vector<double> &vLoop, int ifix, const std::vector<std::vector<double> > & v2vsv1,
89 std::vector<double>& xx, std::vector<double>& yy, bool regularOrder = true);
90 void fillEdgeValues(TH1D* hin);
91 bool isInterpolation(const std::vector<double>& vval, double aimval)
92 {return(vval.front() <= aimval && aimval <= vval.back());};
93 bool isInterpolation(std::vector<double>* vval, double aimval)
94 {return(vval->front() <= aimval && aimval <= vval->back());};
95 static double relativeDistance(double x1, double x2); //difference between x1 x2 scaled to x1
96 static double relativeDistance(double x1, double y1, double x2, double y2);
97 double estimateEfieldLinear(double aimVoltage);
98 void saveTGraph(std::vector<double> vvol, std::vector<double> vflu, std::vector<std::vector<double> > vfluvvol,
99 double aimFlu, double aimVol, std::string_view prepend, bool skipNegative = true);
100};
101
102#endif //> !PIXELDIGITIZATION_EFIELDINTERPOLATOR_H
interpolationMethod
@ interspline
@ interinvdist
@ linearField
const bool debug
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Property< bool > m_initialized
interpolationMethod m_efieldOrigin
virtual StatusCode finalize() override
Gaudi::Property< int > m_sensorDepth
EfieldInterpolator(const std::string &type, const std::string &name, const IInterface *parent)
TH1D * loadEfieldFromDat(const std::string &fname, bool fillEdges=true)
double estimateEfield(std::vector< double > vvol, const std::vector< double > &vflu, const std::vector< std::vector< double > > &vfluvvol, double aimFlu, double aimVol, std::string_view prepend="", bool debug=false)
virtual StatusCode initialize() override
bool initializeFromFile(const std::string &finpath)
std::vector< std::vector< TString > > list_files(const TString &fileList_TCADsamples)
const std::string createInterpolationFromTCADtree(const std::string &fTCAD)
static double relativeDistance(double x1, double x2)
bool initializeFromDirectory(const std::string &fpath)
static double extrapolateLinear(double x1, double y1, double x2, double y2, double xaim)
StatusCode loadTCADlist(const std::string &TCADfileListToLoad)
virtual ~EfieldInterpolator()
const std::string loadTCADfiles(const std::string &targetList="")
int fillXYvectors(const std::vector< double > &vLoop, int ifix, const std::vector< std::vector< double > > &v2vsv1, std::vector< double > &xx, std::vector< double > &yy, bool regularOrder=true)
static void scaleIntegralTo(TH1 *hin, double aimInt, int first=1, int last=-1)
void fillEdgeValues(TH1D *hin)
TH1D * getEfield(double aimFluence, double aimVoltage)
void reliabilityCheck(double aimFluence, const std::vector< double > &fluences, double aimVoltage, const std::vector< double > &voltages)
TH1D * createEfieldProfile(double aimFluence, double aimVoltage)
bool isInterpolation(const std::vector< double > &vval, double aimval)
void saveTGraph(std::vector< double > vvol, std::vector< double > vflu, std::vector< std::vector< double > > vfluvvol, double aimFlu, double aimVol, std::string_view prepend, bool skipNegative=true)
Gaudi::Property< bool > m_useSpline
bool isInterpolation(std::vector< double > *vval, double aimval)
double estimateEfieldInvDistance(const std::vector< double > &vvol, const std::vector< double > &vflu, const std::vector< std::vector< double > > &vfluvvol, double aimFlu, double aimVol, double measure=1.)
double estimateEfieldLinear(double aimVoltage)