GENFIT Rev: NoNumberAvailable
Loading...
Searching...
No Matches
MeasuredStateOnPlane.h
Go to the documentation of this file.
1/* Copyright 2008-2010, Technische Universitaet Muenchen,
2 Authors: Christian Hoeppner & Sebastian Neubert & Johannes Rauch
3
4 This file is part of GENFIT.
5
6 GENFIT is free software: you can redistribute it and/or modify
7 it under the terms of the GNU Lesser General Public License as published
8 by the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 GENFIT is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public License
17 along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18*/
19
23
24#ifndef genfit_MeasuredStateOnPlane_h
25#define genfit_MeasuredStateOnPlane_h
26
27#include "StateOnPlane.h"
28#include "AbsTrackRep.h"
29
30#include <TMatrixDSym.h>
31
32#include <cassert>
33
34namespace genfit {
35
40
41 public:
42
43 MeasuredStateOnPlane(const AbsTrackRep* rep = nullptr);
44 MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const genfit::SharedPlanePtr& plane, const AbsTrackRep* rep);
45 MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const genfit::SharedPlanePtr& plane, const AbsTrackRep* rep, const TVectorD& auxInfo);
47 MeasuredStateOnPlane(const StateOnPlane& state, const TMatrixDSym& cov);
48
50 void swap(MeasuredStateOnPlane& other); // nothrow
51
53 virtual MeasuredStateOnPlane* clone() const override {return new MeasuredStateOnPlane(*this);}
54
55
56 const TMatrixDSym& getCov() const {return cov_;}
57 TMatrixDSym& getCov() {return cov_;}
58
60 void blowUpCov(double blowUpFac, bool resetOffDiagonals = true, double maxVal = -1.);
61
62 void setStateCov(const TVectorD& state, const TMatrixDSym& cov) {setState(state); setCov(cov);}
63 void setStateCovPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane) {setStatePlane(state, plane); setCov(cov);}
64 void setCov(const TMatrixDSym& cov) {if(cov_.GetNrows() == 0) cov_.ResizeTo(cov); cov_ = cov;}
65
66 // Shortcuts to TrackRep functions
67 TMatrixDSym get6DCov() const {return getRep()->get6DCov(*this);};
68 void getPosMomCov(TVector3& pos, TVector3& mom, TMatrixDSym& cov) const {getRep()->getPosMomCov(*this, pos, mom, cov);}
69 void get6DStateCov(TVectorD& stateVec, TMatrixDSym& cov) const {getRep()->get6DStateCov(*this, stateVec, cov);}
70 double getMomVar() const {return getRep()->getMomVar(*this);}
71
72 void setPosMomErr(const TVector3& pos, const TVector3& mom, const TVector3& posErr, const TVector3& momErr) {getRep()->setPosMomErr(*this, pos, mom, posErr, momErr);}
73 void setPosMomCov(const TVector3& pos, const TVector3& mom, const TMatrixDSym& cov6x6) {getRep()->setPosMomCov(*this, pos, mom, cov6x6);}
74 void setPosMomCov(const TVectorD& state6, const TMatrixDSym& cov6x6) {getRep()->setPosMomCov(*this, state6, cov6x6);}
75
76
77 virtual void Print(Option_t* option = "") const override;
78
79 protected:
80
81 TMatrixDSym cov_;
82
83 public:
84 ClassDefOverride(MeasuredStateOnPlane,1)
85
86};
87
88
93
94
96 StateOnPlane::swap(other);
97 this->cov_.ResizeTo(other.cov_);
98 std::swap(this->cov_, other.cov_);
99}
100
102 StateOnPlane(rep), cov_(0,0)
103{
104 if (rep != nullptr) {
105 cov_.ResizeTo(rep->getDim(), rep->getDim());
106 }
107}
108
109inline MeasuredStateOnPlane::MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane, const AbsTrackRep* rep) :
110 StateOnPlane(state, plane, rep), cov_(cov)
111{
112 assert(rep != nullptr);
113 //assert(cov_.GetNcols() == (signed)rep->getDim());
114}
115
116inline MeasuredStateOnPlane::MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane, const AbsTrackRep* rep, const TVectorD& auxInfo) :
117 StateOnPlane(state, plane, rep, auxInfo), cov_(cov)
118{
119 assert(rep != nullptr);
120 //assert(cov_.GetNcols() == (signed)rep->getDim());
121}
122
127
128inline MeasuredStateOnPlane::MeasuredStateOnPlane(const StateOnPlane& state, const TMatrixDSym& cov) :
129 StateOnPlane(state), cov_(cov)
130{
131 //assert(cov_.GetNcols() == (signed)getRep()->getDim());
132}
133
135 swap(other);
136 return *this;
137}
138
139
140} /* End of namespace genfit */
142
143#endif // genfit_MeasuredStateOnPlane_h
Abstract base class for a track representation.
Definition AbsTrackRep.h:66
virtual double getMomVar(const MeasuredStateOnPlane &state) const =0
get the variance of the absolute value of the momentum .
virtual void setPosMomErr(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TVector3 &posErr, const TVector3 &momErr) const =0
Set position and momentum and error of state.
virtual void setPosMomCov(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const =0
Set position, momentum and covariance of state.
virtual unsigned int getDim() const =0
Get the dimension of the state vector used by the track representation.
virtual void get6DStateCov(const MeasuredStateOnPlane &state, TVectorD &stateVec, TMatrixDSym &cov) const
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
virtual TMatrixDSym get6DCov(const MeasuredStateOnPlane &state) const =0
Get the 6D covariance.
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
StateOnPlane with additional covariance matrix.
void getPosMomCov(TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const
MeasuredStateOnPlane & operator=(MeasuredStateOnPlane other)
void setStateCovPlane(const TVectorD &state, const TMatrixDSym &cov, const SharedPlanePtr &plane)
void setPosMomCov(const TVectorD &state6, const TMatrixDSym &cov6x6)
void setPosMomErr(const TVector3 &pos, const TVector3 &mom, const TVector3 &posErr, const TVector3 &momErr)
void setPosMomCov(const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6)
void setCov(const TMatrixDSym &cov)
MeasuredStateOnPlane(const AbsTrackRep *rep=nullptr)
const TMatrixDSym & getCov() const
void setStateCov(const TVectorD &state, const TMatrixDSym &cov)
void swap(MeasuredStateOnPlane &other)
void get6DStateCov(TVectorD &stateVec, TMatrixDSym &cov) const
virtual MeasuredStateOnPlane * clone() const override
void blowUpCov(double blowUpFac, bool resetOffDiagonals=true, double maxVal=-1.)
Blow up covariance matrix with blowUpFac. Per default, off diagonals are reset to 0 and the maximum v...
virtual void Print(Option_t *option="") const override
StateOnPlane(const genfit::StateOnPlane &)=default
void swap(StateOnPlane &other)
void setState(const TVectorD &state)
void setStatePlane(const TVectorD &state, const SharedPlanePtr &plane)
const AbsTrackRep * getRep() const
Defines for I/O streams used for error and debug printing.
MeasuredStateOnPlane calcAverageState(const MeasuredStateOnPlane &forwardState, const MeasuredStateOnPlane &backwardState)
Calculate weighted average between two MeasuredStateOnPlanes.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.