pk3dlocalbasis.hh

Go to the documentation of this file.
00001 #ifndef DUNE_PK3DLOCALBASIS_HH
00002 #define DUNE_PK3DLOCALBASIS_HH
00003 
00004 #include <dune/common/fmatrix.hh>
00005 
00006 #include <dune/localfunctions/common/localbasis.hh>
00007 
00008 namespace Dune
00009 {
00022 template<class D, class R, unsigned int k>
00023 class Pk3DLocalBasis
00024 {
00025 public:
00026     enum {N = (k+1)*(k+2)*(k+3)/6};
00027     enum {O = k};
00028 
00029     typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
00030                                Dune::FieldMatrix<R,1,3> > Traits;
00031 
00033     Pk3DLocalBasis () {}
00034 
00036     unsigned int size () const
00037     {
00038         return N;
00039     }
00040 
00042     inline void evaluateFunction (const typename Traits::DomainType& x,
00043                                   std::vector<typename Traits::RangeType>& out) const
00044     {
00045         out.resize(N);
00046         typename Traits::DomainType kx = x;
00047         kx *= k;
00048         unsigned int n = 0;
00049         unsigned int i[4];
00050         R factor[4];
00051         for (i[2] = 0; i[2] <= k; ++i[2])
00052         {
00053             factor[2] = 1.0;
00054             for (unsigned int j = 0; j < i[2]; ++j)
00055                 factor[2] *= (kx[2]-j) / (i[2]-j);
00056             for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
00057             {
00058                 factor[1] = 1.0;
00059                 for (unsigned int j = 0; j < i[1]; ++j)
00060                     factor[1] *= (kx[1]-j) / (i[1]-j);
00061                 for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
00062                 {
00063                     factor[0] = 1.0;
00064                     for (unsigned int j = 0; j < i[0]; ++j)
00065                         factor[0] *= (kx[0]-j) / (i[0]-j);
00066                     i[3] = k - i[0] - i[1] - i[2];
00067                     D kx3 = k - kx[0] - kx[1] - kx[2];
00068                     factor[3] = 1.0;
00069                     for (unsigned int j = 0; j < i[3]; ++j)
00070                         factor[3] *= (kx3-j) / (i[3]-j);
00071                     out[n++] = factor[0] * factor[1] * factor[2] * factor[3];
00072                 }
00073             }
00074         }
00075     }
00076 
00078     inline void
00079     evaluateJacobian (const typename Traits::DomainType& x,         // position
00080                       std::vector<typename Traits::JacobianType>& out) const      // return value
00081     {
00082         out.resize(N);
00083         typename Traits::DomainType kx = x;
00084         kx *= k;
00085         unsigned int n = 0;
00086         unsigned int i[4];
00087         R factor[4];
00088         for (i[2] = 0; i[2] <= k; ++i[2])
00089         {
00090             factor[2] = 1.0;
00091             for (unsigned int j = 0; j < i[2]; ++j)
00092                 factor[2] *= (kx[2]-j) / (i[2]-j);
00093             for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
00094             {
00095                 factor[1] = 1.0;
00096                 for (unsigned int j = 0; j < i[1]; ++j)
00097                     factor[1] *= (kx[1]-j) / (i[1]-j);
00098                 for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
00099                 {
00100                     factor[0] = 1.0;
00101                     for (unsigned int j = 0; j < i[0]; ++j)
00102                         factor[0] *= (kx[0]-j) / (i[0]-j);
00103                     i[3] = k - i[0] - i[1] - i[2];
00104                     D kx3 = k - kx[0] - kx[1] - kx[2];
00105                     R sum3 = 0.0;
00106                     factor[3] = 1.0;
00107                     for (unsigned int j = 0; j < i[3]; ++j)
00108                         factor[3] /= i[3] - j;
00109                     R prod_all = factor[0] * factor[1] * factor[2] * factor[3];
00110                     for (unsigned int j = 0; j < i[3]; ++j)
00111                     {
00112                         R prod = prod_all;
00113                         for (unsigned int l = 0; l < i[3]; ++l)
00114                             if (j == l)
00115                                 prod *= -R(k);
00116                             else
00117                                 prod *= kx3 - l;
00118                         sum3 += prod;
00119                     }
00120                     for (unsigned int j = 0; j < i[3]; ++j)
00121                         factor[3] *= kx3 - j;
00122                     for (unsigned int m = 0; m < 3; ++m)
00123                     {
00124                         out[n][0][m] = sum3;
00125                         for (unsigned int j = 0; j < i[m]; ++j)
00126                         {
00127                             R prod = factor[3];
00128                             for (unsigned int p = 0; p < 3; ++p)
00129                             {
00130                                 if (m == p)
00131                                     for (unsigned int l = 0; l < i[p]; ++l)
00132                                         if (j == l)
00133                                             prod *= R(k) / (i[p]-l);
00134                                         else
00135                                             prod *= (kx[p]-l) / (i[p]-l);
00136                                 else
00137                                     prod *= factor[p];
00138                             }
00139                             out[n][0][m] += prod;
00140                         }
00141                     }
00142                     n++;
00143                 }
00144             }
00145         }
00146     }
00147 
00149     unsigned int order () const
00150     {
00151         return k;
00152     }
00153 };
00154 
00155 
00156 //Specialization for k=0
00157 template<class D, class R>
00158 class Pk3DLocalBasis<D,R,0>
00159 {
00160 public:
00161     typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
00162                                Dune::FieldMatrix<R,1,3> > Traits;
00163 
00165     enum {N = 1};
00166 
00168     enum {O = 0};
00169 
00170     unsigned int size () const
00171     {
00172         return 1;
00173     }
00174 
00175     inline void evaluateFunction (const typename Traits::DomainType& in,
00176                                   std::vector<typename Traits::RangeType>& out) const
00177     {
00178         out.resize(1);
00179         out[0] = 1;
00180     }
00181 
00182     // evaluate derivative of a single component
00183     inline void
00184     evaluateJacobian (const typename Traits::DomainType& in,         // position
00185                       std::vector<typename Traits::JacobianType>& out) const      // return value
00186     {
00187         out.resize(1);
00188         out[0][0][0] = 0;
00189         out[0][0][1] = 0;
00190         out[0][0][2] = 0;
00191     }
00192 
00193     // local interpolation of a function
00194     template<typename E, typename F, typename C>
00195     void interpolate (const E& e, const F& f, std::vector<C>& out) const
00196     {
00197         typename Traits::DomainType x;
00198         typename Traits::RangeType y;
00199         x[0] = 1.0/4.0;
00200         x[1] = 1.0/4.0;
00201         x[2] = 1.0/4.0;
00202         f.eval_local(e,x,y);
00203         out[0] = y;
00204     }
00205 
00206     unsigned int order () const
00207     {
00208         return 0;
00209     }
00210 };
00211 }
00212 #endif
Generated on Sat Apr 24 11:15:35 2010 for dune-localfunctions by  doxygen 1.6.3