#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]
Definition at line 198 of file AbsShape.cxx.
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);
◆ covarianceMatrix() [2/2]
CovMatrix AbsShape::covarianceMatrix |
( |
unsigned int |
nPoints, |
|
|
bool |
withCorrs = true |
|
) |
| const |
◆ error()
double AbsShape::error |
( |
unsigned int |
i | ) |
const |
|
virtual |
◆ findTimeInterval()
int AbsShape::findTimeInterval |
( |
double |
time | ) |
const |
◆ graph()
TGraphErrors * AbsShape::graph |
( |
bool |
timeInUnitOfSamples = false | ) |
const |
◆ 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.
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; }
162 values.ResizeTo(actualLwb, actualUpb);
163 errors.ResizeTo(actualLwb, actualUpb, actualLwb, actualUpb);
164 return (actualLwb >= 0);
◆ interpolate() [2/2]
int AbsShape::interpolate |
( |
double |
time, |
|
|
double & |
value, |
|
|
double & |
error |
|
) |
| 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.
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; }
181 diffs.ResizeTo(actualLwb, actualUpb);
182 return (actualLwb >= 0);
◆ interpolateDiff() [2/2]
int AbsShape::interpolateDiff |
( |
double |
time, |
|
|
double & |
diff |
|
) |
| const |
◆ invCovarianceMatrix() [1/2]
◆ invCovarianceMatrix() [2/2]
CovMatrix AbsShape::invCovarianceMatrix |
( |
unsigned int |
nPoints, |
|
|
bool |
withCorrs = true |
|
) |
| const |
◆ maxPosition()
int AbsShape::maxPosition |
( |
| ) |
const |
◆ maxValue()
double AbsShape::maxValue |
( |
bool |
withErrors = false | ) |
const |
- Returns
- sample max parameters
Definition at line 30 of file AbsShape.cxx.
◆ minPosition()
int AbsShape::minPosition |
( |
| ) |
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.
242 double dT = (
t1 -
t0)/nPts;
244 for (
unsigned int i = 0;
i < nPts;
i++,
t += dT) {
247 if (
inRange != 0)
return nullptr;
◆ 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.
137 if (lwb < 0) lwb = 0;
140 TVectorD
vals(lwb, upb);
141 for (
int i = lwb;
i <= upb;
i++)
The documentation for this class was generated from the following files: