Skip to content
Snippets Groups Projects
Commit ed01ad24 authored by Oliver Sander's avatar Oliver Sander Committed by sander@PCPOOL.MI.FU-BERLIN.DE
Browse files

make the target space a template parameter

[[Imported from SVN: r4046]]
parent bc319072
No related branches found
No related tags found
No related merge requests found
#ifndef ROD_DIFFERENCE_HH #ifndef GEODESIC_DIFFERENCE_HH
#define ROD_DIFFERENCE_HH #define GEODESIC_DIFFERENCE_HH
#include "rigidbodymotion.hh" template <class TargetSpace>
Dune::BlockVector<typename TargetSpace::TangentVector> computeGeodesicDifference(const std::vector<TargetSpace>& a,
Dune::BlockVector<Dune::FieldVector<double,6> > computeGeodesicDifference(const std::vector<RigidBodyMotion<3> >& a, const std::vector<TargetSpace>& b)
const std::vector<RigidBodyMotion<3> >& b)
{
if (a.size() != b.size())
DUNE_THROW(Dune::Exception, "a and b have to have the same length!");
Dune::BlockVector<Dune::FieldVector<double,6> > result(a.size());
for (size_t i=0; i<result.size(); i++) {
// Subtract centerline position
for (int j=0; j<3; j++)
result[i][j] = a[i].r[j] - b[i].r[j];
// Subtract orientations on the tangent space of 'a'
Dune::FieldVector<double,3> v = Rotation<3,double>::difference(a[i].q, b[i].q);
// Compute difference on T_a SO(3)
for (int j=0; j<3; j++)
result[i][j+3] = v[j];
}
return result;
}
Dune::BlockVector<Dune::FieldVector<double,3> > computeGeodesicDifference(const std::vector<Rotation<3,double> >& a,
const std::vector<Rotation<3,double> >& b)
{ {
if (a.size() != b.size()) if (a.size() != b.size())
DUNE_THROW(Dune::Exception, "a and b have to have the same length!"); DUNE_THROW(Dune::Exception, "a and b have to have the same length!");
Dune::BlockVector<Dune::FieldVector<double,3> > result(a.size()); Dune::BlockVector<typename TargetSpace::TangentVector> result(a.size());
for (size_t i=0; i<result.size(); i++) { for (size_t i=0; i<result.size(); i++) {
// Subtract orientations on the tangent space of 'a' // Subtract orientations on the tangent space of 'a'
result[i] = Rotation<3,double>::difference(a[i], b[i]); result[i] = TargetSpace::difference(a[i], b[i]);
} }
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment