Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:18:21

0001 /* Copyright 2008-2010, Technische Universitaet Muenchen,
0002    Authors: Christian Hoeppner & Sebastian Neubert & Johannes Rauch
0003 
0004    This file is part of GENFIT.
0005 
0006    GENFIT is free software: you can redistribute it and/or modify
0007    it under the terms of the GNU Lesser General Public License as published
0008    by the Free Software Foundation, either version 3 of the License, or
0009    (at your option) any later version.
0010 
0011    GENFIT is distributed in the hope that it will be useful,
0012    but WITHOUT ANY WARRANTY; without even the implied warranty of
0013    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
0014    GNU Lesser General Public License for more details.
0015 
0016    You should have received a copy of the GNU Lesser General Public License
0017    along with GENFIT.  If not, see <http://www.gnu.org/licenses/>.
0018 */
0019 
0020 /** @addtogroup genfit
0021  * @{
0022  */
0023 
0024 #ifndef genfit_MeasuredStateOnPlane_h
0025 #define genfit_MeasuredStateOnPlane_h
0026 
0027 #include "StateOnPlane.h"
0028 #include "AbsTrackRep.h"
0029 
0030 #include <TMatrixDSym.h>
0031 
0032 #include <cassert>
0033 
0034 namespace genfit {
0035 
0036 /**
0037  *  @brief #StateOnPlane with additional covariance matrix.
0038  */
0039 class MeasuredStateOnPlane : public StateOnPlane {
0040 
0041  public:
0042 
0043   MeasuredStateOnPlane(const AbsTrackRep* rep = nullptr);
0044   MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const genfit::SharedPlanePtr& plane, const AbsTrackRep* rep);
0045   MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const genfit::SharedPlanePtr& plane, const AbsTrackRep* rep, const TVectorD& auxInfo);
0046   MeasuredStateOnPlane(const MeasuredStateOnPlane& o);
0047   MeasuredStateOnPlane(const StateOnPlane& state, const TMatrixDSym& cov);
0048 
0049   MeasuredStateOnPlane& operator=(MeasuredStateOnPlane other);
0050   void swap(MeasuredStateOnPlane& other); // nothrow
0051 
0052   virtual ~MeasuredStateOnPlane() {}
0053   virtual MeasuredStateOnPlane* clone() const override {return new MeasuredStateOnPlane(*this);}
0054 
0055 
0056   const TMatrixDSym& getCov() const {return cov_;}
0057   TMatrixDSym& getCov() {return cov_;}
0058 
0059   //! Blow up covariance matrix with blowUpFac. Per default, off diagonals are reset to 0 and the maximum values are limited to maxVal.
0060   void blowUpCov(double blowUpFac, bool resetOffDiagonals = true, double maxVal = -1.);
0061 
0062   void setStateCov(const TVectorD& state, const TMatrixDSym& cov) {setState(state); setCov(cov);}
0063   void setStateCovPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane) {setStatePlane(state, plane); setCov(cov);}
0064   void setCov(const TMatrixDSym& cov) {if(cov_.GetNrows() == 0) cov_.ResizeTo(cov); cov_ = cov;}
0065 
0066   // Shortcuts to TrackRep functions
0067   TMatrixDSym get6DCov() const {return getRep()->get6DCov(*this);};
0068   void getPosMomCov(TVector3& pos, TVector3& mom, TMatrixDSym& cov) const {getRep()->getPosMomCov(*this, pos, mom, cov);}
0069   void get6DStateCov(TVectorD& stateVec, TMatrixDSym& cov) const {getRep()->get6DStateCov(*this, stateVec, cov);}
0070   double getMomVar() const {return getRep()->getMomVar(*this);}
0071 
0072   void setPosMomErr(const TVector3& pos, const TVector3& mom, const TVector3& posErr, const TVector3& momErr) {getRep()->setPosMomErr(*this, pos, mom, posErr, momErr);}
0073   void setPosMomCov(const TVector3& pos, const TVector3& mom, const TMatrixDSym& cov6x6) {getRep()->setPosMomCov(*this, pos, mom, cov6x6);}
0074   void setPosMomCov(const TVectorD& state6, const TMatrixDSym& cov6x6) {getRep()->setPosMomCov(*this, state6, cov6x6);}
0075 
0076 
0077   virtual void Print(Option_t* option = "") const override;
0078 
0079  protected:
0080 
0081   TMatrixDSym cov_;
0082 
0083  public:
0084   ClassDefOverride(MeasuredStateOnPlane,1)
0085 
0086 };
0087 
0088 
0089 /**
0090  * @brief Calculate weighted average between two MeasuredStateOnPlanes
0091  */
0092 MeasuredStateOnPlane calcAverageState(const MeasuredStateOnPlane& forwardState, const MeasuredStateOnPlane& backwardState);
0093 
0094 
0095 inline void MeasuredStateOnPlane::swap(MeasuredStateOnPlane& other) {
0096   StateOnPlane::swap(other);
0097   this->cov_.ResizeTo(other.cov_);
0098   std::swap(this->cov_, other.cov_);
0099 }
0100 
0101 inline MeasuredStateOnPlane::MeasuredStateOnPlane(const AbsTrackRep* rep) :
0102   StateOnPlane(rep), cov_(0,0)
0103 {
0104   if (rep != nullptr) {
0105     cov_.ResizeTo(rep->getDim(), rep->getDim());
0106   }
0107 }
0108 
0109 inline MeasuredStateOnPlane::MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane, const AbsTrackRep* rep) :
0110   StateOnPlane(state, plane, rep), cov_(cov)
0111 {
0112   assert(rep != nullptr);
0113   //assert(cov_.GetNcols() == (signed)rep->getDim());
0114 }
0115 
0116 inline MeasuredStateOnPlane::MeasuredStateOnPlane(const TVectorD& state, const TMatrixDSym& cov, const SharedPlanePtr& plane, const AbsTrackRep* rep, const TVectorD& auxInfo) :
0117   StateOnPlane(state, plane, rep, auxInfo), cov_(cov)
0118 {
0119   assert(rep != nullptr);
0120   //assert(cov_.GetNcols() == (signed)rep->getDim());
0121 }
0122 
0123 inline MeasuredStateOnPlane::MeasuredStateOnPlane(const MeasuredStateOnPlane& o) :
0124   StateOnPlane(o), cov_(o.cov_)
0125 {
0126 }
0127 
0128 inline MeasuredStateOnPlane::MeasuredStateOnPlane(const StateOnPlane& state, const TMatrixDSym& cov) :
0129   StateOnPlane(state), cov_(cov)
0130 {
0131   //assert(cov_.GetNcols() == (signed)getRep()->getDim());
0132 }
0133 
0134 inline MeasuredStateOnPlane& MeasuredStateOnPlane::operator=(MeasuredStateOnPlane other) {
0135   swap(other);
0136   return *this;
0137 }
0138 
0139 
0140 } /* End of namespace genfit */
0141 /** @} */
0142 
0143 #endif // genfit_MeasuredStateOnPlane_h