field.hh

Go to the documentation of this file.
00001 #ifndef DUNE_FIELD_HH
00002 #define DUNE_FIELD_HH
00003 
00004 #if HAVE_ALGLIB
00005 #include <alglib/amp.h>
00006 #endif
00007 
00008 #include <dune/common/gmpfield.hh>
00009 #include <dune/common/fvector.hh>
00010 
00011 namespace Dune
00012 {
00013 
00014   // Unity
00015   // -----
00016 
00027   template< class Field >
00028   struct Unity
00029   {
00030     operator Field () const
00031     {
00032       return Field( 1 );
00033     }
00034   };
00035 
00036   template< class Field >
00037   Field operator+ ( const Unity< Field > &u, const Field &f )
00038   {
00039     return (Field)u + f;
00040   }
00041 
00042   template< class Field >
00043   Field operator- ( const Unity< Field > &u, const Field &f )
00044   {
00045     return (Field)u - f;
00046   }
00047 
00048   template< class Field >
00049   Field operator* ( const Unity< Field > &u, const Field &f )
00050   {
00051     return f;
00052   }
00053 
00054   template< class Field >
00055   Field operator/ ( const Unity< Field > &u, const Field &f )
00056   {
00057     return (Field)u / f;
00058   }
00059 
00060 
00061 
00062   // Zero
00063   // ----
00064 
00076   template< class Field >
00077   struct Zero
00078   {
00079     operator Field () const
00080     {
00081       return Field( 0 );
00082     }
00083     static const Field epsilon()
00084     {
00085       return Field(1e-12);
00086     }
00087   };
00088 #if HAVE_ALGLIB
00089   template< unsigned int precision >
00090   struct Zero< amp::ampf< precision > >
00091   {
00092     typedef amp::ampf< precision > Field;
00093     operator Field () const
00094     {
00095       return Field( 0 );
00096     }
00097     static const Field epsilon()
00098     {
00099       return Field(1e-20);
00100     }
00101   };
00102 #endif
00103 #if HAVE_GMP
00104   template< unsigned int precision >
00105   struct Zero< GMPField< precision > >
00106   {
00107     typedef GMPField< precision > Field;
00108     operator Field () const
00109     {
00110       return Field( 0 );
00111     }
00112     static const Field epsilon()
00113     {
00114       return Field(1e-20);
00115     }
00116   };
00117 #endif
00118 
00119   template< class Field >
00120   inline bool operator == ( const Zero< Field > &, const Field &f )
00121   {
00122     return ( f < Zero<Field>::epsilon() && f > -Zero<Field>::epsilon() );
00123   }
00124 
00125   template< class Field >
00126   inline bool operator == ( const Field &f, const Zero< Field > &z)
00127   {
00128     return ( z == f );
00129   }
00130 
00131   template< class Field >
00132   inline bool operator< ( const Zero< Field > &, const Field &f )
00133   {
00134     return f > Zero<Field>::epsilon();
00135   }
00136 
00137   template< class Field >
00138   inline bool operator< ( const Field &f, const Zero< Field > & )
00139   {
00140     return f < -Zero<Field>::epsilon();
00141   }
00142 
00143   template< class Field >
00144   inline bool operator> ( const Zero< Field > &z, const Field &f )
00145   {
00146     return f < z;
00147   }
00148 
00149   template< class Field >
00150   inline bool operator> ( const Field &f, const Zero< Field > &z )
00151   {
00152     return z < f;
00153   }
00154 
00155 
00156   // field_cast
00157   // ----------
00158 
00171   template< class F2, class F1 >
00172   inline void field_cast ( const F1 &f1, F2 &f2 )
00173   {
00174     f2 = f1;
00175   }
00176 
00177 #if HAVE_ALGLIB
00178   template< unsigned int precision >
00179   inline void field_cast ( const amp::ampf< precision > &f1, double &f2 )
00180   {
00181     f2 = f1.toDouble();
00182   }
00183 
00184   template< unsigned int precision >
00185   inline void field_cast ( const amp::ampf< precision > &f1, long double &f2 )
00186   {
00187     f2 = f1.toDouble();
00188   }
00189 #endif
00190 
00191 #if HAVE_GMP
00192   template< unsigned int precision >
00193   inline void field_cast ( const Dune::GMPField< precision > &f1, double &f2 )
00194   {
00195     f2 = f1.get_d();
00196   }
00197 
00198   template< unsigned int precision >
00199   inline void field_cast ( const Dune::GMPField< precision > &f1, long double &f2 )
00200   {
00201     f2 = f1.get_d();
00202   }
00203 #endif
00204 
00205   template< class F2, class F1, int dim >
00206   inline void field_cast ( const Dune::FieldVector< F1, dim > &f1, Dune::FieldVector< F2, dim > &f2 )
00207   {
00208     for( int d = 0; d < dim; ++d )
00209       field_cast( f1[ d ], f2[ d ] );
00210   }
00211   template< class F2, class F1 >
00212   inline void field_cast ( const Dune::FieldVector< F1, 1 > &f1, F2 &f2 )
00213   {
00214     field_cast( f1[ 0 ], f2 );
00215   }
00216   template< class F2, class F1 >
00217   inline void field_cast ( const F1 &f1, Dune::FieldVector< F2, 1 > &f2 )
00218   {
00219     field_cast( f1, f2[ 0 ] );
00220   }
00221 
00222   template< class F2, class F1, int rdim, int cdim >
00223   inline void field_cast ( const Dune::FieldMatrix< F1, rdim, cdim > &f1, Dune::FieldMatrix< F2, rdim, cdim > &f2 )
00224   {
00225     for( int r = 0; r < rdim; ++r )
00226       field_cast( f1[ r ], f2[ r ] );
00227   }
00228   template< class F2, class F1 >
00229   inline void field_cast ( const Dune::FieldMatrix<F1,1,1> &f1, Dune::FieldMatrix< F2, 1,1 > &f2 )
00230   {
00231     field_cast( f1[ 0 ][ 0 ], f2[ 0 ][ 0 ] );
00232   }
00233   template< class F2, class F1 >
00234   inline void field_cast ( const Dune::FieldMatrix< F1, 1,1 > &f1, F2 &f2 )
00235   {
00236     field_cast( f1[ 0 ][ 0 ], f2 );
00237   }
00238   template< class F2, class F1 >
00239   inline void field_cast ( const F1 &f1, Dune::FieldMatrix< F2, 1,1 > &f2 )
00240   {
00241     field_cast( f1, f2[ 0 ][ 0 ] );
00242   }
00243   template< class F2, class F1 >
00244   inline void field_cast ( const Dune::FieldVector<F1,1> &f1, Dune::FieldMatrix< F2, 1,1 > &f2 )
00245   {
00246     field_cast( f1[ 0 ], f2[ 0 ][ 0 ] );
00247   }
00248   template< class F2, class F1 >
00249   inline void field_cast ( const Dune::FieldMatrix<F1,1,1> &f1, Dune::FieldVector< F2, 1 > &f2 )
00250   {
00251     field_cast( f1[ 0 ][ 0 ], f2[ 0 ] );
00252   }
00253 
00254   template< class F2, class F1 >
00255   inline void field_cast ( const Dune::FieldVector< F1, 1 > &f1, Dune::FieldVector<F2, 1> &f2 )
00256   {
00257     field_cast( f1[ 0 ], f2[ 0 ] );
00258   }
00259 
00260   template< class F2,class V > 
00261   struct FieldCast
00262   {
00263     typedef F2 type;
00264   };
00265   template< class F2,class F1,int dim > 
00266   struct FieldCast< F2, Dune::FieldVector<F1,dim> >
00267   {
00268     typedef Dune::FieldVector<F2,dim> type;
00269   };
00270   template< class F2,class F1,int dim1, int dim2> 
00271   struct FieldCast< F2, Dune::FieldMatrix<F1,dim1,dim2> >
00272   {
00273     typedef Dune::FieldMatrix<F2,dim1,dim2> type;
00274   };
00275   template< class F2,class V > 
00276   inline typename FieldCast<F2,V>::type field_cast ( const V &f1 )
00277   {
00278     typename FieldCast<F2,V>::type f2;
00279     field_cast( f1, f2 );
00280     return f2;
00281   }
00282 
00283 
00284   // Precision
00285   // this is not a perfect solution to obtain the
00286   // precision of a field - definition is not clear
00287   // to be removed
00288   // ---------
00289 
00290   template <class Field>
00291   struct Precision;
00292 
00293   template<>
00294   struct Precision< double >
00295   {
00296     static const unsigned int value = 64;
00297   };
00298 
00299   template<>
00300   struct Precision< long double >
00301   {
00302     static const unsigned int value = 80;
00303   };
00304 
00305   template<>
00306   struct Precision< float >
00307   {
00308     static const unsigned int value = 32;
00309   };
00310 
00311 #if HAVE_ALGLIB
00312   template< unsigned int precision >
00313   struct Precision< amp::ampf< precision > >
00314   {
00315     static const unsigned int value = precision;
00316   };
00317 #endif
00318 #if HAVE_GMP
00319   template< unsigned int precision >
00320   struct Precision< GMPField< precision > >
00321   {
00322     static const unsigned int value = precision;
00323   };
00324 #endif
00325 
00326   // ComputeField
00327   // ------------
00328 
00329   template <class Field,unsigned int sum>
00330   struct ComputeField 
00331   {
00332     typedef Field Type;
00333   };
00334 
00335 #if HAVE_ALGLIB
00336   template< unsigned int precision, unsigned int sum >
00337   struct ComputeField< amp::ampf< precision >, sum >
00338   {
00339     typedef amp::ampf<precision+sum> Type;
00340   };
00341 #endif
00342 #if HAVE_GMP
00343   template< unsigned int precision, unsigned int sum >
00344   struct ComputeField< GMPField< precision >, sum >
00345   {
00346     typedef GMPField<precision+sum> Type;
00347   };
00348 #endif
00349 }
00350 
00351 // to be moved to different location...
00352 namespace std
00353 {
00354 
00355 #if HAVE_ALGLIB
00356   template< unsigned int precision >
00357   inline ostream &
00358   operator<< ( ostream &out, 
00359                const amp::ampf< precision > &value )
00360   {
00361     return out << value.toDec();
00362   }
00363 
00364   template< unsigned int precision >
00365   inline amp::ampf< precision > sqrt ( const amp::ampf< precision > &a )
00366   {
00367     return amp::sqrt( a );
00368   }
00369 
00370   template< unsigned int precision >
00371   inline amp::ampf< precision > abs ( const amp::ampf< precision > &a )
00372   {
00373     return amp::abs( a );
00374   }
00375 #endif // #if HAVE_ALGLIB
00376 
00377 }
00378 
00379 #endif // #ifndef DUNE_FIELD_HH
Generated on Sat Apr 24 11:15:33 2010 for dune-localfunctions by  doxygen 1.6.3