00001
00002 #ifndef DUNE_MONOMLOCALBASIS_HH
00003 #define DUNE_MONOMLOCALBASIS_HH
00004
00005 #include <cassert>
00006
00007 #include <dune/common/fmatrix.hh>
00008
00009 #include <dune/grid/common/referenceelements.hh>
00010
00011 #include"../common/localbasis.hh"
00012
00013 namespace Dune
00014 {
00015 namespace MonomImp {
00019 template<int d, int k>
00020 struct Size {
00021 enum { val = Size<d,k-1>::val+Size<d-1,k>::val };
00022 };
00023 template<int d>
00024 struct Size<d, 0> {
00025 enum { val = 1 };
00026 };
00027 template<int k>
00028 struct Size<0, k> {
00029 enum { val = 1 };
00030 };
00031
00033 template <typename Traits>
00034 class EvalAccess {
00035 std::vector<typename Traits::RangeType> &out;
00036 #ifndef NDEBUG
00037 unsigned int first_unused_index;
00038 #endif
00039
00040 public:
00041 EvalAccess(std::vector<typename Traits::RangeType> &out_)
00042 : out(out_)
00043 #ifndef NDEBUG
00044 , first_unused_index(0)
00045 #endif
00046 { }
00047 #ifndef NDEBUG
00048 ~EvalAccess() {
00049 assert(first_unused_index == out.size());
00050 }
00051 #endif
00052 typename Traits::RangeFieldType &operator[](unsigned int index)
00053 {
00054 assert(index < out.size());
00055 #ifndef NDEBUG
00056 if(first_unused_index <= index)
00057 first_unused_index = index+1;
00058 #endif
00059 return out[index][0];
00060 }
00061 };
00062
00064 template <typename Traits>
00065 class JacobianAccess {
00066 std::vector<typename Traits::JacobianType> &out;
00067 unsigned int row;
00068 #ifndef NDEBUG
00069 unsigned int first_unused_index;
00070 #endif
00071
00072 public:
00073 JacobianAccess(std::vector<typename Traits::JacobianType> &out_,
00074 unsigned int row_)
00075 : out(out_), row(row_)
00076 #ifndef NDEBUG
00077 , first_unused_index(0)
00078 #endif
00079 { }
00080 #ifndef NDEBUG
00081 ~JacobianAccess() {
00082 assert(first_unused_index == out.size());
00083 }
00084 #endif
00085 typename Traits::RangeFieldType &operator[](unsigned int index)
00086 {
00087 assert(index < out.size());
00088 #ifndef NDEBUG
00089 if(first_unused_index <= index)
00090 first_unused_index = index+1;
00091 #endif
00092 return out[index][0][row];
00093 }
00094 };
00095
00108 template <typename Traits, int c>
00109 struct Evaluate
00110 {
00111 enum {
00113 d = Traits::dimDomain - c
00114 };
00121 template <typename Access>
00122 static void eval (
00123 const typename Traits::DomainType &in,
00126 const array<int, Traits::dimDomain> &derivatives,
00129 typename Traits::RangeFieldType prod,
00131 int bound,
00133 int& index,
00135 Access &access)
00136 {
00137
00138 for (int e = bound; e >= 0; --e)
00139 {
00140
00141
00142 int newbound = bound - e;
00143 if(e < derivatives[d])
00144 Evaluate<Traits,c-1>::
00145 eval(in, derivatives, 0, newbound, index, access);
00146 else {
00147 int coeff = 1;
00148 for(int i = e - derivatives[d] + 1; i <= e; ++i)
00149 coeff *= i;
00150
00151 Evaluate<Traits,c-1>::
00152 eval(
00153 in, derivatives,
00154
00155
00156 prod*std::pow(in[d], typename Traits::DomainFieldType(e-derivatives[d]))*coeff,
00157
00158
00159 newbound,
00160
00161
00162 index, access);
00163 }
00164 }
00165 }
00166 };
00167
00172 template <typename Traits>
00173 struct Evaluate<Traits, 1>
00174 {
00175 enum { d = Traits::dimDomain-1 };
00177 template <typename Access>
00178 static void eval (const typename Traits::DomainType &in,
00179 const array<int, Traits::dimDomain> &derivatives,
00180 typename Traits::RangeFieldType prod,
00181 int bound, int& index, Access &access)
00182 {
00183 if(bound < derivatives[d])
00184 prod = 0;
00185 else {
00186 int coeff = 1;
00187 for(int i = bound - derivatives[d] + 1; i <= bound; ++i)
00188 coeff *= i;
00189 prod *= std::pow(in[d], typename Traits::DomainFieldType(bound-derivatives[d]))*coeff;
00190 }
00191 access[index] = prod;
00192 ++index;
00193 }
00194 };
00195
00196 }
00197
00211 template<class D, class R, unsigned int d, unsigned int p>
00212 class MonomLocalBasis
00213 {
00214 enum { static_size = MonomImp::Size<d,p>::val };
00215
00216 public:
00218 typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,Dune::FieldVector<R,1>,
00219 Dune::FieldMatrix<R,1,d>,p> Traits;
00220
00222 unsigned int size () const
00223 {
00224 return static_size;
00225 }
00226
00228 inline void evaluateFunction (const typename Traits::DomainType& in,
00229 std::vector<typename Traits::RangeType>& out) const
00230 {
00231 evaluate<0>(array<int, 0>(), in, out);
00232 }
00233
00235 template<unsigned int k>
00236 inline void evaluate (const array<int,k>& directions,
00237 const typename Traits::DomainType& in,
00238 std::vector<typename Traits::RangeType>& out) const
00239 {
00240 out.resize(size());
00241 int index = 0;
00242 array<int, d> derivatives;
00243 for(unsigned int i = 0; i < d; ++i) derivatives[i] = 0;
00244 for(unsigned int i = 0; i < k; ++i) ++derivatives[directions[i]];
00245 MonomImp::EvalAccess<Traits> access(out);
00246 for(unsigned int lp = 0; lp <= p; ++lp)
00247 MonomImp::Evaluate<Traits, d>::eval(in, derivatives, 1, lp, index,
00248 access);
00249 }
00250
00252 inline void
00253 evaluateJacobian (const typename Traits::DomainType& in,
00254 std::vector<typename Traits::JacobianType>& out) const
00255 {
00256 out.resize(size());
00257 array<int, d> derivatives;
00258 for(unsigned int i = 0; i < d; ++i)
00259 derivatives[i] = 0;
00260 for(unsigned int i = 0; i < d; ++i)
00261 {
00262 derivatives[i] = 1;
00263 int index = 0;
00264 MonomImp::JacobianAccess<Traits> access(out, i);
00265 for(unsigned int lp = 0; lp <= p; ++lp)
00266 MonomImp::Evaluate<Traits, d>::eval(in, derivatives, 1, lp, index, access);
00267 derivatives[i] = 0;
00268 }
00269 }
00270
00272 unsigned int order () const
00273 {
00274 return p;
00275 }
00276 };
00277
00278 }
00279
00280 #endif // DUNE_MONOMLOCALBASIS_HH