common/referenceelements.hh

Go to the documentation of this file.
00001 // $Id: referenceelements.hh 4116 2008-04-14 09:24:23Z christi $
00002 
00003 #ifndef DUNE_REFERENCEELEMENTS_HH
00004 #define DUNE_REFERENCEELEMENTS_HH
00005 
00006 #include <iostream>
00007 #include <dune/common/fvector.hh>
00008 #include <dune/common/exceptions.hh>
00009 #include <dune/common/geometrytype.hh>
00010 #include <dune/grid/common/grid.hh>
00011 #include <dune/common/misc.hh>
00012 
00019 namespace Dune
00020 {
00021   namespace {
00022 
00024     template <class RT, int codim>
00025     FieldVector<typename RT::CoordType, RT::d> 
00026     mapGlobal(const RT& refElem,
00027               const FieldVector<typename RT::CoordType, RT::d-codim>& l, 
00028               int i, int cdim)
00029     {
00030       typedef typename RT::CoordType CoordType;
00031       const int dim = RT::d;
00032       
00033       assert(cdim == codim);
00034       static FieldMatrix<CoordType, dim-codim, dim> mat;
00035       
00036       int zeroLocalIdx = refElem.subEntity(i, cdim, 0, dim);
00037       FieldVector<CoordType, dim> result = 
00038         refElem.position(zeroLocalIdx, dim);
00039       for (int j = 0; j < dim-cdim; ++j) {
00040         int indexLocal = refElem.subEntity(i, cdim, j+1, dim);
00041         mat[j] = refElem.position(indexLocal, dim) - result;
00042       }
00043 
00044       mat.umtv(l, result);
00045 
00046       return result;
00047     }
00048     
00049   } // end unnamed namespace
00050 
00051 
00097   /***********************************************************
00098    * Interface for reference elements 
00099    ***********************************************************/
00100 
00102   template<typename ctype, int dim>
00103   class ReferenceElement
00104   {
00105   public:
00106 
00107     // compile time sizes
00108     enum { d=dim };    // maps from R^d
00109 
00110     // export types
00111     typedef ctype CoordType;
00112 
00114     virtual int size (int c) const = 0;
00115         
00117     virtual int size (int i, int c, int cc) const = 0;
00118         
00120     virtual int subEntity (int i, int c, int ii, int cc) const = 0;
00121 
00123     virtual const FieldVector<ctype,dim>& position (int i, int c) const = 0;
00124  
00126     template <int codim>
00127     FieldVector<ctype, dim> global(const FieldVector<ctype, dim-codim>& local, int i, int cdim) const 
00128     {
00129       return mapGlobal<ReferenceElement<ctype, dim>, codim>(*this, local, i, cdim);
00130     }
00131 
00133     virtual GeometryType type (int i, int c) const = 0;
00134 
00136     virtual double volume () const = 0;
00137 
00139     virtual ~ReferenceElement ()
00140     {}
00141   };
00142 
00144   template<typename Imp>
00145   class ReferenceElementWrapper : 
00146     public ReferenceElement<typename Imp::CoordType,Imp::d>,
00147     private Imp
00148   {
00149   public:
00150 
00151     // compile time sizes
00152     enum { d=Imp::d };    // maps from R^d
00153 
00154     // export types
00155     typedef typename Imp::CoordType CoordType;
00156 
00158     int size (int c) const
00159     {
00160       return Imp::size(c);
00161     }
00162         
00164     int size (int i, int c, int cc) const
00165     {
00166       return Imp::size(i,c,cc);
00167     }
00168         
00170     int subEntity (int i, int c, int ii, int cc) const
00171     {
00172       return Imp::subEntity(i,c,ii,cc);
00173     }
00174 
00176     const FieldVector<CoordType,d>& position (int i, int c) const
00177     {
00178       return Imp::position(i,c);
00179     }
00180 
00182     template <int codim>
00183     FieldVector<CoordType, d> global(const FieldVector<CoordType, d-codim>& l, 
00184                                      int i, int cdim) const 
00185     {
00186       return Imp::global(l, i, cdim);
00187     }
00188 
00190     GeometryType type (int i, int c) const
00191     {
00192       return Imp::type(i,c);
00193     }
00194 
00196     double volume () const
00197     {
00198       return Imp::volume();
00199     }
00200   };
00201 
00202   // output operator for wrapped reference elements
00203   template<typename T>
00204   inline std::ostream& operator<< (std::ostream& s, 
00205                                    const ReferenceElementWrapper<T>& r)
00206   {
00207     enum {dim=T::d};
00208 
00209     std::cout << "REFERENCE ELEMENT " << GeometryName(r.type(0,0))
00210               << " dimension=" << dim 
00211               << " volume=" << r.volume()
00212               << std::endl;
00213 
00214     for (int c=0; c<=dim; c++)
00215       {
00216         std::cout << r.size(c) << " codim " << c << " entitie(s)" << std::endl;
00217         for (int i=0; i<r.size(c); i++)
00218           {
00219             std::cout << "  entity=" << i
00220                       << " codim=" << c
00221                       << " type=" << GeometryName(r.type(i,c))
00222                       << " position=(" << r.position(i,c) << ")"
00223                       << std::endl;
00224 
00225             for (int cc=c+1; cc<=dim; cc++)
00226               {
00227                 std::cout << "    " << r.size(i,c,cc) 
00228                           << " subentities of codim " << cc
00229                           << std::endl;
00230 
00231                 for (int ii=0; ii<r.size(i,c,cc); ii++)
00232                   {
00233                     std::cout << "    index=" << ii
00234                               << " subentity=" << r.subEntity(i,c,ii,cc)
00235                               << " position=(" << r.position(r.subEntity(i,c,ii,cc),cc) << ")"
00236                               << std::endl;
00237                   }
00238               }
00239           }
00240       }
00241 
00242     return s;
00243   }
00244 
00245 
00246   /***********************************************************
00247    * The hypercube in any dimension
00248    ***********************************************************/
00249 
00250 
00252   template<typename ctype, int dim>
00253   class ReferenceCube
00254   {
00255   public:
00256     enum { MAXE = Power_m_p<3,dim>::power }; // maximum number of entities per codim
00257 
00258     // compile time sizes
00259     enum { d=dim };    // maps from R^d
00260 
00261     // export types
00262     typedef ctype CoordType;
00263 
00265     ReferenceCube ()
00266     {
00267       for (int i=0; i<=dim; ++i)
00268         sizes[i] = 0;
00269       for (int i=0; i<MAXE; ++i) 
00270         for (int j=0; j<=dim; ++j) 
00271           for (int k=0; k<=dim; ++k)
00272                         if (j==k)
00273                           subsizes[i][j][k] = 1;
00274                         else
00275                           subsizes[i][j][k] = 0;
00276       for (int i=0; i<MAXE; ++i) 
00277         for (int j=0; j<=dim; ++j) 
00278                   hierarchy[i][j][0][j] = i;
00279       FieldVector<int,dim> direction;
00280       for (int c=dim; c>=0; --c)
00281         generate(0,c,direction);
00282 
00283 
00284     }
00285         
00287     int size (int c) const
00288     {
00289       return sizes[c];
00290     }
00291         
00293     int size (int i, int c, int cc) const
00294     {
00295       return subsizes[i][c][cc];
00296     }
00297         
00299     int subEntity (int i, int c, int ii, int cc) const
00300     {
00301       return hierarchy[i][c][ii][cc];
00302     }
00303 
00305     const FieldVector<ctype,dim>& position (int i, int c) const
00306     {
00307       return pos[i][c];   
00308     }
00309 
00311     template <int codim>
00312     FieldVector<ctype, dim> global(const FieldVector<ctype, dim-codim>& local, int i, int cdim) const 
00313     {
00314       return mapGlobal<ReferenceCube<ctype, dim>, codim>(*this, local, i, cdim);
00315     }
00316 
00318     GeometryType type (int i, int c) const
00319     {
00320         return GeometryType(GeometryType::cube, dim-c);   
00321     }
00322 
00324     double volume () const
00325     {
00326       return 1.0;
00327     }
00328 
00330     const FieldVector<int,dim>& iposition (int i, int c) const
00331     {
00332       return ipos[i][c];          
00333     }
00334 
00335   private:
00336 
00337     class IdMapper {
00338     public:
00339       int& operator() (const FieldVector<int,dim>& x)
00340       {
00341         int index=x[dim-1];
00342         for (int i=dim-2; i>=0; --i) index = 3*index+x[i];
00343         return id[index];
00344       }
00345     private:
00346       int id[1<<(2*dim)];
00347     };
00348 
00349     void generate (int k, int c, FieldVector<int,dim>& direction)
00350     {
00351       if (k<c) 
00352         { 
00353           // select kth direction
00354           for (int i=0; i<dim; ++i)
00355             {
00356               bool done=false;
00357               for (int j=0; j<k; ++j)
00358                 if (i<=direction[j]) {
00359                   done = true;
00360                   break;
00361                 }
00362               if (done) continue;
00363               direction[k] = i; // new direction selected
00364               generate(k+1,c,direction);
00365             }
00366         }
00367       else
00368         {
00369           //              std::cout << "c=" << c 
00370           //                                    << " directions=(" << direction << ")"
00371           //                                    << std::endl;
00372                                 
00373           // c directions have been selected
00374           // for each there are 2 choices, ie 2^c possibilities in total
00375           for (int b=0; b<(1<<c); ++b)
00376             {
00377               // make coordinate in dim-cube
00378               FieldVector<int,dim> x; 
00379               for (int i=0; i<dim; ++i) x[i] = 1;
00380               for (int i=0; i<c; i++)
00381                 if (((1<<i)&b)==0)
00382                   x[direction[i]] = 0;
00383                 else
00384                   x[direction[i]] = 2;
00385 
00386               int entity = sizes[c];
00387               (sizes[c])++;
00388               if (sizes[c]>MAXE)
00389                 DUNE_THROW(GridError, "MAXE in ReferenceCube exceeded");
00390 
00391               // print info
00392               //                          std::cout << " x=(" << x << ")->"
00393               //                                                << "id=" << entity
00394               //                                                << std::endl;
00395 
00396               // store id in map
00397               idmap(x) = entity;
00398 
00399               // assign position
00400               for (int i=0; i<dim; i++)
00401                 pos[entity][c][i] = x[i]*0.5;
00402 
00403               // assign integer position
00404               ipos[entity][c] = x;
00405 
00406               // generate subentities
00407               for (int cc=c+1; cc<=dim; ++cc)
00408                 generatesub(k,cc,direction,x,c);
00409             }
00410         }
00411     }
00412 
00413     void generatesub (int k, int cc, FieldVector<int,dim>& direction,
00414                       FieldVector<int,dim>& e, int c)
00415     {
00416       if (k<cc) 
00417         { 
00418           // select kth direction
00419           for (int i=0; i<dim; ++i)
00420             {
00421               bool done=false;
00422               for (int j=0; j<c; ++j)
00423                 if (i==direction[j]) {
00424                   done = true;
00425                   break;
00426                 }
00427               for (int j=c; j<k; ++j)
00428                 if (i<=direction[j]) {
00429                   done = true;
00430                   break;
00431                 }
00432               if (done) continue;
00433               direction[k] = i; // new direction selected
00434               generatesub(k+1,cc,direction,e,c);
00435             }
00436         }
00437       else
00438         {
00439           //              std::cout << "  cc=" << cc 
00440           //                                    << " directions=(" << direction << ")"
00441           //                                    << std::endl;
00442                                 
00443           // cc-c additional directions have been selected
00444           // for each there are 2 choices, ie 2^(cc-c) possibilities in total
00445           for (int b=0; b<(1<<(cc-c)); ++b)
00446             {
00447               // make coordinate in dim-cube
00448               FieldVector<int,dim> x(e); 
00449               for (int i=0; i<(cc-c); i++)
00450                 if (((1<<i)&b)==0)
00451                   x[direction[i+c]] = 0;
00452                 else
00453                   x[direction[i+c]] = 2;
00454 
00455               int entity = idmap(e);
00456               int subentity = idmap(x);
00457               int index = subsizes[entity][c][cc];
00458               (subsizes[entity][c][cc])++;
00459               if (subsizes[entity][c][cc]>MAXE)
00460                 DUNE_THROW(GridError, "MAXE in ReferenceCube exceeded");
00461 
00462               // print info
00463               //                          std::cout << "  (" << e << ")," << entity
00464               //                                                << " has subentity (" << x << ")," << subentity
00465               //                                                << " at index=" << index
00466               //                                                << std::endl;
00467 
00468               // store id
00469               hierarchy[entity][c][index][cc] = subentity;
00470             }
00471         }
00472     }
00473 
00474     IdMapper idmap;
00475     int sizes[dim+1];
00476     int subsizes[MAXE][dim+1][dim+1];
00477     int hierarchy[MAXE][dim+1][MAXE][dim+1];
00478     FieldVector<ctype,dim> pos[MAXE][dim+1];
00479     FieldVector<int,dim> ipos[MAXE][dim+1];
00480   };
00481 
00482 
00484   template<typename ctype>
00485   class ReferenceCube<ctype,0>
00486   {
00487   public:
00488     enum {MAXE = 1};
00489     enum {d=0};
00490    
00491     typedef ctype CoordType;
00492     
00493     ReferenceCube ()
00494     {
00495         }
00496 
00498     int size (int c) const
00499     {
00500           return 1;
00501     }
00502 
00504     int size (int i, int c, int cc) const
00505     {
00506           return 1;
00507     }
00508 
00510     int subEntity (int i, int c, int ii, int cc) const
00511     {
00512           return 0;
00513     }
00514     
00516     const FieldVector<ctype,d>& position (int i, int c) const
00517     {
00518       return pos;         
00519     }
00520 
00522     template <int codim>
00523     FieldVector<ctype, d> global(const FieldVector<ctype, d-codim>& local, int i, int cdim) const 
00524     {
00525       return pos;
00526     }
00527 
00529     GeometryType type (int i, int c) const
00530     {
00531           return GeometryType(GeometryType::simplex,d-c);         
00532     }
00533     
00535     double volume () const 
00536         { 
00537           return 1; // This must be 1 to make *= volume work!
00538         }
00539 
00540   private:
00541     
00542         FieldVector<ctype,d> pos;
00543   };
00544 
00545 
00547   template<typename ctype, int dim>
00548   class ReferenceCubeContainer
00549   {
00550   public:
00551 
00553     typedef ReferenceCube<ctype,dim> value_type;
00554 
00556     const value_type& operator() (GeometryType type) const
00557     {
00558         if ( type.isCube() )
00559             return cube_;
00560       DUNE_THROW(RangeError, "expected a cube!");
00561     }
00562 
00563   private:
00564     ReferenceCube<ctype,dim> cube_;
00565   };
00566 
00567  
00568 
00569  
00570   /***********************************************************
00571    * the simplex in any dimension (line,triangle,tetrahedron)
00572    ***********************************************************/
00573 
00574 
00575   /*
00576 
00577   see the reference elements:
00578 
00579   y
00580   | 2(0,1)
00581   | |\
00582   | | \
00583   | |  \
00584   | |   \
00585   | |    \ 
00586   | |_ _ _\  
00587   |  0(0,0)  1(1,0)
00588   |_ _ _ _ _ _ _ _ 
00589   x
00590 
00591 
00592   linear triangular 
00593   area=1/2;
00594   ---------------------
00595 
00596 
00597   3 (0,0,1)
00598 
00599   |`  
00600   |\ `
00601   | \  `
00602   |  \   `  
00603   |   \  .' 2 (0,1,0)
00604   |    .' |
00605   |  .' \ | 
00606   |.'_ _ \|  
00607   
00608   0       1
00609   (0,0,0,)  (1,0,0)
00610 
00611   linear tetrahedron
00612   volume = 1/6;
00613   --------------------
00614   */
00615 
00616 
00617   //reference simplex without virtual functions
00619   template<typename ctype, int dim>
00620   class ReferenceSimplex
00621   {
00622   public:
00623     enum {MAXE = 6}; // implement only up to 3D anyway, avoids some warnings
00624     enum {d=dim};
00625    
00626 
00627     typedef ctype CoordType;
00628     
00629     ReferenceSimplex () : volume_( (double) 1.0/Factorial<dim>::factorial )
00630     {
00631           assert(dim<=3);
00632       for (int i=0; i<=dim; ++i)
00633                 sizes[i]=0;
00634       for (int i=0; i<MAXE; ++i) 
00635                 for (int j=0; j<=dim; ++j) 
00636                   for (int k=0; k<=dim; ++k)
00637                         if (j==k)
00638                           subsizes[i][j][k] = 1;
00639                         else
00640                           subsizes[i][j][k] = 0;
00641           for (int i=0; i<MAXE; ++i) 
00642         for (int j=0; j<=dim; ++j) 
00643                   subentityindex[i][j][0][j] = i;
00644      
00645       for (int c=dim; c>=0; --c)
00646                 entity_details (c);
00647     }
00648 
00649     
00650     
00652     int size (int c) const
00653     {
00654       //entity_size(c);
00655       return sizes[c];
00656     }
00658     int size (int i, int c, int cc) const
00659     {
00660       return subsizes[i][c][cc];
00661     }
00663     int subEntity (int i, int c, int ii, int cc) const
00664     {
00665       return subentityindex[i][c][ii][cc];
00666     }
00667     
00669     const FieldVector<ctype,dim>& position (int i, int c) const
00670     {
00671       return pos[i][c];   
00672     }
00673 
00675     template <int codim>
00676     FieldVector<ctype, dim> global(const FieldVector<ctype, dim-codim>& local, int i, int cdim) const 
00677     {
00678       return mapGlobal<ReferenceSimplex<ctype, dim>, codim>(*this, local, i, cdim);
00679     }
00680 
00682     GeometryType type (int i, int c) const
00683     {
00684         return GeometryType(GeometryType::simplex,dim-c);         
00685     }
00686     
00688     double volume () const { return volume_; }
00689 
00691     //   const FieldVector<int,dim>& iposition (int i, int c) const
00692     //     {
00693     //       //return ipos[i][c];         
00694     //     }
00695   private:
00696     
00697     //details of entities and subentities
00698     void entity_details(int c) //(int c) 
00699     {
00700      
00701       sizes[dim]=dim+1; // simplex definition 
00702          
00703       // position of vertices, there are dim+1 vertices
00704       FieldVector<int,dim> x(0);
00705       
00706       // vertex is codim=dim entity
00707       for (int n=0;n<dim;n++)
00708         {
00709           pos[0][dim][n]=x[n];
00710             
00711         }
00712       for(int k=1;k<=dim;++k)
00713         {
00714           for (int j=0;j<dim;++j)
00715             {
00716               x[j]=0;
00717               x[k-1]=1;
00718               pos[k][dim][j]= x[j]; 
00719                 
00720             }
00721         }
00722      
00723       // position of centre of gravity of the element 
00724       // codim=0 for element or cell
00725       sizes[0]=1; // only 1 cell !!
00726       int node;
00727       for(int k=0;k<dim;++k)
00728         { node=0;
00729         pos[sizes[0]-1][0][k]=(pos[0][dim][k])/sizes[dim];
00730         subentityindex[sizes[0]-1][0][0][dim]=0;
00731         node+=1;
00732         for (int j=1;j<sizes[dim];++j)
00733           {
00734             pos[sizes[0]-1][0][k]+=(pos[j][dim][k])/(sizes[dim]);
00735             subentityindex[sizes[0]-1][0][j][dim]=j;
00736             node+=1;
00737           }
00738         }
00739       subsizes[sizes[0]-1][0][dim]=node;
00740 
00741       // the simplex itself has one simplex ;)
00742       subsizes[0][0][0]=1;
00743       
00744       // every vertex entity has one sub entity --> itself 
00745       for (int k=0;k<MAXE;++k) subsizes[k][dim][dim]=1; 
00746 
00747       //++++++++++++++++
00748       if(dim==1)// line
00749         {
00750           // node indices on element
00751           for(int i=0;i<subsizes[0][0][dim];++i)
00752             subentityindex[0][0][i][dim]=i;
00753         }
00754           else if(dim==2) // triangle
00755         {
00756                   sizes[dim-1]=3; // edge
00757             
00758           // hard coding the number of subentities
00759           // triangle has 3 vertices, 3 edges 
00760           subsizes[0][0][dim]=3; 
00761           subsizes[0][0][dim-1]=3; 
00762           // triangle  has 2 vertices on each  edge
00763           for (int k=0;k<3;++k){
00764             subsizes[k][1][dim]=2;
00765                         // triangle  has 1 edge on each edge ;-)
00766                         subsizes[k][1][1]=1;
00767           }
00768           // subentity indices
00769           // node indices on element
00770           for(int i=0;i<subsizes[0][0][dim];++i)
00771             subentityindex[0][0][i][dim]=i;
00772           // edge indices on element
00773           for(int i=0;i<subsizes[0][0][dim-1];++i)
00774             subentityindex[0][0][i][1]=i;
00775           
00776                   // node indices on edge 0
00777           subentityindex[0][1][0][dim]=1;
00778           subentityindex[0][1][1][dim]=2;
00779           // node indices on edge 1
00780           subentityindex[1][1][0][dim]=2;
00781           subentityindex[1][1][1][dim]=0;
00782           // node indices on edge 2
00783                   // 2 = (subsizes[0][0][dim-1])-1  work around icc warning
00784           subentityindex[(subsizes[0][0][dim-1])-1][dim-1][0][dim]=0;
00785           subentityindex[(subsizes[0][0][dim-1])-1][dim-1][1][dim]=1;
00786 
00787           for(int j=0;j<dim;++j)
00788             {
00789 
00790                           //edge 0 (nodes 1,2)
00791                           // 2 = (subsizes[0][0][dim])-1 work around icc warning
00792               pos[0][1][j]=(pos[1][dim][j]+pos[(subsizes[0][0][dim])-1][dim][j])/2;
00793               //edge 1 (nodes 0,2)
00794               pos[1][1][j]=(pos[0][dim][j]+pos[(subsizes[0][0][dim])-1][dim][j])/2;
00795               //edge 2 (nodes 0,1)
00796               pos[(subsizes[0][0][dim-1])-1][1][j]=(pos[0][dim][j]+pos[1][2][j])/2;
00797                           
00798               // //edge 0 (nodes 1,2)
00799 //               pos[0][1][j]=(pos[1][dim][j]+pos[2][dim][j])/2;
00800 //               //edge 1 (nodes 0,2)
00801 //               pos[1][1][j]=(pos[0][dim][j]+pos[2][dim][j])/2;
00802 //               //edge 2 (nodes 0,1)
00803 //               pos[(subsizes[0][0][dim-1])-1][1][j]=(pos[0][dim][j]+pos[1][2][j])/2;
00804             }
00805         }
00806       else if(dim==3)// tetrahedron
00807         {
00808           sizes[1]=4; // face
00809           sizes[dim-1]=6; // edge
00810 
00811           // hard coding the number of subentities
00812           // tetrahedron has 4 vertices, 6 edges and 4 facese on element 
00813           subsizes[0][0][dim]=4; 
00814           subsizes[0][0][dim-1]=6; 
00815           subsizes[0][0][1]=4; 
00816           //  tetrahedron has 3 vertices on each triang. face 
00817           for(int i=0;i<subsizes[0][0][1];++i)
00818             subsizes[i][1][dim]=3;
00819           //  tetrahedron has 3 edges on each triang. face 
00820           for(int i=0;i<subsizes[0][0][1];++i)
00821             subsizes[i][1][dim-1]=3;
00822                   //  tetrahedron has 1 face on each triang. face! 
00823                   for(int i=0;i<subsizes[0][0][1];++i)
00824             subsizes[i][1][1]=1;
00825       
00826           //  tetrahedron has 2 vertices on each edge 
00827           //  we have 6 edges, an we know that
00828           for (int k=0;k<6;++k)
00829             subsizes[k][dim-1][dim]=2;
00830           
00831                   //  tetrahedron has 1 edge on each edge! 
00832           //  we have 6 edges, an we know that
00833           for (int k=0;k<6;++k)
00834             subsizes[k][dim-1][dim-1]=1;
00835           
00836           // subentity indices
00837           // node indices on element
00838           for(int i=0;i<subsizes[0][0][dim];++i)
00839             subentityindex[0][0][i][dim]=i;
00840           // edge indices on element
00841           for(int i=0;i<subsizes[0][0][dim-1];++i)
00842             subentityindex[0][0][i][dim-1]=i; 
00843           // face indices on element
00844           for(int i=0;i<subsizes[0][0][1];++i)
00845             subentityindex[0][0][i][1]=i; 
00846 
00847           // node numbering- normal pointing outward
00848       // reference triangle in dune has a counter-clockwise
00849           // numbering. so  nodes on each triangle face in tetrahedron numbered same way
00850                     int faceindx;
00851                         int edgeindx;
00852           // node indices on face 0
00853                         subentityindex[0][1][0][dim]=1;
00854                         subentityindex[0][1][1][dim]=2;
00855                         subentityindex[0][1][(subsizes[0][1][dim]-1)][dim]=3;
00856           // node indices on face 1
00857           subentityindex[1][1][0][dim]=0;
00858           subentityindex[1][1][1][dim]=3;
00859           subentityindex[1][1][(subsizes[1][1][dim]-1)][dim]=2;
00860           // node indices on face 2
00861                   faceindx= subsizes[0][0][1]-2;
00862           subentityindex[faceindx][1][0][dim]=0;
00863           subentityindex[faceindx][1][1][dim]=1;
00864           subentityindex[faceindx][1][(subsizes[faceindx][1][dim]-1)][dim]=3;
00865           // node indices on face 3
00866                   faceindx= subsizes[0][0][1]-1;
00867           subentityindex[faceindx][1][0][dim]=0;
00868           subentityindex[faceindx][1][1][dim]=2;
00869           subentityindex[faceindx][1][(subsizes[faceindx][1][dim]-1)][dim]=1;
00870                   
00871           // edge indices on face 0
00872           subentityindex[0][1][0][dim-1]=5;
00873           subentityindex[0][1][1][dim-1]=4;
00874                   
00875                   edgeindx=subsizes[0][1][dim-1]-1;
00876           subentityindex[0][1][edgeindx][dim-1]=1;  
00877           // edge indices on face 1
00878           subentityindex[1][1][0][dim-1]=5;
00879           subentityindex[1][1][1][dim-1]=2;
00880                   
00881                   edgeindx=subsizes[1][1][dim-1]-1;
00882           subentityindex[1][1][edgeindx][dim-1]=3;
00883           // edge indices on face 2
00884                   faceindx= subsizes[0][0][1]-2;
00885           subentityindex[faceindx][1][0][dim-1]=4;
00886           subentityindex[faceindx][1][1][dim-1]=3;
00887                   
00888                   edgeindx=subsizes[faceindx][1][dim-1]-1;
00889           subentityindex[faceindx][1][edgeindx][dim-1]=0;
00890           // edge indices on face 3
00891                    faceindx= subsizes[0][0][1]-1;
00892           subentityindex[faceindx][1][0][dim-1]=1;
00893           subentityindex[faceindx][1][1][dim-1]=0;
00894                   
00895                    edgeindx=subsizes[faceindx][1][dim-1]-1;
00896           subentityindex[faceindx][1][edgeindx][dim-1]=2;
00897 
00898           // node indices on edge 0
00899           subentityindex[0][dim-1][0][dim]=0;
00900           subentityindex[0][dim-1][1][dim]=1;
00901           // node indices on edge 1
00902           subentityindex[1][dim-1][0][dim]=1;
00903           subentityindex[1][dim-1][1][dim]=2;
00904           // node indices on edge 2
00905                   
00906           subentityindex[(subsizes[0][0][dim-1]-4)][dim-1][0][dim]=0;
00907           subentityindex[(subsizes[0][0][dim-1]-4)][dim-1][1][dim]=2;
00908           // node indices on edge 3
00909           subentityindex[(subsizes[0][0][dim-1]-3)][dim-1][0][dim]=0;
00910           subentityindex[(subsizes[0][0][dim-1]-3)][dim-1][1][dim]=3;
00911           // node indices on edge 4
00912           subentityindex[(subsizes[0][0][dim-1]-2)][dim-1][0][dim]=1;
00913           subentityindex[(subsizes[0][0][dim-1]-2)][dim-1][1][dim]=3;
00914           // node indices on edge 5
00915           subentityindex[(subsizes[0][0][dim-1]-1)][dim-1][0][dim]=2;
00916           subentityindex[(subsizes[0][0][dim-1]-1)][dim-1][1][dim]=3;
00917 
00918           for(int j=0;j<dim;++j)
00919             {
00920 
00921               //face 0 (nodes 1,2,3)
00922               pos[0][1][j]=(pos[1][dim][j]+pos[2][dim][j]+pos[3][dim][j])/3;
00923               //face 1 (nodes 0,2,3)
00924               pos[1][1][j]=(pos[0][dim][j]+pos[2][dim][j]+pos[3][dim][j])/3;
00925               //face 2 (nodes 0,1,3)
00926               pos[2][1][j]=(pos[0][dim][j]+pos[1][dim][j]+pos[3][dim][j])/3;
00927               //face 3 (nodes 0,1,2)
00928               pos[3][1][j]=(pos[0][dim][j]+pos[1][dim][j]+pos[2][dim][j])/3;
00929 
00930               //edge 0 (nodes 0,1)
00931               pos[0][2][j]=(pos[0][dim][j]+pos[1][dim][j])/2;
00932               //edge 1 (nodes 1,2)
00933               pos[1][2][j]=(pos[1][dim][j]+pos[2][dim][j])/2;
00934               //edge 2 (nodes 0,2)
00935               pos[2][2][j]=(pos[0][dim][j]+pos[2][dim][j])/2;
00936               //edge 3 (nodes 0,3)
00937               pos[3][2][j]=(pos[0][dim][j]+pos[3][dim][j])/2;
00938               //edge 4 (nodes 1,3)
00939               pos[4][2][j]=(pos[1][dim][j]+pos[3][dim][j])/2;
00940               //edge 5 (nodes 2,3)
00941               pos[5][2][j]=(pos[2][dim][j]+pos[3][dim][j])/2;
00942             }
00943         }
00944       else
00945         DUNE_THROW(NotImplemented, "dim not implemented yet");
00946 
00947       //++++++++++++++++
00948            
00949     }
00950 
00951     int sizes[dim+1];
00952     int subsizes[MAXE][dim+1][dim+1];
00953     int subentityindex[MAXE][dim+1][MAXE][dim+1];
00954         FieldVector<ctype,dim> pos[MAXE][dim+1];
00955 
00956     // volume of reference simplex 
00957     const double volume_;
00958   };
00959 
00960 
00962   template<typename ctype>
00963   class ReferenceSimplex<ctype,0>
00964   {
00965   public:
00966     enum {MAXE = 1};
00967     enum {d=0};
00968    
00969     typedef ctype CoordType;
00970     
00971     ReferenceSimplex ()
00972     {
00973         }
00974 
00976     int size (int c) const
00977     {
00978           return 1;
00979     }
00980 
00982     int size (int i, int c, int cc) const
00983     {
00984           return 1;
00985     }
00986 
00988     int subEntity (int i, int c, int ii, int cc) const
00989     {
00990           return 0;
00991     }
00992     
00994     const FieldVector<ctype,d>& position (int i, int c) const
00995     {
00996       return pos;         
00997     }
00998 
01000     template <int codim>
01001     FieldVector<ctype, d> global(const FieldVector<ctype, d-codim>& local, int i, int cdim) const 
01002     {
01003       return pos;
01004     }
01005 
01007     GeometryType type (int i, int c) const
01008     {
01009           return GeometryType(GeometryType::simplex,d-c);         
01010     }
01011     
01013     double volume () const 
01014         { 
01015           return 1; 
01016         }
01017 
01018   private:
01019     
01020         FieldVector<ctype,d> pos;
01021   };
01022 
01023 
01025   template<typename ctype, int dim>
01026   class ReferenceSimplexContainer
01027   {
01028   public:
01029 
01031     typedef ReferenceSimplex<ctype,dim> value_type;
01032 
01034     const value_type& operator() (GeometryType type) const
01035     {
01036         if ( type.isSimplex())
01037         return simplices;
01038       DUNE_THROW(RangeError, "expected a simplex!");
01039     }
01040 
01041   private:
01042     ReferenceSimplex<ctype,dim> simplices;
01043   };
01044 
01045 
01046   //++++++++++++++++++++++++++++++++++++++++
01047   //Reference Prism
01048   //++++++++++++++++++++++++++++++++++++++++
01049 
01050   /*                                   5(0,1,1)  
01051                        ¸.¸
01052                        y            ¸ *  . ` ,
01053                        |           *_________ *4(1,0,1)
01054                        |   3(0,0,1)|     .    |
01055                        |           |     |    |
01056                        |           |     .    |
01057                        |           |     |    |
01058                        |           |   ¸2¸(0,1,0)  
01059                        |           | ¸`    *¸ |
01060                        |           |*________*| 
01061                        |    0(0,0,0)         1(1,0,0)
01062                        |_ _ _ _ _ _ _ _ 
01063                        x
01064 
01065                        volume of triangular prism =1/2 
01066 
01067   */
01068 
01069 
01070   template<typename ctype, int dim>
01071   class ReferencePrism;
01072 
01073   template<typename ctype>
01074   class ReferencePrism<ctype,3>
01075   {
01076   public:
01077     //assert dim=3;// dim must be 3
01078     enum{dim=3};
01079     enum {MAXE=9}; // 9 edges
01080     enum{d=dim};
01081     typedef ctype CoordType;
01082 
01083 
01084     ReferencePrism()
01085     {
01086       for (int i=0; i<=3; ++i)
01087                 sizes[i]=0;
01088 
01089       for (int i=0; i<MAXE; ++i) 
01090                 for (int j=0; j<=dim; ++j) 
01091                   for (int k=0; k<=dim; ++k)
01092                         if (j==k)
01093                           subsizes[i][j][k] = 1;
01094                         else
01095                           subsizes[i][j][k] = 0;
01096           for (int i=0; i<MAXE; ++i) 
01097         for (int j=0; j<=dim; ++j) 
01098                   subentityindex[i][j][0][j] = i;
01099 
01100       for (int c=3; c>=0; --c)
01101                 prism_entities (c);
01102     }
01103 
01105     int size (int c) const
01106     {
01107       return sizes[c];
01108     }
01109 
01111     int size (int i, int c, int cc) const
01112     {
01113       return subsizes[i][c][cc];
01114     }
01115 
01117     int subEntity (int i, int c, int ii, int cc) const
01118     {
01119       return subentityindex[i][c][ii][cc];
01120     }
01121 
01123     const FieldVector<ctype,dim>& position (int i, int c) const
01124     {
01125       return pos[i][c];   
01126     }
01127 
01129     template <int codim>
01130     FieldVector<ctype, dim> global(const FieldVector<ctype, dim-codim>& local, int i, int cdim) const 
01131     {
01132       return mapGlobal<ReferencePrism<ctype, dim>, codim>(*this, local, i, cdim);
01133     }
01134 
01136     GeometryType type (int i, int c) const
01137     {
01138       switch (c) 
01139                 {
01140                 case 3: return GeometryType(GeometryType::cube,0);
01141                 case 2: return GeometryType(GeometryType::cube,1);
01142                 case 0: return GeometryType(GeometryType::prism,3);
01143                 case 1: 
01144                   switch (i)
01145                         {
01146                         case 0: return GeometryType(GeometryType::simplex,2);
01147                         case 1: return GeometryType(GeometryType::cube,2);
01148                         case 2: return GeometryType(GeometryType::cube,2);
01149                         case 3: return GeometryType(GeometryType::cube,2);
01150                         case 4: return GeometryType(GeometryType::simplex,2);
01151                         default:
01152                           DUNE_THROW(RangeError, "i argument out of range");
01153                         }
01154                 default:
01155                   DUNE_THROW(RangeError, "codim argument out of range");
01156                 }
01157     }
01158     
01160     
01161     double volume () const
01162     {
01163       double vol=1.0/2.0;
01164       return vol;
01165       
01166     }
01167   private:
01168 
01169     void prism_entities(int c)
01170     {
01171       // compile time error if dim is not equal to 3
01172       dune_static_assert(dim == 3, "prism required dim==3");
01173       // hard coding the size of entities
01174       sizes[0]=1; // element
01175       sizes[1]=5; // face
01176       sizes[2]=9; // edge
01177       sizes[3]=6; // vertex
01178 
01179       //-----------------------------------------------
01180       // hard coding the number of subentities
01181       // prism has 6 vertices, 9 edges and 5 facese on element 
01182       subsizes[0][0][3]=6; 
01183       subsizes[0][0][2]=9; 
01184       subsizes[0][0][1]=5; 
01185 
01186       // a prism itself  has one prism ;) 
01187       subsizes[0][0][0]=1; 
01188  
01189       // face indices according to that given in 
01190       // http://dune-project.org/doc/appl/refelements.html
01191           
01192       // prism has one face on each face!
01193       for(int ii=0;ii<5;++ii)
01194         subsizes[ii][1][1]=1;
01195       //  prism has 3 vertices on bott triang. face 
01196       subsizes[0][1][3]=3;
01197       //  prism has 4 vertices on front,right and left  rectang. faces
01198       subsizes[1][1][3]=4; 
01199       subsizes[2][1][3]=4; 
01200       subsizes[3][1][3]=4;
01201           //  prism has 3 vertices on top triang. face 
01202       subsizes[4][1][3]=3;
01203           
01204       // prism has 3 edges on a bott triang. face
01205       subsizes[0][1][2]=3; 
01206       // prism has 4 edges on  front,right and left rec. faces 
01207       subsizes[1][1][2]=4; 
01208       subsizes[2][1][2]=4; 
01209       subsizes[3][1][2]=4; 
01210       // prism has 3 edges on a top triang. face
01211       subsizes[4][1][2]=3; 
01212       // prism has 2 vertices on each  edge
01213           // also 1 edge on each edge
01214       for (int k=0;k<9;++k){
01215                 subsizes[k][2][3]=2;
01216                 subsizes[k][2][2]=1;
01217           }
01218           
01219       //------------------------------------------
01220 
01221       // positions of vertex with local index "i", there are 6 vertices
01222       // here c = codim = 3
01223       // -----------------------------------
01224       FieldVector<int,3> x;
01225       x=0;
01226       for (int n=0;n<3;n++)
01227         {
01228           pos[0][3][n]=x[n]; 
01229             
01230         }
01231       for(int k=1;k<=5;++k)
01232         for (int j=0;j<3;++j)
01233           {
01234             if(k>3)
01235               {
01236                 x[j]=0;
01237                 x[k-4]=1;
01238                 x[2]=1;
01239                 pos[k][3][j]= x[j];
01240               }
01241             else
01242               {
01243                 x[j]=0;
01244                 x[k-1]=1;
01245                 pos[k][3][j]= x[j];
01246               } 
01247           }
01248       //---------------------------------------
01249 
01250       // position of centre of gravity of the element 
01251       for(int k=0;k<3;++k)
01252         {
01253           pos[sizes[0]-1][0][k]=(pos[0][3][k])/sizes[3];
01254           for (int j=1;j<sizes[3];++j)
01255             pos[sizes[0]-1][0][k]+=(pos[j][3][k])/(sizes[3]);
01256         }
01257         
01258       // subentity indices - in agreement with UG description
01259       // see ug-cvs/gm/elements.c 
01260       
01261       // node indices on element
01262       for(int i=0;i<subsizes[0][0][3];++i)
01263         subentityindex[0][0][i][3]=i;
01264       // edge indices on element
01265       for(int i=0;i<subsizes[0][0][2];++i)
01266         subentityindex[0][0][i][2]=i; 
01267       // face indices on element
01268       for(int i=0;i<subsizes[0][0][1];++i)
01269         subentityindex[0][0][i][1]=i; 
01270 
01271       // node indices on face 0
01272       subentityindex[0][1][0][3]=0;
01273       subentityindex[0][1][1][3]=2;
01274       subentityindex[0][1][2][3]=1;
01275       // node indices on face 1
01276       subentityindex[1][1][0][3]=0;
01277       subentityindex[1][1][1][3]=1;
01278       subentityindex[1][1][2][3]=4;
01279       subentityindex[1][1][3][3]=3;
01280       // node indices on face 2
01281       subentityindex[2][1][0][3]=1;
01282       subentityindex[2][1][1][3]=2;
01283       subentityindex[2][1][2][3]=5;
01284       subentityindex[2][1][3][3]=4;
01285       // node indices on face 3
01286       subentityindex[3][1][0][3]=2;
01287       subentityindex[3][1][1][3]=0;
01288       subentityindex[3][1][2][3]=3;
01289       subentityindex[3][1][3][3]=5;
01290       // node indices on face 4
01291       for(int i=0;i<subsizes[4][1][3];++i) 
01292         subentityindex[4][1][i][3]=i+3;
01293 
01294       // edge indices on face 0
01295       subentityindex[0][1][0][2]=2;
01296       subentityindex[0][1][1][2]=1;
01297       subentityindex[0][1][2][2]=0;
01298       // edge indices on face 1
01299       subentityindex[1][1][0][2]=0;
01300       subentityindex[1][1][1][2]=4;
01301       subentityindex[1][1][2][2]=6;
01302       subentityindex[1][1][3][2]=3;
01303       // edge indices on face 2
01304       subentityindex[2][1][0][2]=1;
01305       subentityindex[2][1][1][2]=5;
01306       subentityindex[2][1][2][2]=7;
01307       subentityindex[2][1][3][2]=4;
01308       // edge indices on face 3
01309       subentityindex[3][1][0][2]=2;
01310       subentityindex[3][1][1][2]=3;
01311       subentityindex[3][1][2][2]=8;
01312       subentityindex[3][1][3][2]=5;
01313       // edge indices on face 4
01314       subentityindex[4][1][0][2]=6;
01315       subentityindex[4][1][1][2]=7;
01316       subentityindex[4][1][2][2]=8;
01317 
01318       // node indices on edge 0
01319       subentityindex[0][2][0][3]=0;
01320       subentityindex[0][2][1][3]=1;
01321       // node indices on edge 1
01322       subentityindex[1][2][0][3]=1;
01323       subentityindex[1][2][1][3]=2;
01324       // node indices on edge 2
01325       subentityindex[2][2][0][3]=2;
01326       subentityindex[2][2][1][3]=0;
01327       // node indices on edge 3
01328       subentityindex[3][2][0][3]=0;
01329       subentityindex[3][2][1][3]=3;
01330       // node indices on edge 4
01331       subentityindex[4][2][0][3]=1;
01332       subentityindex[4][2][1][3]=4;
01333       // node indices on edge 5
01334       subentityindex[5][2][0][3]=2;
01335       subentityindex[5][2][1][3]=5;
01336       // node indices on edge 6
01337       subentityindex[6][2][0][3]=3;
01338       subentityindex[6][2][1][3]=4;
01339       // node indices on edge 7
01340       subentityindex[7][2][0][3]=4;
01341       subentityindex[7][2][1][3]=5;
01342       // node indices on edge 8
01343       subentityindex[8][2][0][3]=5;
01344       subentityindex[8][2][1][3]=3;
01345       //
01346 
01347       //position of faces and edges
01348       for(int j=0;j<3;++j)
01349         {
01350           //face 0 (nodes 0,1,2)
01351           pos[0][1][j]=(pos[0][3][j]+pos[1][3][j]+pos[2][3][j])/3;
01352           //face 1 (nodes 0,1,3,4)
01353           pos[1][1][j]=(pos[0][3][j]+pos[1][3][j]+pos[3][3][j]+pos[4][3][j])/4;
01354           //face 2 (nodes 1,2,4,5)
01355           pos[2][1][j]=(pos[1][3][j]+pos[2][3][j]+pos[4][3][j]+pos[5][3][j])/4;
01356           //face 3 (nodes 0,2,3,5)
01357           pos[3][1][j]=(pos[0][3][j]+pos[2][3][j]+pos[3][3][j]+pos[5][3][j])/4;
01358           //face 4 (nodes 3,4,5)
01359           pos[4][1][j]=(pos[3][3][j]+pos[4][3][j]+pos[5][3][j])/3;
01360 
01361           //edge 0 (nodes 0,1)
01362           pos[0][2][j]=(pos[0][3][j]+pos[1][3][j])/2;
01363           //edge 1 (nodes 1,2)
01364           pos[1][2][j]=(pos[1][3][j]+pos[2][3][j])/2;
01365           //edge 2 (nodes 2,0)
01366           pos[2][2][j]=(pos[2][3][j]+pos[0][3][j])/2;
01367           //edge 3 (nodes 0,3)
01368           pos[3][2][j]=(pos[0][3][j]+pos[3][3][j])/2;
01369           //edge 4 (nodes 1,4)
01370           pos[4][2][j]=(pos[1][3][j]+pos[4][3][j])/2;
01371           //edge 5 (nodes 2,5)
01372           pos[5][2][j]=(pos[2][3][j]+pos[5][3][j])/2;
01373           //edge 6 (nodes 3,4)
01374           pos[6][2][j]=(pos[3][3][j]+pos[4][3][j])/2;
01375           //edge 7 (nodes 4,5)
01376           pos[7][2][j]=(pos[4][3][j]+pos[5][3][j])/2;
01377           //edge 8 (nodes 5,3)
01378           pos[8][2][j]=(pos[5][3][j]+pos[3][3][j])/2;
01379 
01380         }
01381 
01382 
01383     } 
01384 
01385 
01386     int sizes[dim+1];
01387     int subsizes[MAXE][dim+1][dim+1];
01388     int subentityindex[MAXE][dim+1][MAXE][dim+1];
01389     FieldVector<ctype,dim> pos[MAXE][dim+1];
01390   
01391 
01392   };
01393 
01395   template<typename ctype, int dim>
01396   class ReferencePrismContainer
01397   {
01398   public:
01399 
01401     typedef ReferencePrism<ctype,dim> value_type;
01402 
01404     const value_type& operator() (GeometryType type) const
01405     {
01406         if (type.isPrism())
01407         return pris;
01408       DUNE_THROW(RangeError, "expected a prism!");
01409     }
01410 
01411   private:
01412     ReferencePrism<ctype,dim> pris;
01413   };
01414 
01415 
01416   //++++++++++++++++++++++++++++++++++++
01417   // Reference Pyramid
01418   //++++++++++++++++++++++++++++++++++++
01419 
01420   template<typename ctype, int dim>
01421   class ReferencePyramid;
01422 
01423   template<typename ctype>
01424   class ReferencePyramid<ctype, 3>
01425   {
01426  
01427   public:
01428     enum{dim=3};
01429     enum {MAXE=8}; // 8 edges
01430     enum{d=dim};
01431     typedef ctype CoordType;
01432 
01433 
01434     ReferencePyramid()
01435     {
01436       for (int i=0; i<=3; ++i)
01437                 sizes[i]=0;
01438           
01439       for (int i=0; i<MAXE; ++i) 
01440         for (int j=0; j<=dim; ++j) 
01441           for (int k=0; k<=dim; ++k)
01442                         if (j==k)
01443                           subsizes[i][j][k] = 1;
01444                         else
01445                           subsizes[i][j][k] = 0;
01446           for (int i=0; i<MAXE; ++i) 
01447         for (int j=0; j<=dim; ++j) 
01448                   subentityindex[i][j][0][j] = i;
01449           
01450       for (int c=3; c>=0; --c)
01451                 pyramid_entities (c);
01452     }
01453 
01455     int size (int c) const
01456     {
01457       return sizes[c];
01458     }
01459 
01461     int size (int i, int c, int cc) const
01462     {
01463       return subsizes[i][c][cc];
01464     }
01465 
01467     int subEntity (int i, int c, int ii, int cc) const
01468     {
01469       return subentityindex[i][c][ii][cc];
01470     }
01471 
01473     const FieldVector<ctype,dim>& position (int i, int c) const
01474     {
01475       return pos[i][c];   
01476     }
01477 
01479     template <int codim>
01480     FieldVector<ctype, dim> global(const FieldVector<ctype, dim-codim>& local, int i, int cdim) const 
01481     {
01482       return mapGlobal<ReferencePyramid<ctype, dim>, codim>(*this, local, i, cdim);
01483     }
01484 
01486     GeometryType type (int i, int c) const
01487     {
01488       switch (c) 
01489                 {
01490                 case 3: return GeometryType(GeometryType::cube,0);
01491                 case 2: return GeometryType(GeometryType::cube,1);
01492                 case 0: return GeometryType(GeometryType::pyramid,3);
01493                 case 1: 
01494                   switch (i)
01495                         {
01496                         case 0: return GeometryType(GeometryType::cube,2);
01497                         case 1: 
01498                         case 2: 
01499                         case 3: 
01500                         case 4: return GeometryType(GeometryType::simplex,2);
01501                         default:
01502                           DUNE_THROW(RangeError, "i argument out of range");
01503                         }
01504                 default:
01505                   DUNE_THROW(RangeError, "codim argument out of range");
01506                 }
01507     }
01508     
01510     
01511     double volume () const
01512     {
01513       double vol=1.0/3.0;
01514       return vol;
01515       
01516     }
01517   private:
01518 
01519     void pyramid_entities(int c)
01520   
01521     {
01522       // compile time error if dim is not 3
01523       dune_static_assert(dim == 3, "pyramid required dim==3");
01524       // hard coding the size of entities
01525       sizes[0]=1; // element
01526       sizes[1]=5; // face
01527       sizes[2]=8; // edge
01528       sizes[3]=5; // vertex
01529 
01530       //-----------------------------------------------
01531       // hard coding the number of subentities
01532       // pyramid has 5 vertices, 8 edges and 5 facese on element 
01533       subsizes[0][0][3]=5; 
01534       subsizes[0][0][2]=8; 
01535       subsizes[0][0][1]=5; 
01536       
01537       // a pyramid itself has one pyramid ;) 
01538       subsizes[0][0][0]=1; 
01539       // pyramid has one face on each face!
01540       for(int ii=0;ii<5;++ii)
01541         subsizes[ii][1][1]=1;
01542 
01543       //  pyramid has 4 vertices on bott rect. face 
01544       subsizes[0][1][3]=4; 
01545       //  pyramid has 3 vertices on front,right,back and left triang. faces
01546       subsizes[1][1][3]=3; 
01547       subsizes[2][1][3]=3; 
01548       subsizes[3][1][3]=3;
01549       subsizes[4][1][3]=3; 
01550      
01551       // pyramid has 4 edges on a bott rect. face
01552       subsizes[0][1][2]=4; 
01553       // pyramid has 3 edges on front,right,back and left triang. faces
01554       for(int i=1;i<5;++i)
01555         subsizes[i][1][2]=3; 
01556      
01557       // pyramid has 2 vertices on each  edge
01558       for (int k=0;k<8;++k){
01559         subsizes[k][2][3]=2;
01560         //also one edge on each edge!!
01561         subsizes[k][2][2]=1;
01562           }
01563       //------------------------------------------
01564 
01565       // positions of vertex with local index "i", there are 5 vertices
01566       // here c = codim = 3
01567       // -----------------------------------
01568      
01569      
01570       FieldVector<int,3> x;
01571       x=0;
01572       for (int n=0;n<3;n++)
01573         {
01574           pos[0][3][n]=x[n]; 
01575             
01576         }
01577 
01578       for(int k=1;k<=4;++k)
01579         for (int j=0;j<3;++j)
01580           {
01581             if(k==1)
01582               {
01583                 x[j]=0;
01584                 x[k-1]=1;
01585                 pos[k][3][j]= x[j];
01586               }
01587             else if(k==2)
01588               {
01589                 x[j]=0;
01590                 x[k-2]=1;
01591                 x[1]=1;
01592                 pos[k][3][j]= x[j];
01593               }
01594             else
01595               {
01596                 x[j]=0;
01597                 x[k-2]=1;
01598                 pos[k][3][j]= x[j];
01599               } 
01600           }
01601 
01602       //---------------------------------------
01603 
01604       // position of centre of gravity of the element 
01605       for(int k=0;k<3;++k)
01606         {
01607           pos[sizes[0]-1][0][k]=(pos[0][3][k])/sizes[3];
01608           for (int j=1;j<sizes[3];++j)
01609             pos[sizes[0]-1][0][k]+=(pos[j][3][k])/(sizes[3]);
01610         }
01611         
01612  // subentity indices - in agreement with UG description
01613       // see ug-cvs/gm/elements.c 
01614 
01615       // node indices on element
01616       for(int i=0;i<subsizes[0][0][3];++i)
01617         subentityindex[0][0][i][3]=i;
01618       // edge indices on element
01619       for(int i=0;i<subsizes[0][0][2];++i)
01620         subentityindex[0][0][i][2]=i; 
01621       // face indices on element
01622       for(int i=0;i<subsizes[0][0][1];++i)
01623         subentityindex[0][0][i][1]=i; 
01624 
01625       // node indices on face 0
01626       subentityindex[0][1][0][3]=0;
01627       subentityindex[0][1][1][3]=3;
01628       subentityindex[0][1][2][3]=2;
01629       subentityindex[0][1][3][3]=1;
01630       // node indices on face 1
01631       subentityindex[1][1][0][3]=0;
01632       subentityindex[1][1][1][3]=1;
01633       subentityindex[1][1][2][3]=4;
01634       // node indices on face 2
01635       subentityindex[2][1][0][3]=1;
01636       subentityindex[2][1][1][3]=2;
01637       subentityindex[2][1][2][3]=4;
01638       // node indices on face 3
01639       subentityindex[3][1][0][3]=2;
01640       subentityindex[3][1][1][3]=3;
01641       subentityindex[3][1][2][3]=4;
01642       // node indices on face 4
01643       subentityindex[4][1][0][3]=3;
01644       subentityindex[4][1][1][3]=0;
01645       subentityindex[4][1][2][3]=4;
01646 
01647       // edge indices on face 0
01648       subentityindex[0][1][0][2]=3;
01649       subentityindex[0][1][1][2]=2;
01650       subentityindex[0][1][2][2]=1;
01651       subentityindex[0][1][3][2]=0;
01652       // edge indices on face 1
01653       subentityindex[1][1][0][2]=0;
01654       subentityindex[1][1][1][2]=5;
01655       subentityindex[1][1][2][2]=4;
01656       // edge indices on face 2
01657       subentityindex[2][1][0][2]=1;
01658       subentityindex[2][1][1][2]=6;
01659       subentityindex[2][1][2][2]=5;
01660       // edge indices on face 3
01661       subentityindex[3][1][0][2]=2;
01662       subentityindex[3][1][1][2]=7;
01663       subentityindex[3][1][2][2]=6;
01664       // edge indices on face 4
01665       subentityindex[4][1][0][2]=3;
01666       subentityindex[4][1][1][2]=4;
01667       subentityindex[4][1][2][2]=7;
01668 
01669       // node indices on edge 0
01670       subentityindex[0][2][0][3]=0;
01671       subentityindex[0][2][1][3]=1;
01672       // node indices on edge 1
01673       subentityindex[1][2][0][3]=1;
01674       subentityindex[1][2][1][3]=2;
01675       // node indices on edge 2
01676       subentityindex[2][2][0][3]=2;
01677       subentityindex[2][2][1][3]=3;
01678       // node indices on edge 3
01679       subentityindex[3][2][0][3]=3;
01680       subentityindex[3][2][1][3]=0;
01681       // node indices on edge 4
01682       subentityindex[4][2][0][3]=0;
01683       subentityindex[4][2][1][3]=4;
01684       // node indices on edge 5
01685       subentityindex[5][2][0][3]=1;
01686       subentityindex[5][2][1][3]=4;
01687       // node indices on edge 6
01688       subentityindex[6][2][0][3]=2;
01689       subentityindex[6][2][1][3]=4;
01690       // node indices on edge 7
01691       subentityindex[7][2][0][3]=3;
01692       subentityindex[7][2][1][3]=4;
01693       
01694 
01695       //position of faces and edges
01696       for(int j=0;j<3;++j)
01697         {
01698           //face 0 (nodes 0,1,2,3)
01699           pos[0][1][j]=(pos[0][3][j]+pos[1][3][j]+pos[2][3][j]+pos[3][3][j])/4;
01700           //face 1 (nodes 0,1,4)
01701           pos[1][1][j]=(pos[0][3][j]+pos[1][3][j]+pos[4][3][j])/3;
01702           //face 2 (nodes 1,2,4)
01703           pos[2][1][j]=(pos[1][3][j]+pos[2][3][j]+pos[4][3][j])/3;
01704           //face 3 (nodes 2,3,4)
01705           pos[3][1][j]=(pos[2][3][j]+pos[3][3][j]+pos[4][3][j])/3;
01706           //face 4 (nodes 3,0,4)
01707           pos[4][1][j]=(pos[3][3][j]+pos[0][3][j]+pos[4][3][j])/3;
01708 
01709           //edge 0 (nodes 0,1)
01710           pos[0][2][j]=(pos[0][3][j]+pos[1][3][j])/2;
01711           //edge 1 (nodes 1,2)
01712           pos[1][2][j]=(pos[1][3][j]+pos[2][3][j])/2;
01713           //edge 2 (nodes 2,3)
01714           pos[2][2][j]=(pos[2][3][j]+pos[3][3][j])/2;
01715           //edge 3 (nodes 3,0)
01716           pos[3][2][j]=(pos[3][3][j]+pos[0][3][j])/2;
01717           //edge 4 (nodes 0,4)
01718           pos[4][2][j]=(pos[0][3][j]+pos[4][3][j])/2;
01719           //edge 5 (nodes 1,4)
01720           pos[5][2][j]=(pos[1][3][j]+pos[4][3][j])/2;
01721           //edge 6 (nodes 2,4)
01722           pos[6][2][j]=(pos[2][3][j]+pos[4][3][j])/2;
01723           //edge 7 (nodes 3,4)
01724           pos[7][2][j]=(pos[3][3][j]+pos[4][3][j])/2;
01725           
01726 
01727         }
01728 
01729 
01730     } 
01731 
01732 
01733     int sizes[dim+1];
01734     int subsizes[MAXE][dim+1][dim+1];
01735     int subentityindex[MAXE][dim+1][MAXE][dim+1];
01736     FieldVector<ctype,dim> pos[MAXE][dim+1];
01737   
01738 
01739   };
01740 
01742   template<typename ctype, int dim>
01743   class ReferencePyramidContainer
01744   {
01745   public:
01746 
01748     typedef ReferencePrism<ctype,dim> value_type;
01749 
01751     const value_type& operator() (GeometryType type) const
01752     {
01753         if (type.isPyramid())
01754         return pyram;
01755       DUNE_THROW(RangeError, "expected a pyramid!");
01756     }
01757 
01758   private:
01759     ReferencePyramid<ctype,dim> pyram;
01760   };
01761 
01762 
01763 
01764   /***********************************************************
01765    * The general container and the singletons
01766    ***********************************************************/
01767 
01768 
01770   template<typename ctype, int dim>
01771   class ReferenceElementContainer
01772   {
01773   public:
01774 
01776     typedef ReferenceElement<ctype,dim> value_type;
01777 
01779     const ReferenceElement<ctype,dim>& operator() (GeometryType type) const
01780     {
01781         if ( type.isCube())
01782         return hcube;
01783         else if ( type.isSimplex() )
01784         return simplices;
01785       else
01786         DUNE_THROW(NotImplemented, "type not implemented yet");
01787     }
01788 
01789   private:
01790     ReferenceElementWrapper<ReferenceCube<ctype,dim> > hcube;
01791     ReferenceElementWrapper<ReferenceSimplex<ctype,dim> > simplices;
01792   };
01793 
01795   template<typename ctype>
01796   class ReferenceElementContainer<ctype, 3>
01797   {
01798     enum { dim=3 };
01799   public:
01800 
01802     typedef ReferenceElement<ctype,dim> value_type;
01803 
01805     const ReferenceElement<ctype,dim>& operator() (GeometryType type) const
01806     {
01807         if ( type.isCube() )
01808         return hcube;
01809         else if( type.isSimplex() )
01810         return simplices;
01811         else if ( type.isPrism() )
01812         return pris;
01813         else if( type.isPyramid() )
01814         return pyram;
01815       else
01816         DUNE_THROW(NotImplemented, "type not implemented yet");
01817     }
01818 
01819   private:
01820     ReferenceElementWrapper<ReferenceCube<ctype,dim> > hcube;
01821     ReferenceElementWrapper<ReferenceSimplex<ctype,dim> > simplices;
01822     ReferenceElementWrapper<ReferencePrism<ctype,dim> > pris;
01823     ReferenceElementWrapper<ReferencePyramid<ctype,dim> > pyram;
01824   };
01825 
01827   template<typename ctype, int dim>
01828   struct ReferenceElements {
01829     static ReferenceCubeContainer<ctype,dim> cube;
01830     static ReferenceSimplexContainer<ctype,dim> simplices;
01831     static ReferenceElementContainer<ctype,dim> general;
01832   };
01833 
01835   template<typename ctype>
01836   struct ReferenceElements<ctype,3> {
01837       static ReferenceCubeContainer<ctype,3> cube;
01838       static ReferenceSimplexContainer<ctype,3> simplices;
01839       static ReferencePrismContainer<ctype,3> prism;
01840       static ReferencePyramidContainer<ctype,3> pyramid;
01841       static ReferenceElementContainer<ctype,3> general;
01842   };
01843 
01845 }
01846 #endif

Generated on Sun Nov 15 22:28:43 2009 for dune-grid by  doxygen 1.5.6