dune-localfunctions  2.2.0
brezzidouglasmarini12dlocalbasis.hh
Go to the documentation of this file.
00001 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1TRIANGLELOCALBASIS_HH
00002 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1TRIANGLELOCALBASIS_HH
00003 
00004 #include <vector>
00005 
00006 #include <dune/common/fmatrix.hh>
00007 
00008 #include "../../common/localbasis.hh"
00009 
00010 namespace Dune 
00011 {
00021   template<class D, class R>
00022   class BDM12DLocalBasis 
00023   {
00024 
00025 public:
00026     typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
00027                              Dune::FieldMatrix<R,2,2> > Traits;
00028 
00030     BDM12DLocalBasis ()
00031     {
00032       sign0 = sign1 = sign2 = 1.0;
00033     }
00034 
00040     BDM12DLocalBasis (unsigned int s)
00041     {
00042       sign0 = sign1 = sign2 = 1.0;
00043       if (s & 1)
00044       {
00045         sign0 = -1.0;
00046       }
00047       if (s & 2)
00048       {
00049         sign1 = -1.0;
00050       }
00051       if (s & 4)
00052       {
00053         sign2 = -1.0;
00054       }
00055     }
00056 
00058     unsigned int size () const
00059     {
00060       return 6;
00061     }
00062 
00069     inline void evaluateFunction (const typename Traits::DomainType& in,
00070                                   std::vector<typename Traits::RangeType>& out) const
00071     { 
00072       out.resize(6);
00073       
00074       out[0][0] = sign0*in[0];
00075       out[0][1] = sign0*(in[1] - 1.0);
00076       out[1][0] = sign1*(in[0] - 1.0);
00077       out[1][1] = sign1*in[1];
00078       out[2][0] = sign2*in[0];
00079       out[2][1] = sign2*in[1];
00080       out[3][0] = 3.0*in[0];
00081       out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
00082       out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
00083       out[4][1] = -3.0*in[1];
00084       out[5][0] = -3.0*in[0];
00085       out[5][1] = 3.0*in[1];
00086     }
00087 
00094     inline void evaluateJacobian (const typename Traits::DomainType& in,
00095                                   std::vector<typename Traits::JacobianType>& out) const
00096     {  
00097       out.resize(6);
00098       
00099       out[0][0][0] = sign0;
00100       out[0][0][1] = 0.0; 
00101       out[0][1][0] = 0.0;
00102       out[0][1][1] = sign0;
00103 
00104       out[1][0][0] = sign1;
00105       out[1][0][1] = 0.0;
00106       out[1][1][0] = 0.0;
00107       out[1][1][1] = sign1;
00108 
00109       out[2][0][0] = sign2;
00110       out[2][0][1] = 0.0;
00111       out[2][1][0] = 0.0;
00112       out[2][1][1] = sign2;
00113     
00114       out[3][0][0] = 3.0;
00115       out[3][0][1] = 0.0;
00116       out[3][1][0] = -6.0;
00117       out[3][1][1] = -3.0;
00118   
00119       out[4][0][0] = 3.0;
00120       out[4][0][1] = 6.0;
00121       out[4][1][0] = 0.0;
00122       out[4][1][1] = -3.0;
00123   
00124       out[5][0][0] = -3.0;
00125       out[5][0][1] = 0.0;
00126       out[5][1][0] = 0.0;
00127       out[5][1][1] = 3.0;
00128     }
00129 
00131     unsigned int order () const
00132     {
00133       return 1;
00134     }
00135 
00136   private:
00137     R sign0, sign1, sign2;
00138   };
00139 }
00140 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1TRIANGLELOCALBASIS_HH