hierarchicalprismp2localbasis.hh

Go to the documentation of this file.
00001 // -*- tab-width: 4; indent-tabs-mode: nil -*-
00002 #ifndef DUNE_HIERARCHICAL_PRISM_P2_LOCALBASIS_HH
00003 #define DUNE_HIERARCHICAL_PRISM_P2_LOCALBASIS_HH
00004 
00009 #include <dune/common/fvector.hh>
00010 #include <dune/common/fmatrix.hh>
00011 
00012 #include <dune/localfunctions/common/localbasis.hh>
00013 
00014 namespace Dune
00015 {
00016   template<class D, class R>
00017   class HierarchicalPrismP2LocalBasis
00018   {
00019    public: 
00021    typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,3> > Traits;
00022 
00024 unsigned int size () const
00025 {
00026 return 18;
00027 }
00028 
00030 void evaluateFunction (const typename Traits::DomainType& in,
00031                        std::vector<typename Traits::RangeType> & out) const
00032 {
00033 out.resize(18);
00034 
00035         out[0]=(1.0-in[0]-in[1])*(1.0-in[2]);
00036         out[1]= in[0]*(1-in[2]);
00037         out[2]=in[1]*(1-in[2]);
00038         out[3]=in[2]*(1.0-in[0]-in[1]);
00039         out[4]=in[0]*in[2];
00040         out[5]=in[1]*in[2];
00041 
00042         //edges
00043         out[6]=2*(1.0-in[0]-in[1])*(0.5-in[0]-in[1])*(4*in[2]-4*in[2]*in[2]);   
00044         out[7]=2*in[0]*(-0.5+in[0])*(4*in[2]-4*in[2]*in[2]); 
00045         out[8]=2*in[1]*(-0.5+in[1])*(4*in[2]-4*in[2]*in[2]);
00046         out[9]=4*in[0]*(1-in[0]-in[1])*(1-3*in[2]+2*in[2]*in[2]);   
00047         out[10]=4*in[1]*(1-in[0]-in[1])*(1-3*in[2]+2*in[2]*in[2]);   
00048         out[11]=4*in[0]*in[1]*(1-3*in[2]+2*in[2]*in[2]);   
00049         out[12]=4*in[0]*(1-in[0]-in[1])*(-in[2]+2*in[2]*in[2]);   
00050         out[13]=4*in[1]*(1-in[0]-in[1])*(-in[2]+2*in[2]*in[2]);   
00051         out[14]=4*in[0]*in[1]*(-in[2]+2*in[2]*in[2]);   
00052 
00053         //faces
00054         out[15]=4*in[0]*(1-in[0]-in[1])*(4*in[2]-4*in[2]*in[2]);         
00055         out[16]=4*in[1]*(1-in[0]-in[1])*(4*in[2]-4*in[2]*in[2]);   
00056         out[17]=4*in[0]*in[1]*(4*in[2]-4*in[2]*in[2]);   
00057 }
00058 
00059 
00060 
00062     void evaluateJacobian (const typename Traits::DomainType& in,     //position
00063                            std::vector<typename Traits::JacobianType>& out) const  //return  value
00064 {
00065 out.resize(18);  
00066 
00067         //vertices 
00068         out[0][0][0] = in[2]-1;  
00069         out[0][0][1] = in[2]-1;  
00070         out[0][0][2] = in[0]+in[1]-1;
00071         
00072         out[1][0][0] = 1-in[2];  
00073         out[1][0][1] = 0;         
00074         out[1][0][2] =-in[0]; 
00075         
00076         out[2][0][0] = 0;        
00077         out[2][0][1] = 1-in[2];  
00078         out[2][0][2] = -in[1];
00079         
00080         out[3][0][0] = -in[2];   
00081         out[3][0][1] = -in[2];   
00082         out[3][0][2] = 1-in[0]-in[1];
00083         
00084         out[4][0][0] = in[2];    
00085         out[4][0][1] = 0;         
00086         out[4][0][2] = in[0];
00087         
00088         out[5][0][0] = 0;        
00089         out[5][0][1] = in[2];     
00090         out[5][0][2] = in[1];
00091                 
00092         //edges
00093         out[6][0][0] = (-3+4*in[0]+4*in[1])*(4*in[2]-4*in[2]*in[2]);  
00094         out[6][0][1] = (-3+4*in[0]+4*in[1])*(4*in[2]-4*in[2]*in[2]);  
00095         out[6][0][2] = 2*(1-in[0]-in[1])*(0.5-in[0]-in[1])*(4-8*in[2]);  
00096                 
00097         out[7][0][0] = (-1+4*in[0])*(4*in[2]-4*in[2]*in[2]);  
00098         out[7][0][1] = 0;
00099         out[7][0][2] = 2*in[0]*(-0.5+in[0])*(4-8*in[2]);  
00100         
00101         out[8][0][0] = 0;
00102         out[8][0][1] = (-1+4*in[1])*(4*in[2]-4*in[2]*in[2]);  
00103         out[8][0][2] = 2*in[1]*(-0.5+in[1])*(4-8*in[2]);  
00104         
00105         out[9][0][0] = (4-8*in[0]-4*in[1])*(1-3*in[2]+2*in[2]*in[2]);  
00106         out[9][0][1] = -4*in[0]*(1-3*in[2]+2*in[2]*in[2]);  
00107         out[9][0][2] = 4*in[0]*(1-in[0]-in[1])*(-3+4*in[2]);  
00108         
00109         out[10][0][0] = (-4*in[1])*(1-3*in[2]+2*in[2]*in[2]);  
00110         out[10][0][1] = (4-4*in[0]-8*in[1])*(1-3*in[2]+2*in[2]*in[2]);  
00111         out[10][0][2] = 4*in[1]*(1-in[0]-in[1])*(-3+4*in[2]);  
00112         
00113         out[11][0][0] = 4*in[1]*(1-3*in[2]+2*in[2]*in[2]);  
00114         out[11][0][1] = 4*in[0]*(1-3*in[2]+2*in[2]*in[2]);  
00115         out[11][0][2] = 4*in[0]*in[1]*(-3+4*in[2]);  
00116         
00117         out[12][0][0] = (4-8*in[0]-4*in[1])*(-in[2]+2*in[2]*in[2]);  
00118         out[12][0][1] = (-4*in[0])*(-in[2]+2*in[2]*in[2]);  
00119         out[12][0][2] = 4*in[0]*(1-in[0]-in[1])*(-1+4*in[2]);  
00120         
00121         out[13][0][0] = -4*in[1]*(-in[2]+2*in[2]*in[2]);  
00122         out[13][0][1] = (4-4*in[0]-8*in[1])*(-in[2]+2*in[2]*in[2]);  
00123         out[13][0][2] = 4*in[1]*(1-in[0]-in[1])*(-1+4*in[2]);  
00124         
00125         out[14][0][0] = 4*in[1]*(-in[2]+2*in[2]*in[2]);  
00126         out[14][0][1] = 4*in[0]*(-in[2]+2*in[2]*in[2]);  
00127         out[14][0][2] = 4*in[0]*in[1]*(-1+4*in[2]);  
00128         
00129         //faces
00130         out[15][0][0] = (4-8*in[0]-4*in[1])*(4*in[2]-4*in[2]*in[2]);  
00131         out[15][0][1] = -4*in[0]*(4*in[2]-4*in[2]*in[2]);  
00132         out[15][0][2] = 4*in[0]*(1-in[0]-in[1])*(4-8*in[2]);  
00133         
00134         out[16][0][0] = -4*in[1]*(4*in[2]-4*in[2]*in[2]);  
00135         out[16][0][1] = (4-4*in[0]-8*in[1])*(4*in[2]-4*in[2]*in[2]);  
00136         out[16][0][2] = 4*in[1]*(1-in[0]-in[1])*(4-8*in[2]);  
00137         
00138         out[17][0][0] = 4*in[1]*(4*in[2]-4*in[2]*in[2]);         
00139         out[17][0][1] = 4*in[0]*(4*in[2]-4*in[2]*in[2]);  
00140         out[17][0][2] = 4*in[0]*in[1]*(4-8*in[2]);  
00141          
00142         }
00145 unsigned int order() const
00146 {
00147 return 2;
00148 }
00149 
00150  };
00151 }
00152 #endif
Generated on Sat Apr 24 11:15:33 2010 for dune-localfunctions by  doxygen 1.6.3