Skip to content
Snippets Groups Projects
Commit 24faf399 authored by Carsten Gräser's avatar Carsten Gräser
Browse files

Implement jacobian() and jacobianInverse() in (Cached)MultiLinearGeometry

parent bf5e28be
No related branches found
No related tags found
1 merge request!193Extend jacobian interface of geometry classes
......@@ -201,6 +201,12 @@ namespace Dune
//! type of jacobian inverse transposed
class JacobianInverseTransposed;
/** \brief Type for the Jacobian matrix */
typedef FieldMatrix< ctype, coorddimension, mydimension > Jacobian;
/** \brief Type for the inverse Jacobian matrix */
typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianInverse;
protected:
typedef Dune::ReferenceElements< ctype, mydimension > ReferenceElements;
......@@ -389,6 +395,29 @@ namespace Dune
return geometry.refElement();
}
/** \brief Obtain the Jacobian
*
* \param[in] local local coordinate to evaluate Jacobian in
*
* \returns a copy of the transposed of the Jacobian
*/
Jacobian jacobian (const LocalCoordinate &local) const
{
return jacobianTransposed(local).transposed();
}
/** \brief Obtain the Jacobian's inverse
*
* The Jacobian's inverse is defined as a pseudo-inverse. If we denote
* the Jacobian by \f$J(x)\f$, the following condition holds:
* \f[J^{-1}(x) J(x) = I.\f]
*/
JacobianInverse jacobianInverse (const LocalCoordinate &local) const
{
return jacobianInverseTransposed(local).transposed();
}
protected:
ReferenceElement refElement () const
......@@ -513,6 +542,8 @@ namespace Dune
typedef typename Base::JacobianTransposed JacobianTransposed;
typedef typename Base::JacobianInverseTransposed JacobianInverseTransposed;
typedef typename Base::Jacobian Jacobian;
typedef typename Base::JacobianInverse JacobianInverse;
template< class CornerStorage >
CachedMultiLinearGeometry ( const ReferenceElement &referenceElement, const CornerStorage &cornerStorage )
......@@ -660,6 +691,28 @@ namespace Dune
return Base::jacobianInverseTransposed( local );
}
/** \brief Obtain the Jacobian
*
* \param[in] local local coordinate to evaluate Jacobian in
*
* \returns a copy of the transposed of the Jacobian
*/
Jacobian jacobian (const LocalCoordinate &local) const
{
return jacobianTransposed(local).transposed();
}
/** \brief Obtain the Jacobian's inverse
*
* The Jacobian's inverse is defined as a pseudo-inverse. If we denote
* the Jacobian by \f$J(x)\f$, the following condition holds:
* \f[J^{-1}(x) J(x) = I.\f]
*/
JacobianInverse jacobianInverse (const LocalCoordinate &local) const
{
return jacobianInverseTransposed(local).transposed();
}
protected:
using Base::refElement;
......
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