edges02dlocalbasis.hh
Go to the documentation of this file.00001
00002 #ifndef DUNE_EDGES02DLOCALBASIS_HH
00003 #define DUNE_EDGES02DLOCALBASIS_HH
00004
00005 #include <cmath>
00006
00007 #include <dune/common/fmatrix.hh>
00008
00009 #include <dune/localfunctions/common/localbasis.hh>
00010
00011 namespace Dune
00012 {
00029 template<class D, class R>
00030 class EdgeS02DLocalBasis
00031 {
00032 public:
00034 typedef LocalBasisTraits<
00035 D,2,Dune::FieldVector<D,2>,
00036 R,2,Dune::FieldVector<R,2>,
00037 Dune::FieldMatrix<R,2,2>
00038 > Traits;
00039
00041 EdgeS02DLocalBasis()
00042 {
00043 s[0] = 1; s[1] = 1; s[2] = 1;
00044 }
00045
00050 EdgeS02DLocalBasis(unsigned int orientations)
00051 {
00052 s[0] = 1; s[1] = 1; s[2] = 1;
00053 for(int i = 0; i < 3; ++i)
00054 if(orientations & (1<<i)) s[i] *= -1;
00055 }
00056
00058 unsigned int size () const
00059 {
00060 return 3;
00061 }
00062
00064
00086 template<typename Geometry>
00087 inline void evaluateFunctionGlobal
00088 (const typename Traits::DomainType& in,
00089 std::vector<typename Traits::RangeType>& out,
00090 const Geometry &geometry) const
00091 {
00092 typename Traits::DomainType pos = geometry.global(in);
00093
00094 FieldVector<typename Traits::DomainFieldType, 3> coeff[3];
00095
00096 coefficientsGlobal(coeff, geometry);
00097
00098 out.resize(3);
00099 for(int i = 0; i < 3; ++i) {
00100
00101 out[i][0] = s[i]*( coeff[i][alpha]*pos[1] + coeff[i][a0]);
00102 out[i][1] = s[i]*(-coeff[i][alpha]*pos[0] + coeff[i][a1]);
00103 }
00104 }
00105
00107
00131 template<typename Geometry>
00132 inline void evaluateJacobianGlobal
00133 (const typename Traits::DomainType& in,
00134 std::vector<typename Traits::JacobianType>& out,
00135 const Geometry &geometry) const
00136 {
00137
00138 FieldVector<typename Traits::DomainFieldType, 3> coeff[3];
00139
00140 coefficientsGlobal(coeff, geometry);
00141
00142 out.resize(3);
00143 for(int i = 0; i < 3; ++i) {
00144
00145 out[i][0][0] = 0; out[i][0][1] = s[i]*coeff[i][alpha];
00146 out[i][1][0] = -s[i]*coeff[i][alpha]; out[i][1][1] = 0;
00147 }
00148 }
00149
00151 unsigned int order () const
00152 {
00153 return 1;
00154 }
00155
00156 private:
00157 template<typename Geometry>
00158 void coefficientsGlobal
00159 (FieldVector<typename Traits::DomainFieldType, 3> (&coeff)[3],
00160 const Geometry &geometry) const
00161 {
00162 assert(geometry.type().isTriangle());
00163
00201 typename Traits::DomainType vertex[3];
00202 for(int i = 0; i < 3; ++i)
00203 vertex[i] = geometry.corner(i);
00204
00205 typename Traits::DomainType offset[3];
00206 offset[0] = vertex[1]; offset[0] -= vertex[0];
00207 offset[1] = vertex[2]; offset[1] -= vertex[0];
00208 offset[2] = vertex[2]; offset[2] -= vertex[1];
00209
00210
00211 FieldMatrix<typename Traits::DomainFieldType, 3, 3> M;
00212 M[0][0] = offset[2][0]*vertex[0][1]-offset[2][1]*vertex[0][0];
00213 M[0][1] = offset[2][0];
00214 M[0][2] = offset[2][1];
00215
00216 M[1][0] = offset[1][0]*vertex[1][1]-offset[1][1]*vertex[1][0];
00217 M[1][1] = offset[1][0];
00218 M[1][2] = offset[1][1];
00219
00220 M[2][0] = offset[0][0]*vertex[2][1]-offset[0][1]*vertex[2][0];
00221 M[2][1] = offset[0][0];
00222 M[2][2] = offset[0][1];
00223
00224 M.invert();
00225
00226 FieldVector<typename Traits::DomainFieldType, 3> rhs;
00227
00228
00229 rhs[0] = -offset[0].two_norm();
00230 rhs[1] = offset[0].two_norm();
00231 rhs[2] = 0;
00232 M.mv(rhs, coeff[0]);
00233
00234
00235 rhs[0] = offset[1].two_norm();
00236 rhs[1] = 0;
00237 rhs[2] = offset[1].two_norm();
00238 M.mv(rhs, coeff[1]);
00239
00240
00241 rhs[0] = 0;
00242 rhs[1] = offset[2].two_norm();
00243 rhs[2] = -offset[2].two_norm();
00244 M.mv(rhs, coeff[2]);
00245 }
00246
00247
00248 static const typename FieldVector<typename Traits::DomainFieldType, 3>::size_type
00249 alpha = 0, a0 = 1, a1 = 2;
00250
00252 R s[3];
00253 };
00254 }
00255 #endif // DUNE_EDGES02DLOCALBASIS_HH