ATLAS Offline Software
Loading...
Searching...
No Matches
SimpleShape.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "TH1D.h"
7#include "TMath.h"
9
10
11
12using namespace LArSamples;
13
14
15SimpleShape::SimpleShape(const std::vector<double>& values, const std::vector<double>& errors,
16 double timeInterval, double startTime)
17 : m_values(values), m_errors(errors),
19{
20}
21
22
28
29
30SimpleShape::SimpleShape(const AbsShape& shape, double scale, double shift)
31{
32 for (unsigned int i = 0; i < shape.nPoints(); i++) {
33 m_values.push_back(shape.value(i)*scale);
34 m_errors.push_back(shape.error(i)*scale);
35 }
36
37 m_timeInterval = (shape.time(shape.nPoints() - 1) - shape.time(0))/(shape.nPoints() - 1);
38 m_startTime = shape.time(0) + shift;
39}
40
41
42SimpleShape::SimpleShape(const ShapeInfo& shapeInfo, double scale, double shift, bool samplingTimeOnly)
43{
44 for (unsigned int i = 0; i < shapeInfo.nPoints(); i++) {
45 if (samplingTimeOnly && shapeInfo.phase(i) != 0) continue;
46 m_values.push_back(shapeInfo.value(i)*scale);
47 m_errors.push_back(0);
48 }
49
50 double timeInterval = 1.0*Definitions::samplingInterval/(samplingTimeOnly ? 1 : shapeInfo.nIntervals());
52 m_startTime = shift - shapeInfo.shift();
53}
54
55
56double SimpleShape::time(unsigned int i) const
57{
58 return startTime() + i*timeInterval();
59}
60
61
63{
65 for (unsigned int i = 0; i < nPoints(); i++) {
66 double val;
67 interpolateDiff(time(i), val);
68 diff->set(i, val, 0);
69 }
70 return diff;
71}
72
73
74SimpleShape* SimpleShape::add(const AbsShape& other, double scale, double shift)
75{
76 std::vector<double> values, errors;
77 for (unsigned int k = 0; k < nPoints(); ++k) {
78 double val, err;
79 int pb = other.interpolate(time(k) - shift, val, err);
80 if (pb != 0) val = err = 0;
81 //if (pb == 1) break; //overflow: stop. underflow is OK : shape is assumed to be 0
82 values.push_back(this->value(k) + scale*val);
83 errors.push_back(TMath::Sqrt(TMath::Power(this->error(k),2) + TMath::Power(scale*err, 2)));
84 }
85
86 return new SimpleShape(values, errors, timeInterval(), startTime());
87}
88
89
90bool SimpleShape::add(unsigned int k, double val, double err)
91{
92 m_values[k] = m_values[k] + val;
93 m_errors [k] = TMath::Sqrt(m_errors[k]*m_errors[k] + err*err);
94 return true;
95}
96
97
102
103
104TH1D* SimpleShape::histogram(const char* name, const char* title, bool timeInUnitOfSamples) const
105{
106 double xMin = (timeInUnitOfSamples ? -1.5 : time(0) - 1.5*timeInterval());
107 double xMax = (timeInUnitOfSamples ? nPoints() + 0.5 : time(nPoints()) + 0.5*timeInterval());
108 TH1D* h = new TH1D(name, title, nPoints() + 2, xMin, xMax);
109 h->GetXaxis()->SetTitle(timeInUnitOfSamples ? "Sample Index" : "Time (ns)");
110 h->GetYaxis()->SetTitle("ADC counts");
111 h->GetYaxis()->SetTitleOffset(1.1);
112 h->SetMarkerStyle(20);
113 h->SetMarkerSize(1);
114
115 for (unsigned int i = 0; i < nPoints(); i++) {
116 h->SetBinContent(i + 2, value(i));
117 if (error(i) > 0) h->SetBinError(i + 2, error(i));
118 }
119
120 return h;
121}
122
123
124bool SimpleShape::add(std::unique_ptr<SimpleShape>& s1, const AbsShape& s2)
125{
126 if (!s1) {
127 s1.reset(new SimpleShape(s2));
128 return true;
129 }
130 SimpleShape* sum = s1->add(s2);
131 if (!sum) return false;
132 s1.reset(sum);
133 return true;
134}
135
136
137bool SimpleShape::scaleAndShift(std::unique_ptr<SimpleShape>& s1, double scale, double shift)
138{
139 if (!s1) return false;
140 s1.reset(new SimpleShape(*s1, scale, shift));
141 return true;
142}
virtual unsigned int nPoints() const =0
virtual double time(unsigned int i) const =0
TVectorD values(int lwb, int upb) const
Definition AbsShape.cxx:135
virtual double error(unsigned int i) const
Definition AbsShape.cxx:24
virtual double value(unsigned int i) const =0
int interpolateDiff(double time, double &diff) const
Definition AbsShape.cxx:108
unsigned char phase(unsigned int i) const
Definition ShapeInfo.cxx:66
float shift() const
Definition ShapeInfo.h:42
double value(unsigned int i) const
Definition ShapeInfo.cxx:58
unsigned int nPoints() const
Definition ShapeInfo.h:37
unsigned int nIntervals() const
Definition ShapeInfo.h:41
bool add(unsigned int k, double value, double error)
TH1D * histogram(const char *name="shape", const char *title="", bool timeInUnitOfSamples=false) const
unsigned int nPoints() const
Definition SimpleShape.h:48
SimpleShape(const std::vector< double > &values, const std::vector< double > &errors, double timeInterval=25, double startTime=0)
Constructor.
double timeInterval() const
Definition SimpleShape.h:45
SimpleShape * diff() const
std::vector< double > m_errors
Definition SimpleShape.h:68
double startTime() const
Definition SimpleShape.h:46
static bool scaleAndShift(std::unique_ptr< SimpleShape > &s1, double scale, double shift=0)
SimpleShape * createEmpty() const
double time(unsigned int i) const
std::vector< double > m_values
Definition SimpleShape.h:68
double value(unsigned int i) const
Definition SimpleShape.h:49