3#ifndef DUNE_LOCALFUNCTIONS_UTILITY_FIELD_HH
4#define DUNE_LOCALFUNCTIONS_UTILITY_FIELD_HH
26 template<
class Field >
29 operator Field ()
const
35 template<
class Field >
41 template<
class Field >
42 Field operator- (
const Unity< Field > &u,
const Field &f )
47 template<
class Field >
48 Field operator* (
const Unity< Field > &u,
const Field &f )
53 template<
class Field >
54 Field operator/ (
const Unity< Field > &u,
const Field &f )
75 template<
class Field >
78 operator Field ()
const
82 static const Field epsilon()
89 template<
unsigned int precision >
93 operator Field ()
const
97 static const Field epsilon()
104 template<
class Field >
105 inline bool operator == (
const Zero< Field > &,
const Field &f )
107 return ( f < Zero<Field>::epsilon() && f > -Zero<Field>::epsilon() );
110 template<
class Field >
111 inline bool operator == (
const Field &f,
const Zero< Field > &z)
116 template<
class Field >
117 inline bool operator< (
const Zero< Field > &,
const Field &f )
119 return f > Zero<Field>::epsilon();
122 template<
class Field >
123 inline bool operator< (
const Field &f,
const Zero< Field > & )
125 return f < -Zero<Field>::epsilon();
128 template<
class Field >
129 inline bool operator> (
const Zero< Field > &z,
const Field &f )
134 template<
class Field >
135 inline bool operator> (
const Field &f,
const Zero< Field > &z )
156 template<
class F2,
class F1 >
163 template<
unsigned int precision >
169 template<
unsigned int precision >
176 template<
class F2,
class F1,
int dim >
179 for(
int d = 0; d < dim; ++d )
182 template<
class F2,
class F1 >
187 template<
class F2,
class F1 >
193 template<
class F2,
class F1,
int rdim,
int cdim >
196 for(
int r = 0; r < rdim; ++r )
199 template<
class F2,
class F1 >
204 template<
class F2,
class F1 >
209 template<
class F2,
class F1 >
214 template<
class F2,
class F1 >
219 template<
class F2,
class F1 >
225 template<
class F2,
class F1 >
231 template<
class F2,
class V >
236 template<
class F2,
class F1,
int dim >
237 struct FieldCast< F2,
Dune::FieldVector<F1,dim> >
241 template<
class F2,
class F1,
int dim1,
int dim2>
242 struct FieldCast< F2,
Dune::FieldMatrix<F1,dim1,dim2> >
246 template<
class F2,
class V >
247 inline typename FieldCast<F2,V>::type
field_cast (
const V &f1 )
249 typename FieldCast<F2,V>::type f2;
261 template <
class Field>
265 struct Precision< double >
267 static const unsigned int value = 64;
271 struct Precision< long double >
273 static const unsigned int value = 80;
277 struct Precision< float >
279 static const unsigned int value = 32;
283 template<
unsigned int precision >
284 struct Precision< GMPField< precision > >
286 static const unsigned int value = precision;
293 template <
class Field,
unsigned int sum>
300 template<
unsigned int precision,
unsigned int sum >
301 struct ComputeField< GMPField< precision >, sum >
303 typedef GMPField<precision+sum> Type;
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Number class for high precision floating point number using the GMP library mpf_class implementation.
Definition: gmpfield.hh:29
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
Wrapper for the GNU multiprecision (GMP) library.
EnableIfInterOperable< T1, T2, bool >::type operator<(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:629
EnableIfInterOperable< T1, T2, bool >::type operator>(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:675
EnableIfInterOperable< T1, T2, bool >::type operator==(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for equality.
Definition: iteratorfacades.hh:233
Dune namespace.
Definition: alignedallocator.hh:10
void field_cast(const F1 &f1, F2 &f2)
a helper class to cast from one field to another
Definition: field.hh:157
A class representing the unit of a given Field.
Definition: field.hh:28
A class representing the zero of a given Field.
Definition: field.hh:77