edges03dlocalbasis.hh
Go to the documentation of this file.00001
00002 #ifndef DUNE_EDGES03DLOCALBASIS_HH
00003 #define DUNE_EDGES03DLOCALBASIS_HH
00004
00005 #include <cmath>
00006
00007 #include <dune/common/fmatrix.hh>
00008
00009 #include <dune/grid/common/genericreferenceelements.hh>
00010
00011 #include <dune/localfunctions/common/localbasis.hh>
00012
00013 namespace Dune
00014 {
00030 template<class D, class R>
00031 class EdgeS03DLocalBasis
00032 {
00033 public:
00035 typedef LocalBasisTraits<
00036 D,3,Dune::FieldVector<D,3>,
00037 R,3,Dune::FieldVector<R,3>,
00038 Dune::FieldMatrix<R,3,3>
00039 > Traits;
00040
00042 EdgeS03DLocalBasis()
00043 {
00044 for(int i = 0; i < 6; ++i)
00045 s[i] = 1;
00046 }
00047
00052 EdgeS03DLocalBasis(unsigned int orientations)
00053 {
00054 for(int i = 0; i < 6; ++i)
00055 s[i] = 1;
00056 for(int i = 0; i < 6; ++i)
00057 if(orientations & (1<<i)) s[i] *= -1;
00058 }
00059
00061 unsigned int size () const
00062 {
00063 return 6;
00064 }
00065
00067
00089 template<typename Geometry>
00090 inline void evaluateFunctionGlobal
00091 (const typename Traits::DomainType& in,
00092 std::vector<typename Traits::RangeType>& out,
00093 const Geometry &geometry) const
00094 {
00095 typename Traits::DomainType pos = geometry.global(in);
00096
00097 FieldVector<typename Traits::DomainFieldType, 6> coeff[6];
00098
00099 coefficientsGlobal(coeff, geometry);
00100
00101 out.resize(6);
00102 for(int i = 0; i < 6; ++i) {
00103
00104 out[i][0] = s[i]*( coeff[i][A01]*pos[1]+coeff[i][A02]*pos[2]+coeff[i][a0]);
00105 out[i][1] = s[i]*(-coeff[i][A01]*pos[0]+coeff[i][A12]*pos[2]+coeff[i][a1]);
00106 out[i][2] = s[i]*(-coeff[i][A02]*pos[0]-coeff[i][A12]*pos[1]+coeff[i][a2]);
00107 }
00108 }
00109
00111
00133 template<typename Geometry>
00134 inline void evaluateJacobianGlobal
00135 (const typename Traits::DomainType& in,
00136 std::vector<typename Traits::JacobianType>& out,
00137 const Geometry &geometry) const
00138 {
00139
00140 FieldVector<typename Traits::DomainFieldType, 6> coeff[6];
00141
00142 coefficientsGlobal(coeff, geometry);
00143
00144 out.resize(6);
00145 for(int i = 0; i < 6; ++i) {
00146
00147 out[i][0][0] = 0;
00148 out[i][0][1] = s[i]*coeff[i][A01];
00149 out[i][0][2] = s[i]*coeff[i][A02];
00150
00151 out[i][1][0] = -s[i]*coeff[i][A01];
00152 out[i][1][1] = 0;
00153 out[i][1][2] = s[i]*coeff[i][A12];
00154
00155 out[i][2][0] = -s[i]*coeff[i][A02];
00156 out[i][2][1] = -s[i]*coeff[i][A12];
00157 out[i][2][2] = 0;
00158 }
00159 }
00160
00162 unsigned int order () const
00163 {
00164 return 1;
00165 }
00166
00167 private:
00168 template<typename Geometry>
00169 void coefficientsGlobal
00170 (FieldVector<typename Traits::DomainFieldType, 6> (&coeff)[6],
00171 const Geometry &geometry) const
00172 {
00173 static const GenericReferenceElement<typename Traits::DomainFieldType, 3> &refElem
00174 = GenericReferenceElements<typename Traits::DomainFieldType, 3>::simplex();
00175
00176 assert(geometry.type().isTetrahedron());
00177
00178 typename Traits::DomainType vertex[4];
00179 for(int i = 0; i < 4; ++i)
00180 vertex[i] = geometry.corner(i);
00181
00182
00183 FieldMatrix<typename Traits::DomainFieldType, 6, 6> M;
00184 typename Traits::DomainType distance[6];
00185 for(int j = 0; j < 6; ++j) {
00186 int v0 = refElem.subEntity(j, 2, 0, 3);
00187 int v1 = refElem.subEntity(j, 2, 1, 3);
00188 if(v0 > v1) std::swap(v0, v1);
00189
00190 distance[j] = vertex[v1]; distance[j] -= vertex[v0];
00191
00192 M[j][0] = vertex[v1][0]*vertex[v0][1] - vertex[v1][1]*vertex[v0][0];
00193 M[j][1] = vertex[v1][0]*vertex[v0][2] - vertex[v1][2]*vertex[v0][0];
00194 M[j][2] = vertex[v1][1]*vertex[v0][2] - vertex[v1][2]*vertex[v0][1];
00195 M[j][3] = distance[j][0];
00196 M[j][4] = distance[j][1];
00197 M[j][5] = distance[j][2];
00198 }
00199
00200 M.invert();
00201
00202 for(int i = 0; i < 6; ++i) {
00203 FieldVector<typename Traits::DomainFieldType, 6> rhs(0);
00204 rhs[i] = distance[i].two_norm();
00205
00206 M.mv(rhs, coeff[i]);
00207 }
00208 }
00209
00210
00211 enum { A01, A02, A12, a0, a1, a2 };
00212
00213 R s[6];
00214 };
00215 }
00216 #endif // DUNE_EDGES03DLOCALBASIS_HH