Go to the documentation of this file.
   35#ifndef OPM_VOLUMES_HEADER 
   36#define OPM_VOLUMES_HEADER 
   43#include <dune/common/version.hh> 
   44#include <dune/common/math.hh> 
   45#include <dune/common/fvector.hh> 
   58    FieldVector<T, 3>  cross( const FieldVector<T, 3>& a,  const FieldVector<T, 3>&  b)  
   60        FieldVector<T, 3> res;  
   61        res[0] = a[1]* b[2] - a[2]* b[1];  
   62        res[1] = a[2]* b[0] - a[0]* b[2];  
   63        res[2] = a[0]* b[1] - a[1]* b[0];  
   72    template < class Vector>  
   75        return std::inner_product(a.begin(), a.end(),  b.begin(),  typename Vector::field_type());  
   81    template< typename T,  template < typename,  int>  class Point>  
   84        return a[0][0] * a[1][1] - a[1][0] * a[0][1];  
   91    template< typename T,  template < typename,  int>  class Point>  
   95            a[0][0] * (a[1][1] * a[2][2] - a[2][1] * a[1][2]) -  
   96            a[0][1] * (a[1][0] * a[2][2] - a[2][0] * a[1][2]) +  
   97            a[0][2] * (a[1][0] * a[2][1] - a[2][0] * a[1][1]);  
  103    template< typename T,  template < typename,  int>  class Point,  int Dim>  
  106        Point<T, Dim> tmp[Dim];  
  107        for ( int i = 0; i < Dim; ++i) {  
  108            tmp[i] = a[i+1] - a[i];  
  118    template < typename T,  template < typename,  int>  class Point>  
  119    inline T  area( const Point<T, 2>* c)  
  125    template < typename T,  template < typename,  int>  class Point>  
  126    inline T  area( const Point<T, 3>* c)  
  129        Point<T, 3> d0 = c[1] - c[0];  
  130        Point<T, 3> d1 = c[2] - c[0];  
  131        Point<T, 3> crossprod =  cross(d0,d1);  
  132        return 0.5 * crossprod.two_norm();  
  137    template < typename T,  template < typename,  int>  class Point>  
  144    template < typename T,  template < typename,  int>  class Point>  
  148        Point<T, 3> d0 = c[1] - c[0];  
  149        Point<T, 3> d1 = c[2] - c[0];  
  150        Point<T, 3> crossprod =  cross(d0, d1);  
  151        if ( inner(crossprod, normal) > 0) {  
  152            return 0.5 * crossprod.two_norm();  
  154            return -0.5 * crossprod.two_norm();  
The namespace Dune is the main namespace for all Dune code. Definition: common/CartesianIndexMapper.hpp:10  
T area(const Point< T, 2 > *c) Definition: Volumes.hpp:119  
FieldVector< T, 3 > cross(const FieldVector< T, 3 > &a, const FieldVector< T, 3 > &b) Definition: Volumes.hpp:58  
T volume(const Point< T, 3 > *c) Computes the volume of a 3D simplex (embedded i 3D space). Definition: Volumes.hpp:138  
T determinantOf(const Point< T, 2 > *a) Definition: Volumes.hpp:82  
T simplex_volume(const Point< T, Dim > *a) Definition: Volumes.hpp:104  
Vector::field_type inner(const Vector &a, const Vector &b) Definition: Volumes.hpp:73  
T signed_area(const Point< T, 3 > *c, const Point< T, 3 > &normal) Definition: Volumes.hpp:145  
Dune::FieldVector< double, 3 > Vector Definition: cpgrid/GridHelpers.hpp:385  
double b(const Dune::FieldVector< ct, dimworld > &, double) Definition: transportproblem2.hh:25  
  
  
 
    
     |