KASKADE 7 development version
Classes | Namespaces | Functions
fixdune.hh File Reference

This file contains various utility functions that augment the basic functionality of Dune. More...

#include <type_traits>
#include "dune/grid/config.h"
#include "dune/common/fmatrix.hh"
#include "dune/common/fvector.hh"
#include "dune/common/version.hh"
#include "dune/geometry/type.hh"
#include "dune/geometry/referenceelements.hh"
#include "dune/istl/matrix.hh"
#include "dune/istl/bvector.hh"
#include "dune/istl/bcrsmatrix.hh"
#include "utilities/deprecation.hh"
#include "utilities/scalar.hh"

Go to the source code of this file.

Classes

struct  Transpose< Matrix, transposed >
 Returns the matrix type or its transposed type, depending on the given transpose flag. More...
 
struct  Transpose< Dune::FieldMatrix< Scalar, n, m >, transposed >
 

Namespaces

namespace  Dune
 
namespace  Kaskade
 
 

Functions

template<class T , int n, class S , class enable = typename std::enable_if_t<std::is_arithmetic_v<S>>>
Dune::FieldVector< T, n > Dune::operator* (Dune::FieldVector< T, n > x, S s)
 Scalar-vector multiplication \( (s,x) \mapsto sx \). More...
 
template<class T , int n, class S , class enable = typename std::enable_if_t<std::is_arithmetic_v<S>>>
Dune::FieldVector< T, n > Dune::operator* (S s, Dune::FieldVector< T, n > x)
 Scalar-vector multiplication \( (x,s) \mapsto sx \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::operator+ (Dune::FieldVector< T, n > x, Dune::FieldVector< T, n > const &y)
 Vector addition \( (x,y) \mapsto x+y \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::operator+ (Dune::FieldVector< T, n > x, T const &y)
 Vector addition \( (x,y) \mapsto x+y \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::operator- (Dune::FieldVector< T, n > x)
 Vector negation \( x \mapsto - x \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::max (Dune::FieldVector< T, n > x, Dune::FieldVector< T, n > const &y)
 Componentwise maximum. More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::min (Dune::FieldVector< T, n > x, Dune::FieldVector< T, n > const &y)
 Componentwise minimum. More...
 
template<class T , int n, int m>
Dune::FieldVector< T, n *m > Dune::asVector (Dune::FieldMatrix< T, n, m > const &x)
 Converts a matrix of size n x m to a vector of size n*m by concatenating all the columns. More...
 
template<int n, int m, class T >
Dune::FieldMatrix< T, n, m > Dune::asMatrix (Dune::FieldVector< T, n *m > const &x)
 Converts a vector of size nm to a matrix of size n x m by filling columns successively. More...
 
template<class T , int n, int m>
Dune::FieldMatrix< T, n, m > Dune::outerProduct (Dune::FieldVector< T, n > const &x, Dune::FieldVector< T, m > const &y)
 outer vector product \( (x,y) \mapsto x y^T \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > Dune::normalize (Dune::FieldVector< T, n > x)
 returns the normalized vector More...
 
template<class T , int n, int m>
Dune::FieldVector< T, n+m > Dune::vertcat (Dune::FieldVector< T, n > const &x, Dune::FieldVector< T, m > const &y)
 Concatenation of vectors. More...
 
template<class T >
Dune::vectorProduct (Dune::FieldVector< T, 2 > const &x, Dune::FieldVector< T, 2 > const &y)
 vector product \( (x,y) \mapsto x \times y \in \mathbb{R} \). More...
 
template<class T >
Dune::FieldMatrix< T, 1, 2 > Dune::vectorProductMatrix (Dune::FieldVector< T, 2 > const &x)
 the matrix \( A(x) \in \mathbb{R}^{1\times 2}\) satisfying \( A(x) b = x\times b \) for all \( b \) More...
 
template<class T >
Dune::FieldVector< T, 3 > Dune::vectorProduct (Dune::FieldVector< T, 3 > const &x, Dune::FieldVector< T, 3 > const &y)
 vector product \( (x,y) \mapsto x \times y \). More...
 
template<class T >
Dune::FieldMatrix< T, 3, 3 > Dune::vectorProductMatrix (Dune::FieldVector< T, 3 > const &x)
 the skew-symmetric matrix \( A(x)\in \mathbb{R}^{3\times 3} \) satisfying \( A(x) b = x\times b \) for all \( b \) More...
 
template<class T , class S , int n, int m>
Dune::FieldVector< S, n > Dune::operator* (Dune::FieldMatrix< T, n, m > const &A, Dune::FieldVector< S, m > const &x)
 Matrix-vector multiplication \( (A,x) \mapsto Ax \). More...
 
template<class T , int n, int m>
Dune::FieldMatrix< T, n, m > Dune::operator- (Dune::FieldMatrix< T, n, m > const &A)
 Matrix negation \( A \mapsto - A \). More...
 
template<class T , int n, int m>
auto Dune::normalForm (Dune::FieldMatrix< T, n, m > const &A)
 
template<class Scalar , int n>
auto Kaskade::horzcat (Dune::FieldVector< Scalar, n > const &c1)
 
template<class Scalar , int n, class ... Rest>
auto Kaskade::horzcat (Dune::FieldVector< Scalar, n > const &a, Rest ... rest)
 Concatenates an arbitrary number of column vectors and matrices with coinciding number of rows to a larger matrix. More...
 
template<class Scalar , int n, int m>
auto Kaskade::horzcat (Dune::FieldMatrix< Scalar, n, m > const &a)
 
template<class Scalar , int n, int m, class ... Rest>
auto Kaskade::horzcat (Dune::FieldMatrix< Scalar, n, m > const &a, Rest ... rest)
 
template<class Scalar , int n>
Dune::FieldMatrix< Scalar, n, 3 > Kaskade::horzcat (Dune::FieldVector< Scalar, n > const &col1, Dune::FieldVector< Scalar, n > const &col2, Dune::FieldVector< Scalar, n > const &col3)
 
template<class Scalar , int n, int m>
Dune::FieldVector< Scalar, n > Kaskade::column (Dune::FieldMatrix< Scalar, n, m > const &A, int j)
 Extracts the j-th column of the given matrix A. More...
 
template<class Scalar , int n>
Dune::FieldVector< Scalar, n *n > Kaskade::vectorize (Dune::FieldMatrix< Scalar, n, n > const &A)
 matrix entries as vector, concatenating columns More...
 
template<int n, class Scalar >
Dune::FieldMatrix< Scalar, n, n > Kaskade::unvectorize (Dune::FieldVector< Scalar, n *n > const &v)
 forming quadratic matrix columns of vector segments More...
 
template<class T , int n, int m>
Kaskade::contraction (Dune::FieldMatrix< T, n, m > const &A, Dune::FieldMatrix< T, n, m > const &B)
 Matrix contraction (Frobenius product) \( (A,B) \mapsto A:B = \sum_{i,j} A_{ij} B_{ij} \). More...
 
template<class T , int n>
Dune::trace (Dune::FieldMatrix< T, n, n > const &A)
 Matrix trace \( A \mapsto \mathrm{tr}\,A = \sum_{i} A_{ii} \). More...
 
template<class T , int n>
determinant (Dune::FieldMatrix< T, n, n > const &A)
 Matrix determinant \( A \mapsto \mathrm{det}\,A \). More...
 
template<class T , int n>
Dune::FieldVector< T, n > diag (Dune::FieldMatrix< T, n, n > const &A)
 Matrix diagonal as a vector. More...
 
template<class T , int n>
Dune::FieldMatrix< T, n, n > unitMatrix ()
 Returns the identity matrix of size n \( \mapsto I \). More...
 
template<class T >
Dune::FieldMatrix< T, 3, 3 > skewMatrix (Dune::FieldVector< T, 3 > v)
 Returns the skew symmetric matrix related to the cross product with (3d) vector v. More...
 
template<class T , int n>
Dune::FieldVector< T, n > unitVector (int i)
 Returns the unit vector \( e_i \) with \( (e_i)_k = \delta_{ik} \). More...
 
template<class MatrixZ , class Matrix >
void MatMult (MatrixZ &z, Matrix const &x, Matrix const &y)
 Computes Z = X*Y. X and Y need to be compatible, i.e. X.M()==Y.N() has to hold. No aliasing may occur, i.e. &Z != &X and &Z != &Y must hold. Z is reshaped as needed. More...
 
template<class IndexSet , class Cell >
size_t subIndex (IndexSet const &is, Cell const &cell, int codim, int subentity)
 
template<class Cell >
int count (Cell const &cell, int codim)
 Computes the number of subentities of codimension c of a given entity. More...
 
template<class LocalCoordinate >
LocalCoordinate::field_type checkInside (Dune::GeometryType const &gt, LocalCoordinate const &x)
 Checks whether a point is inside or outside the reference element, and how much. More...
 
template<class Cell >
std::pair< Dune::FieldVector< double, Cell::dimension >, Dune::FieldVector< double, Cell::dimension > > boundingBox (Cell const &cell)
 Computes the bounding box of a cell. More...
 
template<class Scalar , int n, int m, class Allocator1 , class Allocator2 >
void transferBlockVector (Dune::BlockVector< Dune::FieldVector< Scalar, n >, Allocator1 > const &from, Dune::BlockVector< Dune::FieldVector< Scalar, m >, Allocator2 > &to)
 copy between two block vectors of different element lengths More...
 
template<class Scalar , int n, int m, class Allocator >
void bcrsPrint (Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Allocator > const &A, std::ostream &os=std::cout)
 Pretty prints the sparse matrix A to the given output stream. More...
 
template<class Scalar , int n, int m, class Allocator >
std::ostream & operator<< (std::ostream &os, Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Allocator > const &A)
 Pretty prints the sparse matrix A to the given output stream. More...
 

Detailed Description

This file contains various utility functions that augment the basic functionality of Dune.

Definition in file fixdune.hh.

Function Documentation

◆ count()

template<class Cell >
int count ( Cell const &  cell,
int  codim 
)

◆ subIndex()

template<class IndexSet , class Cell >
size_t subIndex ( IndexSet const &  is,
Cell const &  cell,
int  codim,
int  subentity 
)

Computes the index of the k-th subentity of codimension c in the index set is. This is available in Dune only with c given statically.

TODO: use IndexSet::IndexType instead of size_t for Dune 1.1

Definition at line 704 of file fixdune.hh.

◆ transferBlockVector()

template<class Scalar , int n, int m, class Allocator1 , class Allocator2 >
void transferBlockVector ( Dune::BlockVector< Dune::FieldVector< Scalar, n >, Allocator1 > const &  from,
Dune::BlockVector< Dune::FieldVector< Scalar, m >, Allocator2 > &  to 
)

copy between two block vectors of different element lengths

Definition at line 786 of file fixdune.hh.