pk2dlocalbasis.hh
Go to the documentation of this file.00001 #ifndef DUNE_PK2DLOCALBASIS_HH
00002 #define DUNE_PK2DLOCALBASIS_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 Pk2DLocalBasis
00024 {
00025 public:
00026
00028 enum {N = (k+1)*(k+2)/2};
00029
00033 enum {O = k};
00034
00035 typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,1,Dune::FieldVector<R,1>,
00036 Dune::FieldMatrix<R,1,2> > Traits;
00037
00039 Pk2DLocalBasis ()
00040 {
00041 for (unsigned int i=0; i<=k; i++)
00042 pos[i] = (1.0*i)/k;
00043 }
00044
00046 unsigned int size () const
00047 {
00048 return N;
00049 }
00050
00052 inline void evaluateFunction (const typename Traits::DomainType& x,
00053 std::vector<typename Traits::RangeType>& out) const
00054 {
00055 out.resize(N);
00056 int n=0;
00057 for (unsigned int j=0; j<=k; j++)
00058 for (unsigned int i=0; i<=k-j; i++)
00059 {
00060 out[n] = 1.0;
00061 for (unsigned int alpha=0; alpha<i; alpha++)
00062 out[n] *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
00063 for (unsigned int beta=0; beta<j; beta++)
00064 out[n] *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
00065 for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
00066 out[n] *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
00067 n++;
00068 }
00069 }
00070
00072 inline void
00073 evaluateJacobian (const typename Traits::DomainType& x,
00074 std::vector<typename Traits::JacobianType>& out) const
00075 {
00076 out.resize(N);
00077 int n=0;
00078 for (unsigned int j=0; j<=k; j++)
00079 for (unsigned int i=0; i<=k-j; i++)
00080 {
00081
00082 out[n][0][0] = 0.0;
00083 R factor=1.0;
00084 for (unsigned int beta=0; beta<j; beta++)
00085 factor *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
00086 for (unsigned int a=0; a<i; a++)
00087 {
00088 R product=factor;
00089 for (unsigned int alpha=0; alpha<i; alpha++)
00090 if (alpha==a)
00091 product *= 1.0/(pos[i]-pos[alpha]);
00092 else
00093 product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
00094 for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
00095 product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
00096 out[n][0][0] += product;
00097 }
00098 for (unsigned int c=i+j+1; c<=k; c++)
00099 {
00100 R product=factor;
00101 for (unsigned int alpha=0; alpha<i; alpha++)
00102 product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
00103 for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
00104 if (gamma==c)
00105 product *= -1.0/(pos[gamma]-pos[i]-pos[j]);
00106 else
00107 product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
00108 out[n][0][0] += product;
00109 }
00110
00111
00112 out[n][0][1] = 0.0;
00113 factor = 1.0;
00114 for (unsigned int alpha=0; alpha<i; alpha++)
00115 factor *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
00116 for (unsigned int b=0; b<j; b++)
00117 {
00118 R product=factor;
00119 for (unsigned int beta=0; beta<j; beta++)
00120 if (beta==b)
00121 product *= 1.0/(pos[j]-pos[beta]);
00122 else
00123 product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
00124 for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
00125 product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
00126 out[n][0][1] += product;
00127 }
00128 for (unsigned int c=i+j+1; c<=k; c++)
00129 {
00130 R product=factor;
00131 for (unsigned int beta=0; beta<j; beta++)
00132 product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
00133 for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
00134 if (gamma==c)
00135 product *= -1.0/(pos[gamma]-pos[i]-pos[j]);
00136 else
00137 product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
00138 out[n][0][1] += product;
00139 }
00140
00141 n++;
00142 }
00143
00144
00145
00146 }
00147
00149 unsigned int order () const
00150 {
00151 return k;
00152 }
00153
00154 private:
00155 R pos[k+1];
00156 };
00157
00158
00159
00160 template<class D, class R>
00161 class Pk2DLocalBasis<D,R,0>
00162 {
00163 public:
00164 typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,1,Dune::FieldVector<R,1>,
00165 Dune::FieldMatrix<R,1,2> > Traits;
00166
00168 enum {N = 1};
00169
00171 enum {O = 0};
00172
00173 unsigned int size () const
00174 {
00175 return 1;
00176 }
00177
00178 inline void evaluateFunction (const typename Traits::DomainType& in,
00179 std::vector<typename Traits::RangeType>& out) const
00180 {
00181 out.resize(1);
00182 out[0] = 1;
00183 }
00184
00185
00186 inline void
00187 evaluateJacobian (const typename Traits::DomainType& in,
00188 std::vector<typename Traits::JacobianType>& out) const
00189 {
00190 out.resize(1);
00191 out[0][0][0] = 0; out[0][0][1] = 0;
00192 }
00193
00194
00195 template<typename E, typename F, typename C>
00196 void interpolate (const E& e, const F& f, std::vector<C>& out) const
00197 {
00198 typename Traits::DomainType x;
00199 typename Traits::RangeType y;
00200 x[0] = 1.0/3.0; x[1] = 1.0/3.0; f.eval_local(e,x,y); out[0] = y;
00201 }
00202
00203 unsigned int order () const
00204 {
00205 return 0;
00206 }
00207 };
00208 }
00209 #endif