CMS 3D CMS Logo

Public Member Functions | Private Types | Private Member Functions | Private Attributes

BasicMultiTrajectoryState Class Reference

#include <BasicMultiTrajectoryState.h>

Inheritance diagram for BasicMultiTrajectoryState:
BasicTrajectoryState ReferenceCountedPoolAllocated BlockWipedPoolAllocated

List of all members.

Public Member Functions

 BasicMultiTrajectoryState (const std::vector< TSOS > &tsvec)
 BasicMultiTrajectoryState ()
virtual bool canUpdateLocalParameters () const
virtual BasicMultiTrajectoryStateclone () const
virtual std::vector
< TrajectoryStateOnSurface
components () const
void rescaleError (double factor)
virtual void update (const LocalTrajectoryParameters &p, const LocalTrajectoryError &err, const Surface &aSurface, const MagneticField *field, const SurfaceSide side, double weight)
virtual void update (const LocalTrajectoryParameters &p, const Surface &aSurface, const MagneticField *field, const SurfaceSide side)

Private Types

typedef TrajectoryStateOnSurface TSOS

Private Member Functions

void combine ()

Private Attributes

std::vector< TSOStheStates

Detailed Description

Class which combines a set of components of a Gaussian mixture into a single component. Given all the components of a mixture, it calculates the mean and covariance matrix of the entire mixture. This combiner class can also be used in the process of transforming a Gaussian mixture into another Gaussian mixture with a smaller number of components. The relevant formulas can be found in R. Fruhwirth, Computer Physics Communications 100 (1997), 1.

Definition at line 17 of file BasicMultiTrajectoryState.h.


Member Typedef Documentation

Definition at line 19 of file BasicMultiTrajectoryState.h.


Constructor & Destructor Documentation

BasicMultiTrajectoryState::BasicMultiTrajectoryState ( const std::vector< TSOS > &  tsvec)

Definition at line 8 of file BasicMultiTrajectoryState.cc.

References combine(), Exception, i, theStates, and unlikely.

                                                                                  :
  BasicTrajectoryState(tsvec.front().surface()) {
 
  // only accept planes!!
  const BoundPlane* bp = dynamic_cast<const BoundPlane*>(&tsvec.begin()->surface());
  if unlikely( bp==0 )
               throw cms::Exception("LogicError") << "MultiTrajectoryState constructed on cylinder";
   
  theStates.reserve(tsvec.size());
  for (std::vector<TSOS>::const_iterator i=tsvec.begin(); i!=tsvec.end(); i++) {
    if unlikely(!i->isValid()) {
      throw cms::Exception("LogicError") << "MultiTrajectoryState constructed with invalid state";
    }
    if unlikely(i->hasError() != tsvec.front().hasError()) {
      throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with and without errors";
    }
    if unlikely( &i->surface() != &tsvec.front().surface()) {
      throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with different surfaces";
    }
    if unlikely( i->surfaceSide() != tsvec.front().surfaceSide()) {
      throw cms::Exception("LogicError") 
        << "MultiTrajectoryState mixes states defined before and after material";
    }
    if unlikely( i->localParameters().pzSign()*tsvec.front().localParameters().pzSign()<0. ) {
      throw cms::Exception("LogicError") 
        << "MultiTrajectoryState mixes states with different signs of local p_z";
    }

    theStates.push_back( *i);
  }
  combine();
}
BasicMultiTrajectoryState::BasicMultiTrajectoryState ( ) [inline]

Definition at line 25 of file BasicMultiTrajectoryState.h.

Referenced by clone().

{}

Member Function Documentation

virtual bool BasicMultiTrajectoryState::canUpdateLocalParameters ( ) const [inline, virtual]

Reimplemented from BasicTrajectoryState.

Definition at line 45 of file BasicMultiTrajectoryState.h.

{ return false; }
virtual BasicMultiTrajectoryState* BasicMultiTrajectoryState::clone ( void  ) const [inline, virtual]

Implements BasicTrajectoryState.

Definition at line 36 of file BasicMultiTrajectoryState.h.

References BasicMultiTrajectoryState().

                                                   {
    return new BasicMultiTrajectoryState(*this);
  }
void BasicMultiTrajectoryState::combine ( ) [private]

Definition at line 57 of file BasicMultiTrajectoryState.cc.

References diffTreeTool::diff, timingPdfMaker::mean, alignCSCRings::s, theStates, unlikely, update(), and BasicTrajectoryState::weight().

Referenced by BasicMultiTrajectoryState(), and rescaleError().

                                    {
  const std::vector<TrajectoryStateOnSurface>& tsos = theStates;

  if unlikely(tsos.empty()) {
    edm::LogError("MultiTrajectoryStateCombiner") 
      << "Trying to collapse empty set of trajectory states!";
    return;
  }

  double pzSign = tsos.front().localParameters().pzSign();
  for (std::vector<TrajectoryStateOnSurface>::const_iterator it = tsos.begin(); 
       it != tsos.end(); it++) {
    if unlikely(it->localParameters().pzSign() != pzSign) {
      edm::LogError("MultiTrajectoryStateCombiner") 
        << "Trying to collapse trajectory states with different signs on p_z!";
      return;
    }
  }
  
  if unlikely(tsos.size() == 1) {
    BasicTrajectoryState::update(tsos.front().localParameters(), 
                                 tsos.front().localError(), 
                                 tsos.front().surface(), 
                                 tsos.front().magneticField(),
                                 tsos.front().surfaceSide(), 
                                 tsos.front().weight()
                                 );
    return;
  }
  
  double sumw = 0.;
  //int dim = tsos.front().localParameters().vector().num_row();
  AlgebraicVector5 mean;
  AlgebraicSymMatrix55 covarPart1, covarPart2;
  for (std::vector<TrajectoryStateOnSurface>::const_iterator it1 = tsos.begin(); 
       it1 != tsos.end(); it1++) {
    double weight = it1->weight();
    AlgebraicVector5 param = it1->localParameters().vector();
    sumw += weight;
    mean += weight * param;
    covarPart1 += weight * it1->localError().matrix();
    for (std::vector<TrajectoryStateOnSurface>::const_iterator it2 = it1 + 1; 
         it2 != tsos.end(); it2++) {
      AlgebraicVector5 diff = param - it2->localParameters().vector();
      AlgebraicSymMatrix11 s = AlgebraicMatrixID(); //stupid trick to make CLHEP work decently
      covarPart2 += (weight * it2->weight()) * 
                                ROOT::Math::Similarity(AlgebraicMatrix51(diff.Array(), 5), s);
                        //FIXME: we can surely write this thing in a better way
    }   
  }
  double sumwI = 1.0/sumw;
  mean *= sumwI;
  covarPart1 *= sumwI; covarPart2 *= (sumwI*sumwI);
  AlgebraicSymMatrix55 covar = covarPart1 + covarPart2;

  BasicTrajectoryState::update(LocalTrajectoryParameters(mean, pzSign), 
                               LocalTrajectoryError(covar), 
                               tsos.front().surface(), 
                               tsos.front().magneticField(),
                               tsos.front().surfaceSide(), 
                               sumw);
}
virtual std::vector<TrajectoryStateOnSurface> BasicMultiTrajectoryState::components ( ) const [inline, virtual]

Reimplemented from BasicTrajectoryState.

Definition at line 40 of file BasicMultiTrajectoryState.h.

References theStates.

                                                                 {
    return theStates;
  }
void BasicMultiTrajectoryState::rescaleError ( double  factor)

Rescaling the error of the mixture with a given factor. Please note that this rescaling is imposed on each of the components of the mixture and does therefore not exactly correspond to rescaling theCombinedState with the same factor.

Reimplemented from BasicTrajectoryState.

Definition at line 43 of file BasicMultiTrajectoryState.cc.

References combine(), theStates, and unlikely.

                                                          {

  if unlikely(theStates.empty()) {
    edm::LogError("BasicMultiTrajectoryState") << "Trying to rescale errors of empty MultiTrajectoryState!";
    return;
  }
  
  for (std::vector<TSOS>::iterator it = theStates.begin(); it != theStates.end(); it++) {
    it->rescaleError(factor);
  }
  combine();
}
void BasicMultiTrajectoryState::update ( const LocalTrajectoryParameters p,
const LocalTrajectoryError err,
const Surface aSurface,
const MagneticField field,
const SurfaceSide  side,
double  weight 
) [virtual]

Reimplemented from BasicTrajectoryState.

Definition at line 133 of file BasicMultiTrajectoryState.cc.

References Exception.

{
  throw cms::Exception("LogicError", 
                       "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, LocalTrajectoryError, ...) called even if canUpdateLocalParameters() is false");
}
void BasicMultiTrajectoryState::update ( const LocalTrajectoryParameters p,
const Surface aSurface,
const MagneticField field,
const SurfaceSide  side 
) [virtual]

Reimplemented from BasicTrajectoryState.

Definition at line 122 of file BasicMultiTrajectoryState.cc.

References Exception.

Referenced by combine().

{
  throw cms::Exception("LogicError", 
                       "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, Surface, ...) called even if canUpdateLocalParameters() is false");
}

Member Data Documentation