#include <AbsShape.h>
|
| virtual | ~AbsShape () |
| virtual unsigned int | nPoints () const =0 |
| virtual double | value (unsigned int i) const =0 |
| virtual double | covariance (unsigned int i, unsigned int j) const =0 |
| virtual double | time (unsigned int i) const =0 |
| virtual double | error (unsigned int i) const |
| TVectorD | values (int lwb, int upb) const |
| int | findTimeInterval (double time) const |
| int | interpolate (double time, double &value, double &error) const |
| int | interpolateDiff (double time, double &diff) const |
| bool | interpolate (const AbsShape &other, TVectorD &values, CovMatrix &errors, int lwb=-1, int upb=-1) const |
| bool | interpolateDiff (const AbsShape &other, TVectorD &diffs, int lwb=-1, int upb=-1) const |
| CovMatrix | covarianceMatrix (int lwb=-1, int upb=-1, const CovMatrix &refErr=CovMatrix(), bool withCorrs=true) const |
| CovMatrix | invCovarianceMatrix (int lwb=-1, int upb=-1, const CovMatrix &refErr=CovMatrix(), bool withCorrs=true) const |
| CovMatrix | covarianceMatrix (unsigned int nPoints, bool withCorrs=true) const |
| CovMatrix | invCovarianceMatrix (unsigned int nPoints, bool withCorrs=true) const |
| double | maxValue (bool withErrors=false) const |
| double | minValue (bool withErrors=false) const |
| int | maxPosition () const |
| int | minPosition () const |
| TGraphErrors * | graph (bool timeInUnitOfSamples=false) const |
| SimpleShape * | resample (unsigned int nPts) const |
Definition at line 28 of file AbsShape.h.
◆ ~AbsShape()
| virtual LArSamples::AbsShape::~AbsShape |
( |
| ) |
|
|
inlinevirtual |
◆ AbsShape()
| LArSamples::AbsShape::AbsShape |
( |
| ) |
|
|
inlineprotected |
◆ covariance()
| virtual double LArSamples::AbsShape::covariance |
( |
unsigned int | i, |
|
|
unsigned int | j ) const |
|
pure virtual |
◆ covarianceMatrix() [1/2]
| CovMatrix AbsShape::covarianceMatrix |
( |
int | lwb = -1, |
|
|
int | upb = -1, |
|
|
const CovMatrix & | refErr = CovMatrix(), |
|
|
bool | withCorrs = true ) const |
Definition at line 198 of file AbsShape.cxx.
201{
202 if (lwb < 0 || upb < 0) { lwb = 0; upb =
nPoints() - 1; }
204 for (
int i = lwb;
i <= upb;
i++) {
205 for (int j = lwb; j <= upb; j++) {
207 if (refErr.GetNrows() > 0) cov += refErr(i, j);
209 }
210 }
212}
virtual unsigned int nPoints() const =0
virtual double covariance(unsigned int i, unsigned int j) const =0
TMatrixTSym< double > CovMatrix
◆ covarianceMatrix() [2/2]
| CovMatrix AbsShape::covarianceMatrix |
( |
unsigned int | nPoints, |
|
|
bool | withCorrs = true ) const |
Definition at line 223 of file AbsShape.cxx.
225{
227}
CovMatrix covarianceMatrix(int lwb=-1, int upb=-1, const CovMatrix &refErr=CovMatrix(), bool withCorrs=true) const
◆ error()
| double AbsShape::error |
( |
unsigned int | i | ) |
const |
|
virtual |
◆ findTimeInterval()
| int AbsShape::findTimeInterval |
( |
double | time | ) |
const |
Definition at line 71 of file AbsShape.cxx.
72{
73 if (t <
time(0))
return -1;
75
76 for (
unsigned int i = 0;
i <
nPoints() - 1;
i++)
77 if (
time(i) <= t && t <
time(i+1))
return i;
78
80}
virtual double time(unsigned int i) const =0
◆ graph()
| TGraphErrors * AbsShape::graph |
( |
bool | timeInUnitOfSamples = false | ) |
const |
Definition at line 186 of file AbsShape.cxx.
187{
190 for (
unsigned int i = 0;
i <
nPoints();
i++) {
193 }
195}
const PlainObject unit() const
This is a plugin that makes Eigen look like CLHEP & defines some convenience methods.
virtual double error(unsigned int i) const
virtual double value(unsigned int i) const =0
TGraphErrors * graph(bool timeInUnitOfSamples=false) const
static const unsigned int samplingInterval
◆ interpolate() [1/2]
| bool AbsShape::interpolate |
( |
const AbsShape & | other, |
|
|
TVectorD & | values, |
|
|
CovMatrix & | errors, |
|
|
int | lwb = -1, |
|
|
int | upb = -1 ) const |
Definition at line 147 of file AbsShape.cxx.
148{
149 if (lwb < 0) lwb = 0;
150 if (upb < 0) upb =
other.nPoints() - 1;
151 values.ResizeTo(lwb, upb);
152 errors.ResizeTo(lwb, upb, lwb, upb);
153 int actualLwb = -1, actualUpb = upb;
154 for (
int i = lwb;
i <= upb;
i++) {
157 if (stat == 0 && actualLwb < 0) actualLwb =
i;
158 if (stat == +1) { actualUpb =
i - 1;
break; }
161 }
162 values.ResizeTo(actualLwb, actualUpb);
163 errors.ResizeTo(actualLwb, actualUpb, actualLwb, actualUpb);
164 return (actualLwb >= 0);
165}
TVectorD values(int lwb, int upb) const
int interpolate(double time, double &value, double &error) const
◆ interpolate() [2/2]
| int AbsShape::interpolate |
( |
double | time, |
|
|
double & | value, |
|
|
double & | error ) const |
Definition at line 83 of file AbsShape.cxx.
84{
86
89
90 if (i1 < 0) {
92 i1 = 0;
93 }
94
98 }
99
103
105}
int findTimeInterval(double time) const
◆ interpolateDiff() [1/2]
| bool AbsShape::interpolateDiff |
( |
const AbsShape & | other, |
|
|
TVectorD & | diffs, |
|
|
int | lwb = -1, |
|
|
int | upb = -1 ) const |
Definition at line 168 of file AbsShape.cxx.
169{
170 if (lwb < 0) lwb = 0;
171 if (upb < 0) upb =
other.nPoints() - 1;
172 diffs.ResizeTo(lwb, upb);
173 int actualLwb = -1, actualUpb = upb;
174 for (
int i = lwb;
i <= upb;
i++) {
177 if (stat == 0 && actualLwb < 0) actualLwb =
i;
178 if (stat == +1) { actualUpb =
i - 1;
break; }
180 }
181 diffs.ResizeTo(actualLwb, actualUpb);
182 return (actualLwb >= 0);
183}
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
int interpolateDiff(double time, double &diff) const
◆ interpolateDiff() [2/2]
| int AbsShape::interpolateDiff |
( |
double | time, |
|
|
double & | diff ) const |
Definition at line 108 of file AbsShape.cxx.
109{
112 if (i1 < 0) {
113 i1 = 0;
115 }
119 }
120 if (i1 > 1 && i1 <
int(
nPoints()) - 2)
122 else
123 if (i1 > 0)
126 else
128 else
130
132}
◆ invCovarianceMatrix() [1/2]
| CovMatrix AbsShape::invCovarianceMatrix |
( |
int | lwb = -1, |
|
|
int | upb = -1, |
|
|
const CovMatrix & | refErr = CovMatrix(), |
|
|
bool | withCorrs = true ) const |
◆ invCovarianceMatrix() [2/2]
| CovMatrix AbsShape::invCovarianceMatrix |
( |
unsigned int | nPoints, |
|
|
bool | withCorrs = true ) const |
◆ maxPosition()
| int AbsShape::maxPosition |
( |
| ) |
const |
Definition at line 52 of file AbsShape.cxx.
53{
56 for (
unsigned int i = 0;
i <
nPoints();
i++)
59}
double maxValue(bool withErrors=false) const
◆ maxValue()
| double AbsShape::maxValue |
( |
bool | withErrors = false | ) |
const |
- Returns
- sample max parameters
Definition at line 30 of file AbsShape.cxx.
31{
33 for (
unsigned int i = 0;
i <
nPoints();
i++) {
36 }
38}
◆ minPosition()
| int AbsShape::minPosition |
( |
| ) |
const |
Definition at line 61 of file AbsShape.cxx.
62{
65 for (
unsigned int i = 0;
i <
nPoints();
i++)
68}
double minValue(bool withErrors=false) const
◆ minValue()
| double AbsShape::minValue |
( |
bool | withErrors = false | ) |
const |
◆ nPoints()
| virtual unsigned int LArSamples::AbsShape::nPoints |
( |
| ) |
const |
|
pure virtual |
◆ resample()
| SimpleShape * AbsShape::resample |
( |
unsigned int | nPts | ) |
const |
Definition at line 237 of file AbsShape.cxx.
238{
242 double dT = (
t1 -
t0)/nPts;
244 for (
unsigned int i = 0;
i < nPts;
i++,
t += dT) {
247 if (
inRange != 0)
return nullptr;
250 }
251
252 return new SimpleShape(
values, errors, dT,
t0);
253}
bool inRange(const double *boundaries, const double value, const double tolerance=0.02)
std::vector< ALFA_RawDataCollection_p1 > t1
◆ time()
| virtual double LArSamples::AbsShape::time |
( |
unsigned int | i | ) |
const |
|
pure virtual |
◆ value()
| virtual double LArSamples::AbsShape::value |
( |
unsigned int | i | ) |
const |
|
pure virtual |
◆ values()
| TVectorD AbsShape::values |
( |
int | lwb, |
|
|
int | upb ) const |
Definition at line 135 of file AbsShape.cxx.
136{
137 if (lwb < 0) lwb = 0;
140 TVectorD
vals(lwb, upb);
141 for (
int i = lwb;
i <= upb;
i++)
144}
The documentation for this class was generated from the following files: