Dune Core Modules (2.8.0)

field.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_LOCALFUNCTIONS_UTILITY_FIELD_HH
4#define DUNE_LOCALFUNCTIONS_UTILITY_FIELD_HH
5
9
10namespace Dune
11{
12
13 // Unity
14 // -----
15
26 template< class Field >
27 struct Unity
28 {
29 operator Field () const
30 {
31 return Field( 1 );
32 }
33 };
34
35 template< class Field >
36 Field operator+ ( const Unity< Field > &u, const Field &f )
37 {
38 return (Field)u + f;
39 }
40
41 template< class Field >
42 Field operator- ( const Unity< Field > &u, const Field &f )
43 {
44 return (Field)u - f;
45 }
46
47 template< class Field >
48 Field operator* ( const Unity< Field > &u, const Field &f )
49 {
50 return f;
51 }
52
53 template< class Field >
54 Field operator/ ( const Unity< Field > &u, const Field &f )
55 {
56 return (Field)u / f;
57 }
58
59
60
61 // Zero
62 // ----
63
75 template< class Field >
76 struct Zero
77 {
78 operator Field () const
79 {
80 return Field( 0 );
81 }
82 static const Field epsilon()
83 {
84 return Field(1e-12);
85 }
86 };
87
88#if HAVE_GMP
89 template< unsigned int precision >
90 struct Zero< GMPField< precision > >
91 {
92 typedef GMPField< precision > Field;
93 operator Field () const
94 {
95 return Field( 0 );
96 }
97 static const Field epsilon()
98 {
99 return Field(1e-20);
100 }
101 };
102#endif
103
104 template< class Field >
105 inline bool operator == ( const Zero< Field > &, const Field &f )
106 {
107 return ( f < Zero<Field>::epsilon() && f > -Zero<Field>::epsilon() );
108 }
109
110 template< class Field >
111 inline bool operator == ( const Field &f, const Zero< Field > &z)
112 {
113 return ( z == f );
114 }
115
116 template< class Field >
117 inline bool operator< ( const Zero< Field > &, const Field &f )
118 {
119 return f > Zero<Field>::epsilon();
120 }
121
122 template< class Field >
123 inline bool operator< ( const Field &f, const Zero< Field > & )
124 {
125 return f < -Zero<Field>::epsilon();
126 }
127
128 template< class Field >
129 inline bool operator> ( const Zero< Field > &z, const Field &f )
130 {
131 return f < z;
132 }
133
134 template< class Field >
135 inline bool operator> ( const Field &f, const Zero< Field > &z )
136 {
137 return z < f;
138 }
139
140
141 // field_cast
142 // ----------
143
156 template< class F2, class F1 >
157 inline void field_cast ( const F1 &f1, F2 &f2 )
158 {
159 f2 = f1;
160 }
161
162#if HAVE_GMP
163 template< unsigned int precision >
164 inline void field_cast ( const Dune::GMPField< precision > &f1, double &f2 )
165 {
166 f2 = f1.get_d();
167 }
168
169 template< unsigned int precision >
170 inline void field_cast ( const Dune::GMPField< precision > &f1, long double &f2 )
171 {
172 f2 = f1.get_d();
173 }
174#endif
175
176 template< class F2, class F1, int dim >
178 {
179 for( int d = 0; d < dim; ++d )
180 field_cast( f1[ d ], f2[ d ] );
181 }
182 template< class F2, class F1 >
183 inline void field_cast ( const Dune::FieldVector< F1, 1 > &f1, F2 &f2 )
184 {
185 field_cast( f1[ 0 ], f2 );
186 }
187 template< class F2, class F1 >
188 inline void field_cast ( const F1 &f1, Dune::FieldVector< F2, 1 > &f2 )
189 {
190 field_cast( f1, f2[ 0 ] );
191 }
192
193 template< class F2, class F1, int rdim, int cdim >
195 {
196 for( int r = 0; r < rdim; ++r )
197 field_cast( f1[ r ], f2[ r ] );
198 }
199 template< class F2, class F1 >
201 {
202 field_cast( f1[ 0 ][ 0 ], f2[ 0 ][ 0 ] );
203 }
204 template< class F2, class F1 >
205 inline void field_cast ( const Dune::FieldMatrix< F1, 1,1 > &f1, F2 &f2 )
206 {
207 field_cast( f1[ 0 ][ 0 ], f2 );
208 }
209 template< class F2, class F1 >
210 inline void field_cast ( const F1 &f1, Dune::FieldMatrix< F2, 1,1 > &f2 )
211 {
212 field_cast( f1, f2[ 0 ][ 0 ] );
213 }
214 template< class F2, class F1 >
216 {
217 field_cast( f1[ 0 ], f2[ 0 ][ 0 ] );
218 }
219 template< class F2, class F1 >
221 {
222 field_cast( f1[ 0 ][ 0 ], f2[ 0 ] );
223 }
224
225 template< class F2, class F1 >
227 {
228 field_cast( f1[ 0 ], f2[ 0 ] );
229 }
230
231 template< class F2,class V >
232 struct FieldCast
233 {
234 typedef F2 type;
235 };
236 template< class F2,class F1,int dim >
237 struct FieldCast< F2, Dune::FieldVector<F1,dim> >
238 {
239 typedef Dune::FieldVector<F2,dim> type;
240 };
241 template< class F2,class F1,int dim1, int dim2>
242 struct FieldCast< F2, Dune::FieldMatrix<F1,dim1,dim2> >
243 {
245 };
246 template< class F2,class V >
247 inline typename FieldCast<F2,V>::type field_cast ( const V &f1 )
248 {
249 typename FieldCast<F2,V>::type f2;
250 field_cast( f1, f2 );
251 return f2;
252 }
253
254
255 // Precision
256 // this is not a perfect solution to obtain the
257 // precision of a field - definition is not clear
258 // to be removed
259 // ---------
260
261 template <class Field>
262 struct Precision;
263
264 template<>
265 struct Precision< double >
266 {
267 static const unsigned int value = 64;
268 };
269
270 template<>
271 struct Precision< long double >
272 {
273 static const unsigned int value = 80;
274 };
275
276 template<>
277 struct Precision< float >
278 {
279 static const unsigned int value = 32;
280 };
281
282#if HAVE_GMP
283 template< unsigned int precision >
284 struct Precision< GMPField< precision > >
285 {
286 static const unsigned int value = precision;
287 };
288#endif
289
290 // ComputeField
291 // ------------
292
293 template <class Field,unsigned int sum>
294 struct ComputeField
295 {
296 typedef Field Type;
297 };
298
299#if HAVE_GMP
300 template< unsigned int precision, unsigned int sum >
301 struct ComputeField< GMPField< precision >, sum >
302 {
303 typedef GMPField<precision+sum> Type;
304 };
305#endif
306} // namespace Dune
307
308#endif // #ifndef DUNE_LOCALFUNCTIONS_UTILITY_FIELD_HH
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Number class for high precision floating point number using the GMP library mpf_class implementation.
Definition: gmpfield.hh:31
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:269
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:635
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:681
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:235
Dune namespace.
Definition: alignedallocator.hh:11
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)