yaspgrid.hh

Go to the documentation of this file.
00001 #ifndef DUNE_YASPGRID_HH
00002 #define DUNE_YASPGRID_HH
00003 
00004 #include<iostream>
00005 #include<vector>
00006 #include<map>
00007 #include<algorithm>
00008 #include<stack>
00009 
00010 #include <dune/grid/common/grid.hh>     // the grid base classes
00011 #include <dune/grid/yaspgrid/grids.hh>  // the yaspgrid base classes
00012 #include <dune/grid/common/capabilities.hh> // the capabilities
00013 #include <dune/common/misc.hh>
00014 #include <dune/common/helpertemplates.hh>
00015 #include <dune/common/bigunsignedint.hh>
00016 #include <dune/common/typetraits.hh>
00017 #include <dune/common/collectivecommunication.hh>
00018 #include <dune/grid/common/indexidset.hh>
00019 #include <dune/grid/common/datahandleif.hh>
00020 
00021 
00022 #if HAVE_MPI
00023 #include <dune/common/mpicollectivecommunication.hh>
00024 #endif
00025 
00033 namespace Dune {
00034 
00035 //************************************************************************
00039 typedef double yaspgrid_ctype;
00040 static const yaspgrid_ctype yasptolerance=1E-13; // tolerance in coordinate computations
00041 
00042 /* some sizes for building global ids
00043  */
00044   const int yaspgrid_dim_bits = 24;   // bits for encoding each dimension
00045   const int yaspgrid_level_bits = 6;  // bits for encoding level number
00046   const int yaspgrid_codim_bits = 4;  // bits for encoding codimension
00047 
00048 
00049 //************************************************************************
00050 // forward declaration of templates
00051 
00052 template<int dim, int dimworld>               class YaspGrid;
00053 template<int mydim, int cdim, class GridImp>  class YaspGeometry;
00054 template<int codim, int dim, class GridImp>   class YaspEntity;
00055 template<int codim, class GridImp>            class YaspEntityPointer;
00056 template<int codim, PartitionIteratorType pitype, class GridImp> class YaspLevelIterator;
00057 template<class GridImp>            class YaspIntersectionIterator;
00058 template<class GridImp>            class YaspHierarchicIterator;
00059 template<class GridImp>            class YaspLevelIndexSet;
00060 template<class GridImp>            class YaspLeafIndexSet;
00061 template<class GridImp>            class YaspGlobalIdSet;
00062 
00063 //========================================================================
00071 //========================================================================
00072 
00073 template<int mydim, int cdim, class GridImp>
00074 class YaspSpecialGeometry : public Geometry<mydim, cdim, GridImp, YaspGeometry>
00075 {
00077   typedef typename GridImp::ctype ctype;
00078 public:
00079   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, int& m) :
00080     Geometry<mydim, cdim, GridImp, YaspGeometry>(YaspGeometry<mydim, cdim, GridImp>(p,h,m))
00081     {}
00082   YaspSpecialGeometry() :
00083     Geometry<mydim, cdim, GridImp, YaspGeometry>(YaspGeometry<mydim, cdim, GridImp>(false))
00084     {};
00085 };
00086 
00087 template<int mydim, class GridImp>
00088 class YaspSpecialGeometry<mydim,mydim,GridImp> : public Geometry<mydim, mydim, GridImp, YaspGeometry>
00089 {
00091   typedef typename GridImp::ctype ctype;
00092 public:
00093   YaspSpecialGeometry(const FieldVector<ctype, mydim>& p, const FieldVector<ctype, mydim>& h) :
00094     Geometry<mydim, mydim, GridImp, YaspGeometry>(YaspGeometry<mydim, mydim, GridImp>(p,h))
00095     {}
00096   YaspSpecialGeometry() :
00097     Geometry<mydim, mydim, GridImp, YaspGeometry>(YaspGeometry<mydim, mydim, GridImp>(false))
00098     {};
00099 };
00100 
00101 template<int cdim, class GridImp>
00102 class YaspSpecialGeometry<0,cdim,GridImp> : public Geometry<0, cdim, GridImp, YaspGeometry>
00103 {
00105   typedef typename GridImp::ctype ctype;
00106 public:
00107   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p) :
00108     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(p))
00109     {}
00110   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, int& m) :
00111     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(p))
00112     {}
00113   YaspSpecialGeometry() :
00114     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(false))
00115     {};
00116 };
00117 
00118 //========================================================================
00119 // The transformation describing the refinement rule
00120 
00121 template<int dim, class GridImp>
00122 class YaspFatherRelativeLocalElement {
00123 public:
00124   static FieldVector<yaspgrid_ctype, dim> midpoint;  // data neded for the refelem below
00125   static FieldVector<yaspgrid_ctype, dim> extension; // data needed for the refelem below
00126   static YaspSpecialGeometry<dim,dim,GridImp> geo;
00127   static YaspSpecialGeometry<dim,dim,GridImp>& getson (int i)
00128   {
00129     for (int k=0; k<dim; k++)
00130       if (i&(1<<k))
00131         midpoint[k] = 0.75;
00132       else
00133         midpoint[k] = 0.25;
00134     return geo;
00135   }
00136 };
00137 
00138 // initialize static variable with bool constructor (which makes reference elements)
00139 template<int dim, class GridImp>
00140 YaspSpecialGeometry<dim,dim,GridImp> 
00141 YaspFatherRelativeLocalElement<dim,GridImp>::geo(YaspFatherRelativeLocalElement<dim,GridImp>::midpoint,
00142                                                  YaspFatherRelativeLocalElement<dim,GridImp>::extension);
00143 template<int dim, class GridImp>
00144 FieldVector<yaspgrid_ctype,dim> YaspFatherRelativeLocalElement<dim,GridImp>::midpoint(0.25);
00145 
00146 template<int dim, class GridImp>
00147 FieldVector<yaspgrid_ctype,dim> YaspFatherRelativeLocalElement<dim,GridImp>::extension(0.5);
00148 
00150 template<int mydim,int cdim, class GridImp>
00151 class YaspGeometry : public GeometryDefaultImplementation<mydim,cdim,GridImp,YaspGeometry>
00152 {
00153 public:
00155   typedef typename GridImp::ctype ctype;
00157   typedef Geometry<mydim,mydim,GridImp,Dune::YaspGeometry> ReferenceGeometry;
00158   
00160   GeometryType type () const
00161   {
00162       return GeometryType(GeometryType::cube,mydim);
00163   }
00164 
00166   int corners () const
00167   {
00168         return 1<<mydim;
00169   }
00170 
00172   const FieldVector<ctype, cdim>& operator[] (int i) const
00173   {
00174         assert( i >= 0 && i < (int) coord_.N() ); 
00175         FieldVector<ctype, cdim>& c = coord_[i];
00176         int bit=0;
00177         for (int k=0; k<cdim; k++) // run over all directions in world
00178           {
00179                 if (k==missing)
00180                   {
00181                         c[k] = midpoint[k];
00182                         continue;
00183                   }
00184                 //k is not the missing direction
00185                 if (i&(1<<bit))  // check whether bit is set or not
00186                   c[k] = midpoint[k]+0.5*extension[k]; // bit is 1 in i
00187                 else
00188                   c[k] = midpoint[k]-0.5*extension[k]; // bit is 0 in i
00189                 bit++; // we have processed a direction
00190           }
00191 
00192         return c;
00193   }
00194 
00196   FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
00197   {
00198       FieldVector<ctype, cdim> g;
00199       int bit=0;
00200       for (int k=0; k<cdim; k++)
00201           if (k==missing)
00202               g[k] = midpoint[k];
00203           else
00204           {
00205               g[k] = midpoint[k] + (local[bit]-0.5)*extension[k];
00206               bit++;
00207           }
00208       return g;
00209   }
00210 
00212   FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& global) const
00213   {
00214       FieldVector<ctype, mydim> l; // result
00215       int bit=0;
00216       for (int k=0; k<cdim; k++)
00217           if (k!=missing)
00218           {
00219               l[bit] = (global[k]-midpoint[k])/extension[k] + 0.5;
00220               bit++;
00221           }
00222       return l;
00223   }
00224 
00226   ctype volume () const
00227   {
00228     ctype volume=1.0;
00229     for (int k=0; k<cdim; k++)
00230       if (k!=missing) volume *= extension[k];
00231     return volume;
00232   }
00233 
00236   ctype integrationElement (const FieldVector<ctype, mydim>& local) const
00237   {
00238     return volume();
00239   }
00240 
00242   FieldMatrix<ctype,cdim,mydim>& jacobianInverseTransposed (const FieldVector<ctype, mydim>& local) const
00243   {
00244       Jinv = 0.0;
00245       int k=0;
00246       for (int i=0; i<cdim; ++i)
00247       {
00248           if (i != missing)
00249           {
00250               Jinv[i][k] = 1.0/extension[i]; // set diagonal element
00251               k++;
00252           }
00253       }
00254       return Jinv;
00255   }
00256 
00258   bool checkInside (const FieldVector<ctype, mydim>& local) const
00259   {
00260         for (int i=0; i<mydim; i++)
00261           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00262         return true;
00263   }
00264 
00266   YaspGeometry (const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, int& m)
00267         : midpoint(p), extension(h), missing(m)
00268   {
00269         if (cdim!=mydim+1)
00270           DUNE_THROW(GridError, "general YaspGeometry assumes cdim=mydim+1");
00271   }
00272 
00274   YaspGeometry (const YaspGeometry& other) 
00275         : midpoint(other.midpoint), 
00276           extension(other.extension), 
00277           missing(other.missing)
00278   {
00279   }
00280 
00282   void print (std::ostream& s) const
00283   {
00284         s << "YaspGeometry<"<<mydim<<","<<cdim<< "> ";
00285         s << "midpoint";
00286         for (int i=0; i<cdim; i++)
00287                   s << " " << midpoint[i];
00288                 s << " extension";
00289         for (int i=0; i<cdim; i++)
00290                   s << " " << extension[i];
00291                 s << " missing is " << missing;
00292   }
00293 private:
00294   // the element is fully defined by its midpoint the extension
00295   // in each direction and the missing direction.
00296   // References are used because this information
00297   // is known outside the element in many cases.
00298   // Note cdim==cdim+1
00299 
00300   // IMPORTANT midpoint and extension are references,
00301   // YaspGeometry can't be copied
00302   const FieldVector<ctype, cdim> & midpoint;    // the midpoint
00303   const FieldVector<ctype, cdim> & extension;   // the extension
00304   int& missing;                         // the missing, i.e. constant direction
00305 
00306   // In addition we need memory in order to return references.
00307   // Possibly we should change this in the interface ...
00308   mutable FieldMatrix<ctype, cdim, mydim> Jinv;   // the jacobian inverse
00309   mutable FieldMatrix<ctype,  Power_m_p<2,mydim>::power, cdim> coord_;   // the coordinates 
00310 
00311   const YaspGeometry<mydim,cdim,GridImp>&
00312   operator = (const YaspGeometry<mydim,cdim,GridImp>& g);
00313   
00314 };
00315 
00316 
00317 
00319 template<int mydim, class GridImp>
00320 class YaspGeometry<mydim,mydim,GridImp> : public GeometryDefaultImplementation<mydim,mydim,GridImp,YaspGeometry>
00321 {
00322 public:
00323   typedef typename GridImp::ctype ctype;
00325   typedef Geometry<mydim,mydim,GridImp,Dune::YaspGeometry> ReferenceGeometry;
00326   
00328   GeometryType type () const
00329   {
00330       return GeometryType(GeometryType::cube,mydim);
00331   }
00332 
00334   int corners () const
00335   {
00336         return 1<<mydim;
00337   }
00338 
00340   const FieldVector<ctype, mydim>& operator[] (int i) const
00341   {
00342         assert( i >= 0 && i < (int) coord_.N() ); 
00343         FieldVector<ctype, mydim>& c = coord_[i];
00344         for (int k=0; k<mydim; k++)
00345           if (i&(1<<k))
00346                 c[k] = midpoint[k]+0.5*extension[k]; // kth bit is 1 in i
00347           else
00348                 c[k] = midpoint[k]-0.5*extension[k]; // kth bit is 0 in i
00349         return c;
00350   }
00351 
00353   FieldVector<ctype, mydim> global (const FieldVector<ctype, mydim>& local) const
00354   {
00355       FieldVector<ctype,mydim> g;
00356         for (int k=0; k<mydim; k++)
00357           g[k] = midpoint[k] + (local[k]-0.5)*extension[k];
00358         return g;
00359   }
00360 
00362   FieldVector<ctype, mydim> local (const FieldVector<ctype,mydim>& global) const
00363   {
00364       FieldVector<ctype, mydim> l; // result
00365         for (int k=0; k<mydim; k++)
00366           l[k] = (global[k]-midpoint[k])/extension[k] + 0.5;
00367         return l;
00368   }
00369 
00372   ctype integrationElement (const FieldVector<ctype, mydim>& local) const
00373   {
00374     return volume();
00375   }
00376 
00378   ctype volume () const
00379   {
00380     ctype vol=1.0;
00381     for (int k=0; k<mydim; k++) vol *= extension[k];
00382     return vol;
00383   }
00384 
00386   FieldMatrix<ctype,mydim,mydim>& jacobianInverseTransposed (const FieldVector<ctype, mydim>& local) const
00387   {
00388         for (int i=0; i<mydim; ++i)
00389           {
00390                 Jinv[i] = 0.0;                // set column to zero
00391                 Jinv[i][i] = 1.0/extension[i]; // set diagonal element
00392           }
00393         return Jinv;
00394   }
00395 
00397   bool checkInside (const FieldVector<ctype, mydim>& local) const
00398   {
00399         for (int i=0; i<mydim; i++)
00400           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00401         return true;
00402   }
00403 
00405   YaspGeometry (const FieldVector<ctype, mydim>& p, const FieldVector<ctype, mydim>& h)
00406         : midpoint(p), extension(h)
00407   {}
00408 
00410   YaspGeometry (const YaspGeometry& other) 
00411         : midpoint(other.midpoint), 
00412           extension(other.extension) 
00413   {
00414   }
00415 
00417   void print (std::ostream& s) const
00418   {
00419         s << "YaspGeometry<"<<mydim<<","<<mydim<< "> ";
00420         s << "midpoint";
00421         for (int i=0; i<mydim; i++)
00422                   s << " " << midpoint[i];
00423                 s << " extension";
00424         for (int i=0; i<mydim; i++)
00425                   s << " " << extension[i];
00426   }
00427 
00428 private:
00429   // the element is fully defined by midpoint and the extension
00430   // in each direction. References are used because this information
00431   // is known outside the element in many cases.
00432   // Note mydim==cdim
00433 
00434   // IMPORTANT midpoint and extension are references,
00435   // YaspGeometry can't be copied
00436   const FieldVector<ctype, mydim> & midpoint;   // the midpoint
00437   const FieldVector<ctype, mydim> & extension;  // the extension
00438 
00439   // In addition we need memory in order to return references.
00440   // Possibly we should change this in the interface ...
00441   mutable FieldMatrix<ctype, mydim, mydim> Jinv;   // the jacobian inverse
00442   mutable FieldMatrix<ctype, Power_m_p<2,mydim>::power, mydim> coord_;   // the coordinates 
00443 
00444   // disable copy
00445   const YaspGeometry<mydim,mydim,GridImp>&
00446   operator = (const YaspGeometry<mydim,mydim,GridImp>& g);
00447 };
00448 
00450 template<int cdim, class GridImp>
00451 class YaspGeometry<0,cdim,GridImp> : public GeometryDefaultImplementation<0,cdim,GridImp,YaspGeometry>
00452 {
00453 public:
00454   typedef typename GridImp::ctype ctype;
00455   
00457   GeometryType type () const
00458   {
00459       return GeometryType(GeometryType::cube,0);
00460   }
00461 
00463   int corners () const
00464   {
00465         return 1;
00466   }
00467 
00469   const FieldVector<ctype, cdim>& operator[] (int i) const
00470   {
00471         return position;
00472   }
00473 
00476   ctype integrationElement (const FieldVector<ctype, 0>& local) const
00477   {
00478         return 1.0;
00479   }
00480 
00482   FieldMatrix<ctype,cdim,0>& jacobianInverseTransposed (const FieldVector<ctype, 0>& local) const
00483   {
00484       static FieldMatrix<ctype,cdim,0> Jinv(0.0);
00485       return Jinv;
00486   }
00487 
00489   YaspGeometry (const FieldVector<ctype, cdim>& p) : position(p)
00490   {}
00491 
00493   void print (std::ostream& s) const
00494   {
00495         s << "YaspGeometry<"<<0<<","<<cdim<< "> ";
00496         s << "position " << position;
00497   }
00498 
00499 private:
00500   // IMPORTANT position is a reference,
00501   // YaspGeometry can't be copied
00502   const FieldVector<ctype, cdim> & position; 
00503 
00504   const YaspGeometry<0,cdim,GridImp>&
00505   operator = (const YaspGeometry<0,cdim,GridImp>& g);
00506 };
00507 
00508 // operator<< for all YaspGeometrys
00509 template <int mydim, int cdim, class GridImp>
00510 inline
00511 std::ostream& operator<< (std::ostream& s, YaspGeometry<mydim,cdim,GridImp>& e)
00512 {
00513   e.print(s);
00514   return s;
00515 }
00516 
00517 //========================================================================
00525 //========================================================================
00526 
00527 template<int codim, int dim, class GridImp>
00528 class YaspSpecialEntity :
00529   public GridImp::template Codim<codim>::Entity
00530 {
00531 public:
00532   typedef typename GridImp::ctype ctype;
00533   
00534   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00535   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00536 
00537   YaspSpecialEntity(const YGLI& g, const TSI& it) :
00538     GridImp::template Codim<codim>::Entity (YaspEntity<codim, dim, GridImp>(g,it))
00539     {};
00540   YaspSpecialEntity(const YaspEntity<codim, dim, GridImp>& e) :
00541     GridImp::template Codim<codim>::Entity (e)
00542     {};
00543   const TSI& transformingsubiterator () const
00544   {
00545         return this->realEntity.transformingsubiterator();
00546   }
00547   const YGLI& gridlevel () const
00548   {
00549         return this->realEntity.gridlevel();
00550   }
00551 };
00552 
00553 template<int codim, int dim, class GridImp>
00554 class YaspEntity 
00555 :  public EntityDefaultImplementation <codim,dim,GridImp,YaspEntity>
00556 {
00557 public:
00558   typedef typename GridImp::ctype ctype;
00559   
00560   typedef typename GridImp::template Codim<codim>::Geometry Geometry;
00562   int level () const
00563   {
00564         DUNE_THROW(GridError, "YaspEntity not implemented");
00565   }
00566 
00568   int index () const
00569   {
00570         DUNE_THROW(GridError, "YaspEntity not implemented");
00571   }
00572 
00574   const Geometry& geometry () const
00575   {
00576         DUNE_THROW(GridError, "YaspEntity not implemented");
00577   }
00578 
00580   PartitionType partitionType () const
00581   {
00582         DUNE_THROW(GridError, "YaspEntity not implemented");
00583   }
00584 
00585   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00586   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00587   YaspEntity (const YGLI& g, const TSI& it)
00588   {
00589         DUNE_THROW(GridError, "YaspEntity not implemented");
00590   }
00591 
00592   // IndexSets needs access to the private index methods
00593   friend class Dune::YaspLevelIndexSet<GridImp>;
00594   friend class Dune::YaspLeafIndexSet<GridImp>;
00595   friend class Dune::YaspGlobalIdSet<GridImp>;
00596   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00597 
00599   PersistentIndexType persistentIndex () const 
00600   { 
00601     DUNE_THROW(GridError, "YaspEntity not implemented");
00602   }
00603 
00605   int compressedIndex () const 
00606   {
00607     DUNE_THROW(GridError, "YaspEntity not implemented");
00608   }
00609 
00611   int compressedLeafIndex () const 
00612   {
00613     DUNE_THROW(GridError, "YaspEntity not implemented");
00614   }
00615 };
00616 
00617 
00618 // specialization for codim=0
00619 template<int dim, class GridImp>
00620 class YaspEntity<0,dim,GridImp> 
00621 : public EntityDefaultImplementation <0,dim,GridImp,YaspEntity>
00622 {
00623   enum { dimworld = GridImp::dimensionworld };
00624 public:
00625   typedef typename GridImp::ctype ctype;
00626   
00627   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00628   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00629 
00630   typedef YaspSpecialGeometry<dim-0,dim,GridImp> SpecialGeometry;
00631   
00632   typedef typename GridImp::template Codim<0>::Geometry Geometry;
00633   template <int cd>
00634   struct Codim
00635   {
00636     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
00637   };
00638   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
00639   typedef typename GridImp::template Codim<0>::LevelIntersectionIterator IntersectionIterator;
00640   typedef typename GridImp::template Codim<0>::LevelIntersectionIterator LevelIntersectionIterator;
00641   typedef typename GridImp::template Codim<0>::LeafIntersectionIterator  LeafIntersectionIterator;
00642   typedef typename GridImp::template Codim<0>::HierarchicIterator HierarchicIterator;
00643 
00645   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00646 
00648   typedef typename YGrid<dim,ctype>::iTupel iTupel;
00649 
00650   // constructor
00651   YaspEntity (const YGLI& g, const TSI& it)
00652     : _it(it), _g(g), _geometry(it.position(),it.meshsize())
00653   {
00654   }
00655 
00657   int level () const {return _g.level();}
00658 
00660   int index () const { return _it.superindex();} // superindex works also for iteration over subgrids
00661 
00663   int globalIndex () const {
00664     return _g.cell_global().index(_it.coord());
00665   }
00666 
00668   PartitionType partitionType () const
00669   {
00670         if (_g.cell_interior().inside(_it.coord())) return InteriorEntity;
00671         if (_g.cell_overlap().inside(_it.coord())) return OverlapEntity;
00672         return GhostEntity;
00673   }
00674 
00676   const Geometry& geometry () const { return _geometry; }
00677 
00680   template<int cc> int count () const
00681   {
00682     if (cc==dim) return 1<<dim;
00683     if (cc==1) return 2*dim;
00684     if (cc==dim-1) return dim*(1<<(dim-1));
00685     if (cc==0) return 1;
00686     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00687   }
00688 
00691   template<int cc>
00692   typename Codim<cc>::EntityPointer entity (int i) const
00693   {
00694       dune_static_assert( cc == dim || cc == 0 ,
00695           "YaspGrid only supports Entities with codim=dim and codim=0");
00696       // coordinates of the cell == coordinates of lower left corner
00697       if (cc==dim)
00698           {
00699                 iTupel coord = _it.coord();
00700 
00701                 // get corner from there
00702                 for (int k=0; k<dim; k++)
00703                   if (i&(1<<k)) (coord[k])++;
00704 
00705                 return YaspEntityPointer<cc,GridImp>(_g,_g.vertex_overlapfront().tsubbegin(coord));
00706           }
00707       if (cc==0)
00708           {
00709                 return YaspEntityPointer<cc,GridImp>(_g,_it);
00710           }
00711       DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00712   }
00713 
00715   EntityPointer father () const
00716   {
00717         // check if coarse level exists
00718         if (_g.level()<=0)
00719           DUNE_THROW(GridError, "tried to call father on level 0");
00720 
00721         // yes, get iterator to it
00722         YGLI cg = _g.coarser();
00723 
00724         // coordinates of the cell
00725         iTupel coord = _it.coord();
00726 
00727         // get coordinates on next coarser level
00728         for (int k=0; k<dim; k++) coord[k] = coord[k]/2;
00729 
00730         return YaspEntityPointer<0,GridImp>(cg,cg.cell_overlap().tsubbegin(coord));
00731   }
00732 
00742   const Geometry& geometryInFather () const
00743   {
00744     // determine which son we are
00745     int son = 0;
00746     for (int k=0; k<dim; k++)
00747       if (_it.coord(k)%2)
00748         son += (1<<k);
00749     
00750     // configure one of the 2^dim transformations
00751     return YaspFatherRelativeLocalElement<dim,GridImp>::getson(son);
00752   }
00753 
00754   const TSI& transformingsubiterator () const
00755   {
00756         return _it;
00757   }
00758 
00759   const YGLI& gridlevel () const
00760   {
00761         return _g;
00762   }
00763 
00764   bool isLeaf() const
00765   {
00766     return (_g.level() == _g.mg()->maxlevel());
00767   }
00768   
00770   IntersectionIterator ibegin () const
00771   {
00772     return YaspIntersectionIterator<GridImp>(*this,false);
00773   }
00774 
00776   LeafIntersectionIterator ileafbegin () const
00777   {
00778     // only if entity is leaf this iterator delivers intersections
00779     return YaspIntersectionIterator<GridImp>(*this, ! isLeaf() );
00780   }
00781 
00783   LevelIntersectionIterator ilevelbegin () const
00784   {
00785     return ibegin(); 
00786   }
00787 
00789   IntersectionIterator iend () const
00790   {
00791     return YaspIntersectionIterator<GridImp>(*this,true);
00792   }
00793 
00795   LeafIntersectionIterator ileafend () const
00796   {
00797     return iend();
00798   }
00799 
00801   LevelIntersectionIterator ilevelend () const
00802   {
00803     return iend();
00804   }
00805 
00810   HierarchicIterator hbegin (int maxlevel) const
00811   {
00812         return YaspHierarchicIterator<GridImp>(_g,_it,maxlevel);
00813   }
00814 
00816   HierarchicIterator hend (int maxlevel) const
00817   {
00818         return YaspHierarchicIterator<GridImp>(_g,_it,_g.level());
00819   }
00820 
00821 private:
00822   // IndexSets needs access to the private index methods
00823   friend class Dune::YaspLevelIndexSet<GridImp>;
00824   friend class Dune::YaspLeafIndexSet<GridImp>;
00825   friend class Dune::YaspGlobalIdSet<GridImp>;
00826 
00828   PersistentIndexType persistentIndex () const 
00829   { 
00830     // get size of global grid
00831     const iTupel& size =  _g.cell_global().size();
00832 
00833     // get coordinate correction for periodic boundaries
00834     int coord[dim];
00835     for (int i=0; i<dim; i++)
00836       {
00837         coord[i] = _it.coord(i);
00838         if (coord[i]<0) coord[i] += size[i];
00839         if (coord[i]>=size[i]) coord[i] -= size[i];
00840       }
00841 
00842     // encode codim
00843     PersistentIndexType id(0);
00844 
00845     // encode level
00846     id = id << yaspgrid_level_bits;
00847     id = id+PersistentIndexType(_g.level());
00848     
00849 
00850     // encode coordinates
00851     for (int i=dim-1; i>=0; i--)
00852       {
00853         id = id << yaspgrid_dim_bits;
00854         id = id+PersistentIndexType(coord[i]);
00855       }
00856     
00857     return id;
00858   }
00859 
00861   int compressedIndex () const 
00862   {
00863     return _it.superindex();
00864   }
00865 
00867   int compressedLeafIndex () const 
00868   {
00869     return _it.superindex();
00870   }
00871 
00873   template<int cc>
00874   PersistentIndexType subPersistentIndex (int i) const
00875   {
00876     if (cc==0)
00877       return persistentIndex();
00878 
00879     // get position of cell, note that global origin is zero
00880     // adjust for periodic boundaries
00881     int coord[dim];
00882     for (int k=0; k<dim; k++)
00883       {
00884         coord[k] = _it.coord(k);
00885         if (coord[k]<0) coord[k] += _g.cell_global().size(k);
00886         if (coord[k]>=_g.cell_global().size(k)) coord[k] -= _g.cell_global().size(k);
00887       }
00888     
00889     if (cc==dim)
00890       {
00891         // transform to vertex coordinates
00892         for (int k=0; k<dim; k++)
00893           if (i&(1<<k)) (coord[k])++;
00894 
00895         // determine min number of trailing zeroes
00896         int trailing = 1000;
00897         for (int i=0; i<dim; i++)
00898           {
00899             // count trailing zeros
00900             int zeros = 0;
00901             for (int j=0; j<_g.level(); j++)
00902               if (coord[i]&(1<<j))
00903                 break;
00904               else
00905                 zeros++;
00906             trailing = std::min(trailing,zeros);
00907           } 
00908 
00909         // determine the level of this vertex
00910         int level = _g.level()-trailing;
00911 
00912         // encode codim
00913         PersistentIndexType id(dim);
00914         
00915         // encode level
00916         id = id << yaspgrid_level_bits;
00917         id = id+PersistentIndexType(level);
00918         
00919         // encode coordinates
00920         for (int i=dim-1; i>=0; i--)
00921           {
00922             id = id << yaspgrid_dim_bits;
00923             id = id+PersistentIndexType(coord[i]>>trailing);
00924           }
00925         
00926         return id;
00927       }
00928 
00929     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
00930       {
00931         // Idea: Use the doubled grid to assign coordinates to faces
00932 
00933         // ivar is the direction that varies
00934         int ivar=i/2; 
00935 
00936         // compute position from cell position
00937         for (int k=0; k<dim; k++)
00938           coord[k] = coord[k]*2 + 1; // the doubled grid
00939         if (i%2) 
00940           coord[ivar] += 1;
00941         else
00942           coord[ivar] -= 1;
00943 
00944         // encode codim
00945         PersistentIndexType id(1);
00946         
00947         // encode level
00948         id = id << yaspgrid_level_bits;
00949         id = id+PersistentIndexType(_g.level());
00950         
00951         // encode coordinates
00952         for (int i=dim-1; i>=0; i--)
00953           {
00954             id = id << yaspgrid_dim_bits;
00955             id = id+PersistentIndexType(coord[i]);
00956           }
00957         
00958         return id;
00959       }
00960 
00961     if (cc==dim-1) // edges, exist only for dim>2
00962       {
00963         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
00964 
00965         // number of entities per direction
00966         int m=1<<(dim-1);
00967 
00968         // ifix is the direction that is fixed
00969         int ifix=(dim-1)-(i/m); 
00970 
00971         // compute position from cell position
00972         int bit=1;
00973         for (int k=0; k<dim; k++)
00974           {
00975             coord[k] = coord[k]*2+1; // cell position in doubled grid
00976             if (k==ifix) continue;
00977             if ((i%m)&bit) coord[k] += 1; else coord[k] -= 1;
00978             bit *= 2;
00979           } 
00980 
00981         // encode codim
00982         PersistentIndexType id(dim-1);
00983         
00984         // encode level
00985         id = id << yaspgrid_level_bits;
00986         id = id+PersistentIndexType(_g.level());
00987         
00988         // encode coordinates
00989         for (int i=dim-1; i>=0; i--)
00990           {
00991             id = id << yaspgrid_dim_bits;
00992             id = id+PersistentIndexType(coord[i]);
00993           }
00994         
00995         return id;
00996       }
00997 
00998     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00999   }
01000 
01002   template<int cc>
01003   int subCompressedIndex (int i) const
01004   {
01005     if (cc==0)
01006       return compressedIndex();
01007     
01008     // get cell position relative to origin of local cell grid
01009     iTupel coord;
01010     for (int k=0; k<dim; ++k) 
01011       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
01012 
01013     if (cc==dim) // vertices
01014       {
01015         // transform cell coordinate to corner coordinate
01016         for (int k=0; k<dim; k++)
01017           if (i&(1<<k)) (coord[k])++;
01018 
01019         // do lexicographic numbering
01020         int index = coord[dim-1];
01021         for (int k=dim-2; k>=0; --k)
01022           index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
01023         return index;
01024       }
01025 
01026     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01027       {
01028         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
01029 
01030         // ivar is the direction that varies
01031         int ivar=i/2; 
01032 
01033         // compute position from cell position
01034         if (i%2) coord[ivar] += 1; 
01035 
01036         // do lexicographic numbering
01037         int index = coord[dim-1];
01038         for (int k=dim-2; k>=0; --k)
01039           if (k==ivar)
01040             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01041           else
01042             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01043 
01044         // add size of all subsets for smaller directions
01045         for (int j=0; j<ivar; j++)
01046           {
01047             int n=_g.cell_overlap().size(j)+1;
01048             for (int l=0; l<dim; l++)
01049               if (l!=j) n *= _g.cell_overlap().size(l);
01050             index += n;
01051           }
01052 
01053         return index;
01054       }
01055 
01056     if (cc==dim-1) // edges, exist only for dim>2
01057       {
01058         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01059 
01060         // number of entities per direction
01061         int m=1<<(dim-1);
01062 
01063         // ifix is the direction that is fixed
01064         int ifix=(dim-1)-(i/m); 
01065 
01066         // compute position from cell position
01067         int bit=1;
01068         for (int k=0; k<dim; k++)
01069           {
01070             if (k==ifix) continue;
01071             if ((i%m)&bit) coord[k] += 1;
01072             bit *= 2;
01073           } 
01074 
01075         // do lexicographic numbering
01076         int index = coord[dim-1];
01077         for (int k=dim-2; k>=0; --k)
01078           if (k!=ifix)
01079             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01080           else
01081             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01082 
01083         // add size of all subsets for smaller directions
01084         for (int j=dim-1; j>ifix; j--)
01085           {
01086             int n=_g.cell_overlap().size(j);
01087             for (int l=0; l<dim; l++)
01088               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01089             index += n;
01090           }
01091 
01092         return index;
01093       }
01094 
01095     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01096   }
01097 
01099   template<int cc>
01100   int subCompressedLeafIndex (int i) const
01101   {
01102     if (cc==0)
01103       return compressedIndex();
01104     
01105     // get cell position relative to origin of local cell grid
01106     iTupel coord;
01107     for (int k=0; k<dim; ++k) 
01108       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
01109 
01110     if (cc==dim) // vertices
01111       {
01112         // transform cell coordinate to corner coordinate
01113         for (int k=0; k<dim; k++)
01114           if (i&(1<<k)) (coord[k])++;
01115 
01116         // move coordinates up to maxlevel
01117         for (int k=0; k<dim; k++)
01118           coord[k] = coord[k]<<(_g.mg()->maxlevel()-_g.level());
01119 
01120         // do lexicographic numbering
01121         int index = coord[dim-1];
01122         for (int k=dim-2; k>=0; --k)
01123           index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01124         return index;
01125       }
01126 
01127     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01128       {
01129         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
01130 
01131         // ivar is the direction that varies
01132         int ivar=i/2; 
01133 
01134         // compute position from cell position
01135         if (i%2) coord[ivar] += 1; 
01136 
01137         // do lexicographic numbering
01138         int index = coord[dim-1];
01139         for (int k=dim-2; k>=0; --k)
01140           if (k==ivar)
01141             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01142           else
01143             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01144 
01145         // add size of all subsets for smaller directions
01146         for (int j=0; j<ivar; j++)
01147           {
01148             int n=_g.cell_overlap().size(j)+1;
01149             for (int l=0; l<dim; l++)
01150               if (l!=j) n *= _g.cell_overlap().size(l);
01151             index += n;
01152           }
01153 
01154         return index;
01155       }
01156 
01157     if (cc==dim-1) // edges, exist only for dim>2
01158       {
01159         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01160 
01161         // number of entities per direction
01162         int m=1<<(dim-1);
01163 
01164         // ifix is the direction that is fixed
01165         int ifix=(dim-1)-(i/m); 
01166 
01167         // compute position from cell position
01168         int bit=1;
01169         for (int k=0; k<dim; k++)
01170           {
01171             if (k==ifix) continue;
01172             if ((i%m)&bit) coord[k] += 1;
01173             bit *= 2;
01174           } 
01175 
01176         // do lexicographic numbering
01177         int index = coord[dim-1];
01178         for (int k=dim-2; k>=0; --k)
01179           if (k!=ifix)
01180             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01181           else
01182             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01183 
01184         // add size of all subsets for smaller directions
01185         for (int j=dim-1; j>ifix; j--)
01186           {
01187             int n=_g.cell_overlap().size(j);
01188             for (int l=0; l<dim; l++)
01189               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01190             index += n;
01191           }
01192 
01193         return index;
01194       }
01195 
01196     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01197   }
01198   
01199   const TSI& _it;           // position in the grid level
01200   const YGLI& _g;           // access to grid level
01201   SpecialGeometry _geometry; // the element geometry
01202 };
01203 
01204 
01205 // specialization for codim=dim (vertex)
01206 template<int dim, class GridImp>
01207 class YaspEntity<dim,dim,GridImp> 
01208 : public EntityDefaultImplementation <dim,dim,GridImp,YaspEntity>
01209 {
01210   enum { dimworld = GridImp::dimensionworld };
01211 public:
01212   typedef typename GridImp::ctype ctype;
01213   
01214   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01215   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01216 
01217   typedef YaspSpecialGeometry<dim-dim,dim,GridImp> SpecialGeometry;
01218   
01219   typedef typename GridImp::template Codim<dim>::Geometry Geometry;
01220   template <int cd>
01221   struct Codim
01222   {
01223     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
01224   };
01225   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01226 
01228   typedef typename GridImp::PersistentIndexType PersistentIndexType;
01229 
01231   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01232 
01233   // constructor
01234   YaspEntity (const YGLI& g, const TSI& it)
01235         : _it(it), _g(g), _geometry(it.position())
01236   {  }
01237 
01239   int level () const {return _g.level();}
01240 
01242   int index () const {return _it.superindex();}
01243 
01245   int globalIndex () const { return _g.cell_global().index(_it.coord()); }
01246 
01247 
01249   const Geometry& geometry () const { return _geometry; }
01250 
01252   PartitionType partitionType () const
01253   {
01254         if (_g.vertex_interior().inside(_it.coord())) return InteriorEntity;
01255         if (_g.vertex_interiorborder().inside(_it.coord())) return BorderEntity;
01256         if (_g.vertex_overlap().inside(_it.coord())) return OverlapEntity;
01257         if (_g.vertex_overlapfront().inside(_it.coord())) return FrontEntity;
01258         return GhostEntity;
01259   }
01260 
01261 private:
01262   // IndexSets needs access to the private index methods
01263   friend class Dune::YaspLevelIndexSet<GridImp>;
01264   friend class Dune::YaspLeafIndexSet<GridImp>;
01265   friend class Dune::YaspGlobalIdSet<GridImp>;
01266 
01268   PersistentIndexType persistentIndex () const 
01269   { 
01270     // get coordinate and size of global grid
01271     const iTupel& size =  _g.vertex_global().size();
01272     int coord[dim];
01273 
01274     // correction for periodic boundaries
01275     for (int i=0; i<dim; i++)
01276       {
01277         coord[i] = _it.coord(i);
01278         if (coord[i]<0) coord[i] += size[i];
01279         if (coord[i]>=size[i]) coord[i] -= size[i];
01280       }
01281 
01282     // determine min number of trailing zeroes
01283     int trailing = 1000;
01284     for (int i=0; i<dim; i++)
01285       {
01286         // count trailing zeros
01287         int zeros = 0;
01288         for (int j=0; j<_g.level(); j++)
01289           if (coord[i]&(1<<j))
01290             break;
01291           else
01292             zeros++;
01293         trailing = std::min(trailing,zeros);
01294       } 
01295 
01296     // determine the level of this vertex
01297     int level = _g.level()-trailing;
01298 
01299     // encode codim
01300     PersistentIndexType id(dim);
01301 
01302     // encode level
01303     id = id << yaspgrid_level_bits;
01304     id = id+PersistentIndexType(level);
01305     
01306     // encode coordinates
01307     for (int i=dim-1; i>=0; i--)
01308       {
01309         id = id << yaspgrid_dim_bits;
01310         id = id+PersistentIndexType(coord[i]>>trailing);
01311       }
01312     
01313     return id;
01314   }
01315 
01317   int compressedIndex () const { return _it.superindex();}
01318 
01320   int compressedLeafIndex () const 
01321   { 
01322     if (_g.level()==_g.mg()->maxlevel())
01323       return _it.superindex();
01324 
01325     // not on leaf level, interpolate to finest grid
01326     int coord[dim];
01327     for (int i=0; i<dim; i++) coord[i] = _it.coord(i)-(_g).vertex_overlap().origin(i);
01328 
01329     // move coordinates up to maxlevel (multiply by 2 for each level
01330     for (int k=0; k<dim; k++)
01331       coord[k] = coord[k]*(1<<(_g.mg()->maxlevel()-_g.level()));
01332 
01333     // do lexicographic numbering
01334     int index = coord[dim-1];
01335     for (int k=dim-2; k>=0; --k)
01336       index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01337     return index;
01338   }
01339 
01340 public:  
01341   const TSI& transformingsubiterator() const { return _it; }
01342   const YGLI& gridlevel() const { return _g; }
01343 protected:
01344   const TSI& _it;                // position in the grid level
01345   const YGLI& _g;                // access to grid level
01346   SpecialGeometry  _geometry;     // the element geometry
01347   // temporary object
01348   mutable FieldVector<ctype, dim> loc;   // always computed before being returned
01349 };
01350 
01351 
01352 //========================================================================
01357 //========================================================================
01358 
01359 template<class GridImp>
01360 class YaspIntersectionIterator
01361 {
01362   enum { dim=GridImp::dimension };
01363   enum { dimworld=GridImp::dimensionworld };  
01364   typedef typename GridImp::ctype ctype;
01365   YaspIntersectionIterator();
01366 public:
01367   // types used from grids
01368   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01369   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01370   typedef typename GridImp::template Codim<0>::Entity Entity;
01371   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01372   typedef typename GridImp::template Codim<1>::Geometry Geometry;
01373   typedef typename GridImp::template Codim<1>::LocalGeometry LocalGeometry;
01374   typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01375   typedef YaspSpecialGeometry<dim-1,dimworld,GridImp> SpecialGeometry;
01376   typedef YaspSpecialGeometry<dim-1,dim,GridImp> SpecialLocalGeometry;
01377   typedef Dune::Intersection<const GridImp, Dune::YaspIntersectionIterator> Intersection;
01378   
01380   void increment()
01381   {
01382         // check end
01383         if (_count==2*dim) return; // end iterator reached, we are done
01384         // update count, check end
01385         _count++;
01386 
01387         // update intersection iterator from current position
01388         if (_face==0) // direction remains valid
01389           {
01390                 _face = 1; // 0->1, _dir remains
01391 
01392                 // move transforming iterator
01393                 _outside.transformingsubiterator().move(_dir,2); // move two cells in positive direction
01394 
01395                 // make up faces
01396                 _pos_self_local[_dir] = 1.0;
01397                 _pos_nb_local[_dir] = 0.0;
01398                 _pos_world[_dir] += _inside.transformingsubiterator().meshsize(_dir);
01399 
01400                 // make up unit outer normal direction
01401                 _normal[_dir] = 1.0;
01402           }
01403         else // change direction
01404           {
01405                 // move transforming iterator
01406                 _outside.transformingsubiterator().move(_dir,-1); // move one cell back
01407                 if (_count==2*dim) return;
01408                 
01409                 // make up faces
01410                 _pos_self_local[_dir] = 0.5;
01411                 _pos_nb_local[_dir] = 0.5;
01412                 _pos_world[_dir] = _inside.transformingsubiterator().position(_dir);
01413 
01414                 // make up unit outer normal direction
01415                 _normal[_dir] = 0.0;
01416 
01417                 _face = 0;
01418                 _dir += 1;
01419 
01420                 // move transforming iterator
01421                 _outside.transformingsubiterator().move(_dir,-1); // move one cell in negative direction
01422 
01423                 // make up faces
01424                 _pos_self_local[_dir] = 0.0;
01425                 _pos_nb_local[_dir] = 1.0;
01426                 _pos_world[_dir] -= 0.5*_inside.transformingsubiterator().meshsize(_dir);
01427 
01428                 // make up unit outer normal direction
01429                 _normal[_dir] = -1.0;
01430           }
01431   }
01432 
01434   bool equals (const YaspIntersectionIterator& other) const
01435   {
01436       return _inside.equals(other._inside) && _outside.equals(other._outside);
01437   }
01438 
01440   const Intersection & dereference() const
01441   {
01442     return reinterpret_cast<const Intersection&>(*this);
01443   }
01444   
01449   bool boundary () const
01450   {
01451     // The transforming iterator can be safely moved beyond the boundary.
01452     // So we only have to compare against the cell_global grid
01453       return (_outside.transformingsubiterator().coord(_dir) < _inside.gridlevel().cell_global().min(_dir)
01454             ||
01455             _outside.transformingsubiterator().coord(_dir) > _inside.gridlevel().cell_global().max(_dir));
01456   }
01457 
01459   bool neighbor () const
01460   {
01461     return
01462       (_outside.transformingsubiterator().coord(_dir) >= _inside.gridlevel().cell_overlap().min(_dir)
01463        &&
01464        _outside.transformingsubiterator().coord(_dir) <= _inside.gridlevel().cell_overlap().max(_dir));
01465     for (int d = 0; d < dim; d++)
01466     {
01467       if (_outside.transformingsubiterator().coord(_dir)<_inside.gridlevel().cell_overlap().min(_dir)
01468           ||
01469           _outside.transformingsubiterator().coord(_dir) > _inside.gridlevel().cell_overlap().max(_dir))
01470         return false;
01471     }
01472     return true;
01473   }
01474 
01476   bool conforming () const
01477   {
01478     return true;
01479   }
01480 
01483   EntityPointer inside() const
01484     {
01485       return _inside;
01486     }
01487   
01490   EntityPointer outside() const
01491     {
01492         return _outside;
01493     }
01494 
01497   int boundaryId() const
01498     {
01499       if(boundary()) return numberInSelf()+1;
01500       return 0;
01501     }
01502   
01504   FieldVector<ctype, dimworld> outerNormal (const FieldVector<ctype, dim-1>& local) const
01505   {
01506         return _normal;
01507   }
01508 
01509 
01511   FieldVector<ctype, dimworld> unitOuterNormal (const FieldVector<ctype, dim-1>& local) const
01512   {
01513         return _normal;
01514   }
01515 
01519   FieldVector<ctype, dimworld> integrationOuterNormal (const FieldVector<ctype, dim-1>& local) const
01520   {
01521         FieldVector<ctype, dimworld> n = _normal;
01522         n *= intersectionGlobal().volume();
01523         return n;
01524   }
01525  
01529   const LocalGeometry& intersectionSelfLocal () const
01530   {
01531         return _is_self_local;
01532   }
01533 
01537   const LocalGeometry& intersectionNeighborLocal () const
01538   {
01539         return _is_nb_local;
01540   }
01541 
01545   const Geometry& intersectionGlobal () const
01546   {
01547         return _is_global;
01548   }
01549 
01551   GeometryType type () const
01552   {
01553     return intersectionSelfLocal().type();
01554   }
01555 
01557   int numberInSelf () const
01558   {
01559         return _count;
01560   }
01561 
01563   int numberInNeighbor () const
01564   {
01565         return _count + 1-2*_face;
01566   }
01567 
01569   YaspIntersectionIterator (const YaspEntity<0,dim,GridImp>& myself, bool toend)
01570       :   _inside(myself.gridlevel(),
01571                   myself.transformingsubiterator()),
01572           _outside(myself.gridlevel(),
01573                    myself.transformingsubiterator()),
01574           // initialize to first neighbor
01575           _count(0),
01576           _dir(0),
01577           _face(0),
01578           _pos_self_local(0.5),
01579           _pos_nb_local(0.5),
01580           _pos_world(myself.transformingsubiterator().position()),
01581           _ext_local(1.0),
01582           _normal(0.0),
01583           _is_self_local(_pos_self_local,_ext_local,_dir),
01584           _is_nb_local(_pos_nb_local,_ext_local,_dir),
01585           _is_global(_pos_world,_inside.transformingsubiterator().meshsize(),_dir)
01586   {
01587         // making an end iterator?
01588         if (toend)
01589           {
01590                 // initialize end iterator
01591                 _count = 2*dim;
01592                 return;
01593           }
01594 
01595         // move transforming iterator
01596         _outside.transformingsubiterator().move(_dir,-1);
01597 
01598         // make up faces
01599         _pos_self_local[0] = 0.0;
01600         _pos_nb_local[0] = 1.0;
01601         _pos_world[0] -= 0.5*_inside.transformingsubiterator().meshsize(0);
01602 
01603         // make up unit outer normal direction
01604         _normal[0] = -1.0;
01605   }
01606 
01608   YaspIntersectionIterator (const YaspIntersectionIterator& it)
01609       :   _inside(it._inside),
01610           _outside(it._outside),
01611           _count(it._count),
01612           _dir(it._dir),
01613           _face(it._face),
01614           _pos_self_local(it._pos_self_local),
01615           _pos_nb_local(it._pos_nb_local),
01616           _pos_world(it._pos_world),
01617           _ext_local(it._ext_local),
01618           _normal(it._normal),
01619           // Important: _is_* must be recreated -- not copied!
01620           _is_self_local(_pos_self_local,_ext_local,_dir),
01621           _is_nb_local(_pos_nb_local,_ext_local,_dir),
01622           _is_global(_pos_world,_inside.transformingsubiterator().meshsize(),_dir)
01623   {}
01624 
01626   YaspIntersectionIterator & operator = (const YaspIntersectionIterator& it)
01627     {
01628       /* Assert same Iterator Context */
01629       //_inside = it._inside;
01630       if (! _inside.equals(it._inside))
01631         DUNE_THROW(GridError, "assignment of YaspIntersectionIterator "
01632                    << "with different inside Entity");
01633       
01634       /* Assign EntityPointers */
01635       _outside = it._outside;
01636       
01637       /* Assign current position */
01638       _count = it._count;
01639       _dir = it._dir;
01640       _face = it._face;
01641       _pos_self_local = it._pos_self_local;
01642       _pos_nb_local = it._pos_nb_local;
01643       _pos_world = it._pos_world;
01644       _ext_local = it._ext_local;
01645       _normal = it._normal;
01646 
01647       return *this;
01648   }
01649 
01650 private:
01651   /* EntityPointers (get automatically updated) */
01652   const YaspEntityPointer<0,GridImp> _inside;  
01653 //  const YaspEntity<0,dim,GridImp>&   _myself;  //!< reference to myself
01654   YaspEntityPointer<0,GridImp> _outside;       
01655   /* current position */
01656   int _count;                             
01657   int _dir;                               
01658   int _face;                              
01659   /* precalculated data for current position */
01660   FieldVector<ctype, dim> _pos_self_local;
01661   FieldVector<ctype, dim> _pos_nb_local;  
01662   FieldVector<ctype, dimworld>_pos_world;     
01663   FieldVector<ctype, dim> _ext_local;     
01664   /* normal vector */
01665   FieldVector<ctype, dimworld> _normal;   
01666   /* geometry objects (get automatically updated) */
01667   SpecialLocalGeometry    _is_self_local; 
01668   SpecialLocalGeometry    _is_nb_local;   
01669   SpecialGeometry         _is_global;     
01670 };
01671 
01672 
01673 //========================================================================
01677 //========================================================================
01678 
01679 template<class GridImp>
01680 class YaspHierarchicIterator :
01681   public YaspEntityPointer<0,GridImp>
01682 {
01683   enum { dim=GridImp::dimension };
01684   enum { dimworld=GridImp::dimensionworld };  
01685   typedef typename GridImp::ctype ctype;
01686 public:
01687   // types used from grids
01688   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01689   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01690   typedef typename GridImp::template Codim<0>::Entity Entity;
01691   typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01692 
01694   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01695 
01697   YaspHierarchicIterator (const YGLI& g, const TSI& it, int maxlevel) :
01698       YaspEntityPointer<0,GridImp>(g,it)
01699   {
01700       // now iterator points to current cell
01701       StackElem se(this->_g);
01702       se.coord = this->_it.coord();
01703       stack.push(se);
01704 
01705       // determine maximum level
01706       _maxlevel = std::min(maxlevel,this->_g.mg()->maxlevel());
01707 
01708       // if maxlevel not reached then push yourself and sons
01709       if (this->_g.level()<_maxlevel)
01710       {
01711           push_sons();
01712       }
01713 
01714       // and make iterator point to first son if stack is not empty
01715       if (!stack.empty())
01716           pop_tos();
01717   }
01718 
01720   YaspHierarchicIterator (const YaspHierarchicIterator& it) :
01721     YaspEntityPointer<0,GridImp>(it._g, it._it),
01722     _maxlevel(it._maxlevel), stack(it.stack)
01723   {}
01724 
01726   void increment ()
01727   {
01728         // sanity check: do nothing when stack is empty
01729         if (stack.empty()) return;
01730 
01731         // if maxlevel not reached then push sons
01732         if (this->_g.level()<_maxlevel)
01733           push_sons();
01734 
01735         // in any case pop one element
01736         pop_tos();
01737   }
01738 
01739   void print (std::ostream& s) const
01740   {
01741         s << "HIER: " << "level=" << this->_g.level()
01742           << " position=" << this->_it.coord()
01743           << " superindex=" << this->_it.superindex()
01744           << " maxlevel=" << this->_maxlevel
01745           << " stacksize=" << stack.size()
01746           << std::endl;
01747   }
01748 
01749 private:
01750   int _maxlevel;         
01751 
01752   struct StackElem {
01753         YGLI g;       // grid level of the element
01754         iTupel coord; // and the coordinates
01755         StackElem(YGLI gg) : g(gg) {}
01756   };
01757   std::stack<StackElem> stack;      
01758 
01759   // push sons of current element on the stack
01760   void push_sons ()
01761   {
01762         // yes, process all 1<<dim sons
01763         StackElem se(this->_g.finer());
01764         for (int i=0; i<(1<<dim); i++)
01765           {
01766                 for (int k=0; k<dim; k++)
01767                   if (i&(1<<k))
01768                         se.coord[k] = this->_it.coord(k)*2+1;
01769                   else
01770                         se.coord[k] = this->_it.coord(k)*2;
01771                 stack.push(se);
01772           }
01773   }
01774 
01775   // make TOS the current element
01776   void pop_tos ()
01777   {
01778         StackElem se = stack.top();
01779         stack.pop();
01780         this->_g = se.g;
01781         this->_it.reinit(this->_g.cell_overlap(),se.coord);
01782   }
01783 };
01784 
01785 //========================================================================
01795 //========================================================================
01796 template<int codim, class GridImp>
01797 class YaspEntityPointer
01798 {
01800   enum { dim=GridImp::dimension };
01802   enum { dimworld=GridImp::dimensionworld };
01803   typedef typename GridImp::ctype ctype;
01804 public:
01805   typedef typename GridImp::template Codim<codim>::Entity Entity;
01806   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01807   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01808   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
01809   typedef YaspEntityPointer<codim,GridImp> EntityPointerImp;
01810 protected:  
01811   typedef YaspEntity<codim, dim, GridImp> YaspEntityImp; 
01812 
01813 public:
01815   enum { codimension = codim };
01816   
01818   YaspEntityPointer (const YGLI & g, const TSI & it) : _g(g), _it(it), _entity(_g,_it)
01819   {
01820         if (codim>0 && codim<dim)
01821           {
01822                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
01823           }
01824   }
01825 
01827   YaspEntityPointer (const YaspEntityImp& entity) 
01828     : _g(entity.gridlevel()), 
01829       _it(entity.transformingsubiterator()), 
01830       _entity(entity)
01831   {
01832         if (codim>0 && codim<dim)
01833           {
01834                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
01835           }
01836   }
01837 
01839   YaspEntityPointer (const YaspEntityPointer& rhs) : _g(rhs._g), _it(rhs._it), _entity(_g,_it)
01840   {
01841         if (codim>0 && codim<dim)
01842           {
01843                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
01844           }
01845   }
01846 
01848   bool equals (const YaspEntityPointer& rhs) const
01849   {
01850         return (_it==rhs._it && _g == rhs._g);
01851   }
01852 
01854   Entity& dereference() const
01855   {
01856         return _entity;
01857   }
01858 
01860   void compactify() {}
01861 
01863   int level () const {return _g.level();}
01864 
01865   const YaspEntityPointer&
01866   operator = (const YaspEntityPointer& rhs)
01867     {
01868       _g = rhs._g;
01869       _it = rhs._it;
01870       /* _entity = i._entity
01871        * is done implicitely, as the entity is completely
01872        * defined via the interator it belongs to
01873        */
01874       return *this;
01875     }  
01876   
01877   const TSI& transformingsubiterator () const
01878   {
01879         return _it;
01880   }
01881 
01882   const YGLI& gridlevel () const
01883   {
01884         return _g;
01885   }
01886 
01887   TSI& transformingsubiterator ()
01888   {
01889         return _it;
01890   }
01891 
01892   YGLI& gridlevel ()
01893   {
01894         return _g;
01895   }
01896 
01897 protected:
01898   YGLI _g;               // access to grid level
01899   TSI _it;               // position in the grid level
01900   mutable SpecialEntity _entity; 
01901 };
01902 
01903 //========================================================================
01911 //========================================================================
01912 
01913 template<int codim, PartitionIteratorType pitype, class GridImp>
01914 class YaspLevelIterator :
01915   public YaspEntityPointer<codim,GridImp>
01916 {
01918   enum { dim=GridImp::dimension };
01920   enum { dimworld=GridImp::dimensionworld };
01921   typedef typename GridImp::ctype ctype;
01922 public:
01923   typedef typename GridImp::template Codim<codim>::Entity Entity;
01924   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01925   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01926   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
01927 
01929   YaspLevelIterator (const YGLI & g, const TSI & it) :
01930     YaspEntityPointer<codim,GridImp>(g,it) {}
01931 
01933   YaspLevelIterator (const YaspLevelIterator& i) :
01934     YaspEntityPointer<codim,GridImp>(i) {}
01935 
01937   void increment()
01938   {
01939         ++(this->_it);
01940   }
01941 };
01942 
01943 
01944 //========================================================================
01949 //========================================================================
01950 
01951 template <class GridImp> 
01952 struct YaspLevelIndexSetTypes
01953 {
01955   template<int cd>
01956   struct Codim
01957   {
01958     template<PartitionIteratorType pitype>
01959     struct Partition
01960     {
01961       typedef typename GridImp::Traits::template Codim<cd>::template Partition<pitype>::LevelIterator Iterator;
01962     };
01963   };
01964 };
01965 
01966 template<class GridImp>
01967 class YaspLevelIndexSet : public IndexSetDefaultImplementation<GridImp,YaspLevelIndexSet<GridImp>,YaspLevelIndexSetTypes<GridImp> >
01968 {
01969   typedef IndexSetDefaultImplementation<GridImp,YaspLevelIndexSet<GridImp>,YaspLevelIndexSetTypes<GridImp> > Base;
01970 public:
01971 
01973   YaspLevelIndexSet (const GridImp& g, int l) : grid(g), level(l)
01974   {
01975     // contains a single element type;
01976     for (int codim=0; codim<=GridImp::dimension; codim++)
01977       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
01978   }
01979 
01981   template<int cd>
01982   int index (const typename GridImp::Traits::template Codim<cd>::Entity& e) const 
01983   {
01984     return grid.template getRealEntity<cd>(e).compressedIndex(); 
01985   }
01986 
01988   template<int cc>
01989   int subIndex (const typename GridImp::Traits::template Codim<0>::Entity& e, int i) const
01990   {
01991     return grid.template getRealEntity<0>(e).template subCompressedIndex<cc>(i);
01992   }
01993 
01995   int size (GeometryType type) const
01996   {
01997     return grid.size( level, type );
01998   }
01999 
02001   int size (int codim) const
02002   {
02003     return grid.size( level, codim );
02004   }
02005   
02007   template<class EntityType>
02008   bool contains (const EntityType& e) const
02009   {
02010     return e.level() == level;
02011   }
02012   
02014   const std::vector<GeometryType>& geomTypes (int codim) const
02015   {
02016     return mytypes[codim];
02017   }
02018 
02020   template<int cd, PartitionIteratorType pitype>
02021   typename Base::template Codim<cd>::template Partition<pitype>::Iterator begin () const
02022   {
02023     return grid.template lbegin<cd,pitype>(level);
02024   }
02025   
02027   template<int cd, PartitionIteratorType pitype>
02028   typename Base::template Codim<cd>::template Partition<pitype>::Iterator end () const
02029   {
02030     return grid.template lend<cd,pitype>(level);
02031   }
02032 
02033 private:
02034   const GridImp& grid;
02035   int level;
02036   std::vector<GeometryType> mytypes[GridImp::dimension+1];
02037 };
02038 
02039 
02040 // Leaf Index Set
02041 
02042 template <class GridImp> 
02043 struct YaspLeafIndexSetTypes
02044 {
02046   template<int cd>
02047   struct Codim
02048   {
02049     template<PartitionIteratorType pitype>
02050     struct Partition
02051     {
02052       typedef typename GridImp::Traits::template Codim<cd>::template Partition<pitype>::LevelIterator Iterator;
02053     };
02054   };
02055 };
02056 
02057 template<class GridImp>
02058 class YaspLeafIndexSet : public IndexSetDefaultImplementation<GridImp,YaspLeafIndexSet<GridImp>,YaspLeafIndexSetTypes<GridImp> >
02059 {
02060   typedef IndexSetDefaultImplementation<GridImp,YaspLeafIndexSet<GridImp>,YaspLeafIndexSetTypes<GridImp> > Base;
02061 public:
02062 
02064   YaspLeafIndexSet (const GridImp& g) : grid(g)
02065   {
02066     // contains a single element type;
02067     for (int codim=0; codim<=GridImp::dimension; codim++)
02068       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
02069   }
02070 
02072   /*
02073     We use the remove_const to extract the Type from the mutable class,
02074     because the const class is not instantiated yet.
02075   */
02076   template<int cd>
02077   int index (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const 
02078   {
02079     assert(cd==0 || cd==GridImp::dimension);
02080     return grid.template getRealEntity<cd>(e).compressedLeafIndex(); 
02081   }
02082 
02084   /*
02085     We use the remove_const to extract the Type from the mutable class,
02086     because the const class is not instantiated yet.
02087   */
02088   template<int cc>
02089   int subIndex (const typename remove_const<GridImp>::type::Traits::template Codim<0>::Entity& e, int i) const
02090   {
02091     return grid.template getRealEntity<0>(e).template subCompressedLeafIndex<cc>(i);
02092   }
02093 
02095   int size (GeometryType type) const
02096   {
02097     return grid.size( grid.maxLevel(), type );
02098   }
02099 
02101   int size (int codim) const
02102   {
02103     return grid.size( grid.maxLevel(), codim );
02104   }
02105 
02107   template<class EntityType>
02108   bool contains (const EntityType& e) const
02109   {
02110     //return e.isLeaf();
02111     return (e.level() == grid.maxLevel());
02112   }
02113   
02115   const std::vector<GeometryType>& geomTypes (int codim) const
02116   {
02117     return mytypes[codim];
02118   }
02119 
02121   template<int cd, PartitionIteratorType pitype>
02122   typename Base::template Codim<cd>::template Partition<pitype>::Iterator begin () const
02123   {
02124     return grid.template lbegin<cd,pitype>(grid.maxLevel());
02125   }
02126   
02128   template<int cd, PartitionIteratorType pitype>
02129   typename Base::template Codim<cd>::template Partition<pitype>::Iterator end () const
02130   {
02131     return grid.template lend<cd,pitype>(grid.maxLevel());
02132   }
02133 
02134 private:
02135   const GridImp& grid;
02136   enum { ncodim = remove_const<GridImp>::type::dimension+1 };
02137   std::vector<GeometryType> mytypes[ncodim];
02138 };
02139 
02140 
02141 
02142 
02143 //========================================================================
02148 //========================================================================
02149 
02150 template<class GridImp>
02151 class YaspGlobalIdSet : public IdSetDefaultImplementation<GridImp,YaspGlobalIdSet<GridImp>,
02152                                      typename remove_const<GridImp>::type::PersistentIndexType >
02153 /*
02154   We used the remove_const to extract the Type from the mutable class,
02155   because the const class is not instantiated yet.
02156 */
02157 {
02158 public:
02160   typedef typename remove_const<GridImp>::type::PersistentIndexType IdType;
02161 
02163   YaspGlobalIdSet (const GridImp& g) : grid(g) {}
02164 
02166   /*
02167     We use the remove_const to extract the Type from the mutable class,
02168     because the const class is not instantiated yet.
02169   */
02170   template<int cd>
02171   IdType id (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const 
02172   {
02173     return grid.template getRealEntity<cd>(e).persistentIndex();
02174   }
02175 
02177   /*
02178     We use the remove_const to extract the Type from the mutable class,
02179     because the const class is not instantiated yet.
02180   */
02181   template<int cc>
02182   IdType subId (const typename remove_const<GridImp>::type::Traits::template Codim<0>::Entity& e, int i) const
02183   {
02184     return grid.template getRealEntity<0>(e).template subPersistentIndex<cc>(i);
02185   }
02186 
02187 private:
02188   const GridImp& grid;
02189 };
02190 
02191 
02192 template<int dim, int dimworld>
02193 struct YaspGridFamily
02194 {
02195 #if HAVE_MPI
02196   typedef CollectiveCommunication<MPI_Comm> CCType;
02197 #else
02198   typedef CollectiveCommunication<Dune::YaspGrid<dim,dimworld> > CCType;
02199 #endif
02200 
02201   typedef GridTraits<dim,dimworld,Dune::YaspGrid<dim,dimworld>,
02202                      YaspGeometry,YaspEntity,
02203                      YaspEntityPointer,YaspLevelIterator,
02204                      YaspIntersectionIterator, // leaf  intersection 
02205                      YaspIntersectionIterator, // level intersection 
02206                      YaspIntersectionIterator, // leaf  intersection iter 
02207                      YaspIntersectionIterator, // level intersection iter 
02208                      YaspHierarchicIterator,
02209                      YaspLevelIterator,
02210                      YaspLevelIndexSet<const YaspGrid<dim,dimworld> >,
02211                      YaspLevelIndexSetTypes<const YaspGrid<dim,dimworld> >,
02212                      YaspLeafIndexSet<const YaspGrid<dim,dimworld> >,
02213                      YaspLeafIndexSetTypes<const YaspGrid<dim,dimworld> >,
02214                      YaspGlobalIdSet<const YaspGrid<dim,dimworld> >, 
02215                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02216                      YaspGlobalIdSet<const YaspGrid<dim,dimworld> >, 
02217                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02218                      CCType> 
02219   Traits;  
02220 };
02221 
02222 template<int dim, int codim>
02223 struct YaspCommunicateMeta {
02224   template<class G, class DataHandle>
02225   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02226   {
02227     if (data.contains(dim,codim))
02228       g.template communicateCodim<DataHandle,codim>(data,iftype,dir,level);
02229     YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
02230   }
02231 };
02232 
02233 template<int dim>
02234 struct YaspCommunicateMeta<dim,0> {
02235   template<class G, class DataHandle>
02236   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02237   {
02238     if (data.contains(dim,0))
02239       g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
02240   }
02241 };
02242 
02243 
02244 //************************************************************************
02260 template<int dim, int deprecated_dimworld=dim>
02261 class YaspGrid :
02262   public GridDefaultImplementation<dim,dim,yaspgrid_ctype,YaspGridFamily<dim,dim> >,
02263   public MultiYGrid<dim,yaspgrid_ctype>
02264 {
02265   typedef const YaspGrid<dim,dim> GridImp;
02266 
02267   dune_static_assert(dim==deprecated_dimworld, "YaspGrid can only be instantiated with dim==dimworld");
02268 
02269   void init()
02270   {
02271     setsizes();
02272     indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dim> >(*this,0) );
02273     theleafindexset.push_back( new YaspLeafIndexSet<const YaspGrid<dim,dim> >(*this) );
02274     theglobalidset.push_back( new YaspGlobalIdSet<const YaspGrid<dim,dim> >(*this) );
02275   }
02276 public:
02278   typedef yaspgrid_ctype ctype;
02279 
02280   // define the persistent index type
02281   typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits> PersistentIndexType;
02282 
02284   typedef YaspGridFamily<dim,dim> GridFamily;
02285   // the Traits 
02286   typedef typename YaspGridFamily<dim,dim>::Traits Traits;
02287 
02288   // need for friend declarations in entity
02289   typedef YaspLevelIndexSet<YaspGrid<dim,dim> > LevelIndexSetType;
02290   typedef YaspLeafIndexSet<YaspGrid<dim,dim> > LeafIndexSetType;
02291   typedef YaspGlobalIdSet<YaspGrid<dim,dim> > GlobalIdSetType;
02292 
02294   enum { MAXL=64 };
02295 
02297   typedef MultiYGrid<dim,ctype> YMG;
02298   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
02299   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
02300   typedef typename MultiYGrid<dim,ctype>::Intersection IS;
02301   typedef typename std::deque<IS>::const_iterator ISIT;
02302 
02304   std::string name() const { return "YaspGrid"; };
02305 
02313 #if HAVE_MPI
02314   YaspGrid (MPI_Comm comm, Dune::FieldVector<ctype, dim> L, 
02315             Dune::FieldVector<int, dim> s, 
02316             Dune::FieldVector<bool, dim> periodic, int overlap) 
02317     : MultiYGrid<dim,ctype>(comm,L,s,periodic,overlap), ccobj(comm)
02318   {
02319     init();
02320   }
02321 #endif
02322   
02334   YaspGrid (Dune::FieldVector<ctype, dim> L, 
02335             Dune::FieldVector<int, dim> s, 
02336             Dune::FieldVector<bool, dim> periodic, int overlap) 
02337 #if HAVE_MPI
02338     : MultiYGrid<dim,ctype>(MPI_COMM_SELF,L,s,periodic,overlap), ccobj(MPI_COMM_SELF)
02339 #else
02340     : MultiYGrid<dim,ctype>(L,s,periodic,overlap)
02341 #endif
02342   {
02343     init();
02344   }
02345 
02346   ~YaspGrid()
02347   {
02348     deallocatePointers(indexsets);
02349     deallocatePointers(theleafindexset);
02350     deallocatePointers(theglobalidset);
02351   }
02352 
02353   private:
02354   // do not copy this class 
02355   YaspGrid(const YaspGrid&);
02356 
02357   public:
02358   
02362   int maxLevel() const {return MultiYGrid<dim,ctype>::maxlevel();} // delegate
02363 
02365   void globalRefine (int refCount)
02366   {
02367     if (refCount < -maxLevel())
02368       DUNE_THROW(GridError, "Only " << maxLevel() << " levels left. " <<
02369                  "Coarsening " << -refCount << " levels requested!");
02370     bool b=true;
02371     for (int k=refCount; k<0; k++)
02372       {
02373         MultiYGrid<dim,ctype>::coarsen();
02374         setsizes();
02375         indexsets.pop_back();
02376       }
02377     for (int k=0; k<refCount; k++)
02378       {
02379         MultiYGrid<dim,ctype>::refine(b);
02380         setsizes();
02381         indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dim> >(*this,maxLevel()) );
02382       }
02383   }
02384 
02385   // set options for refinement
02386   void refineOptions (bool b)
02387   {
02388     keep_ovlp=b;
02389   }
02390   
02392   void refine (bool b)
02393   {
02394     MultiYGrid<dim,ctype>::refine(b);
02395     setsizes();
02396     indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dim> >(*this,maxLevel()) );
02397   }
02398 
02400   bool adapt ()    
02401   {
02402     globalRefine(1);
02403     return true; 
02404   }
02405   
02407   template<int cd, PartitionIteratorType pitype>
02408   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
02409     {
02410       return levelbegin<cd,pitype>(level);
02411     }
02412 
02414   template<int cd, PartitionIteratorType pitype>
02415   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
02416     {
02417       return levelend<cd,pitype>(level);
02418     }
02419 
02421   template<int cd>
02422   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
02423     {
02424       return levelbegin<cd,All_Partition>(level);
02425     }
02426 
02428   template<int cd>
02429   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
02430     {
02431       return levelend<cd,All_Partition>(level);
02432     }
02433 
02435   template<int cd, PartitionIteratorType pitype>
02436   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const
02437     {
02438       return levelbegin<cd,pitype>(maxLevel());
02439     }
02440 
02442   template<int cd, PartitionIteratorType pitype>
02443   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const
02444     {
02445       return levelend<cd,pitype>(maxLevel());
02446     }
02447 
02449   template<int cd>
02450   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
02451     {
02452       return levelbegin<cd,All_Partition>(maxLevel());
02453     }
02454 
02456   template<int cd>
02457   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
02458     {
02459       return levelend<cd,All_Partition>(maxLevel());
02460     }
02461 
02463   int overlapSize (int level, int codim) const
02464   {
02465         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02466         return g.overlap();
02467   }
02468 
02470   int overlapSize (int codim) const
02471   {
02472         YGLI g = MultiYGrid<dim,ctype>::begin(maxLevel());
02473         return g.overlap();
02474   }
02475 
02477   int ghostSize (int level, int codim) const
02478   {
02479         return 0;
02480   }
02481 
02483   int ghostSize (int codim) const
02484   {
02485         return 0;
02486   }
02487 
02489   int size (int level, int codim) const
02490   {
02491     return sizes[level][codim];
02492   }
02493 
02495   int size (int codim) const
02496   {
02497     return sizes[maxLevel()][codim];
02498   }
02499 
02501   int size (int level, GeometryType type) const
02502   {
02503       return (type.isCube()) ? sizes[level][dim-type.dim()] : 0;
02504   }
02505 
02507   int size (GeometryType type) const
02508   {
02509     return size(maxLevel(),type);
02510   }
02511 
02516   template<class DataHandleImp, class DataType>
02517   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir, int level) const
02518   {
02519     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
02520   }
02521 
02526   template<class DataHandleImp, class DataType>
02527   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir) const
02528   {
02529     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
02530   }
02531 
02536   template<class DataHandle, int codim>
02537   void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
02538   {
02539     // check input
02540     if (!data.contains(dim,codim)) return; // should have been checked outside
02541 
02542     // data types
02543     typedef typename DataHandle::DataType DataType;
02544 
02545     // access to grid level
02546     YGLI g = MultiYGrid<dim,ctype>::begin(level);
02547 
02548     // find send/recv lists or throw error
02549     const std::deque<IS>* sendlist=0;
02550     const std::deque<IS>* recvlist=0;
02551     if (codim==0) // the elements
02552       {
02553         if (iftype==InteriorBorder_InteriorBorder_Interface)
02554           return; // there is nothing to do in this case
02555         if (iftype==InteriorBorder_All_Interface)
02556           {
02557             sendlist = &g.send_cell_interior_overlap();
02558             recvlist = &g.recv_cell_overlap_interior();
02559           }
02560         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface || iftype==All_All_Interface)
02561           {
02562             sendlist = &g.send_cell_overlap_overlap();
02563             recvlist = &g.recv_cell_overlap_overlap();
02564           }
02565       }
02566     if (codim==dim) // the vertices
02567       {
02568         if (iftype==InteriorBorder_InteriorBorder_Interface)
02569           {
02570             sendlist = &g.send_vertex_interiorborder_interiorborder();
02571             recvlist = &g.recv_vertex_interiorborder_interiorborder();
02572           }
02573 
02574         if (iftype==InteriorBorder_All_Interface)
02575           {
02576             sendlist = &g.send_vertex_interiorborder_overlapfront();
02577             recvlist = &g.recv_vertex_overlapfront_interiorborder();
02578           }
02579         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface)
02580           {
02581             sendlist = &g.send_vertex_overlap_overlapfront();
02582             recvlist = &g.recv_vertex_overlapfront_overlap();
02583           }
02584         if (iftype==All_All_Interface)
02585           {
02586             sendlist = &g.send_vertex_overlapfront_overlapfront();
02587             recvlist = &g.recv_vertex_overlapfront_overlapfront();
02588           }
02589       }
02590     if (codim>0 && codim<dim)
02591       {
02592         DUNE_THROW(GridError, "interface communication not implemented");
02593       }
02594 
02595     // change communication direction?
02596     if (dir==BackwardCommunication)
02597       std::swap(sendlist,recvlist);
02598 
02599     int cnt;
02600 
02601     // Size computation (requires communication if variable size)
02602     std::map<int,int> send_size;      // map rank to total number of objects (of type DataType) to be sent
02603     std::map<int,int> recv_size;      // map rank to total number of objects (of type DataType) to be recvd
02604     std::map<int,size_t*> send_sizes; // map rank to array giving number of objects per entity to be sent
02605     std::map<int,size_t*> recv_sizes; // map rank to array giving number of objects per entity to be recvd
02606     if (data.fixedsize(dim,codim))
02607       {
02608         // fixed size: just take a dummy entity, size can be computed without communication
02609         cnt=0;
02610         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02611           {
02612             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02613               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02614             send_size[cnt] = is->grid.totalsize() * data.size(*it);
02615             cnt++;
02616           }
02617         cnt=0;
02618         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02619           {
02620             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02621               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02622             recv_size[cnt] = is->grid.totalsize() * data.size(*it);
02623             cnt++;
02624           }
02625       }
02626     else
02627       {
02628         // variable size case: sender side determines the size
02629         cnt=0;
02630         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02631           {
02632             // allocate send buffer for sizes per entitiy
02633             size_t *buf = new size_t[is->grid.totalsize()];
02634             send_sizes[cnt] = buf;
02635 
02636             // loop over entities and ask for size
02637             int i=0; size_t n=0;
02638             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02639               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02640             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02641               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02642             for ( ; it!=tsubend; ++it)
02643               {
02644                 buf[i] = data.size(*it);
02645                 n += buf[i];
02646                 i++;
02647               }
02648 
02649             // now we know the size for this rank
02650             send_size[cnt] = n;
02651 
02652             // hand over send request to torus class
02653             MultiYGrid<dim,ctype>::torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
02654             cnt++;
02655           }
02656 
02657         // allocate recv buffers for sizes and store receive request
02658         cnt=0;
02659         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02660           {
02661             // allocate recv buffer
02662             size_t *buf = new size_t[is->grid.totalsize()];
02663             recv_sizes[cnt] = buf;
02664             
02665             // hand over recv request to torus class
02666             MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
02667             cnt++;
02668           }
02669 
02670         // exchange all size buffers now
02671         MultiYGrid<dim,ctype>::torus().exchange();
02672 
02673         // release send size buffers
02674         cnt=0;
02675         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02676           {
02677             delete[] send_sizes[cnt];
02678             send_sizes[cnt] = 0;
02679             cnt++;
02680           }
02681 
02682         // process receive size buffers
02683         cnt=0;
02684         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02685           {
02686             // get recv buffer
02687             size_t *buf = recv_sizes[cnt];
02688 
02689             // compute total size
02690             size_t n=0;
02691             for (int i=0; i<is->grid.totalsize(); ++i)
02692               n += buf[i];
02693 
02694             // ... and store it
02695             recv_size[cnt] = n;
02696           }
02697       }
02698 
02699 
02700     // allocate & fill the send buffers & store send request
02701     std::map<int,DataType*> sends; // store pointers to send buffers
02702     cnt=0;
02703     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02704       {
02705 //      std::cout << "[" << this->comm().rank() << "] "
02706 //                << " send " << " dest=" << is->rank
02707 //                << " size=" << send_size[cnt] 
02708 //                << std::endl;
02709 
02710         // allocate send buffer
02711         DataType *buf = new DataType[send_size[cnt]];
02712 
02713         // remember send buffer
02714         sends[cnt] = buf;
02715 
02716         // make a message buffer
02717         MessageBuffer<DataType> mb(buf);
02718 
02719         // fill send buffer; iterate over cells in intersection
02720         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02721           it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02722         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02723           tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02724         for ( ; it!=tsubend; ++it)
02725           data.gather(mb,*it);
02726 
02727         // hand over send request to torus class
02728         MultiYGrid<dim,ctype>::torus().send(is->rank,buf,send_size[cnt]*sizeof(DataType));
02729         cnt++;
02730       }
02731 
02732     // allocate recv buffers and store receive request
02733     std::map<int,DataType*> recvs; // store pointers to send buffers
02734     cnt=0;
02735     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02736       {
02737 //      std::cout << "[" << this->comm().rank() << "] "
02738 //                << " recv " << "  src=" << is->rank
02739 //                << " size=" << recv_size[cnt] 
02740 //                << std::endl;
02741 
02742         // allocate recv buffer
02743         DataType *buf = new DataType[recv_size[cnt]];
02744 
02745         // remember recv buffer
02746         recvs[cnt] = buf;
02747 
02748         // hand over recv request to torus class
02749         MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,recv_size[cnt]*sizeof(DataType));
02750         cnt++;
02751       }
02752 
02753     // exchange all buffers now
02754     MultiYGrid<dim,ctype>::torus().exchange();
02755 
02756     // release send buffers
02757     cnt=0;
02758     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02759     {
02760       delete[] sends[cnt];
02761       sends[cnt] = 0;
02762     }
02763 
02764     // process receive buffers and delete them
02765     cnt=0;
02766     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02767       {
02768         // get recv buffer
02769         DataType *buf = recvs[cnt];
02770 
02771         // make a message buffer
02772         MessageBuffer<DataType> mb(buf);
02773 
02774         // copy data from receive buffer; iterate over cells in intersection
02775         if (data.fixedsize(dim,codim))
02776           {
02777             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02778               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02779             size_t n=data.size(*it);
02780             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02781               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02782             for ( ; it!=tsubend; ++it)
02783               data.scatter(mb,*it,n);
02784           }
02785         else
02786           {
02787             int i=0;
02788             size_t *sbuf = recv_sizes[cnt];
02789             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02790               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02791             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02792               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02793             for ( ; it!=tsubend; ++it)
02794               data.scatter(mb,*it,sbuf[i++]);
02795             delete[] sbuf;
02796           }
02797 
02798         // delete buffer
02799         delete[] buf; // hier krachts !
02800         cnt++;
02801       }
02802   }
02803 
02804   // The new index sets from DDM 11.07.2005
02805   const typename Traits::GlobalIdSet& globalIdSet() const
02806   {
02807     return *(theglobalidset[0]);
02808   }
02809   
02810   const typename Traits::LocalIdSet& localIdSet() const
02811   {
02812     return *(theglobalidset[0]);
02813   }
02814 
02815   const typename Traits::LevelIndexSet& levelIndexSet(int level) const
02816   {
02817     if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
02818     return *(indexsets[level]);
02819   }
02820 
02821   const typename Traits::LeafIndexSet& leafIndexSet() const
02822   {
02823     return *(theleafindexset[0]);
02824   }
02825 
02826 #if HAVE_MPI
02827 
02829   const CollectiveCommunication<MPI_Comm>& comm () const
02830   {
02831     return ccobj;
02832   }
02833 #else
02834 
02836   const CollectiveCommunication<YaspGrid>& comm () const
02837   {
02838     return ccobj;
02839   }
02840 #endif
02841 
02842   YaspIntersectionIterator<const YaspGrid<dim, dim> >&
02843   getRealIntersectionIterator(typename Traits::LevelIntersectionIterator& it)
02844   {
02845     return this->getRealImplementation(it);
02846   }
02847 
02848   const YaspIntersectionIterator<const YaspGrid<dim, dim> >&
02849   getRealIntersectionIterator(const typename Traits::LevelIntersectionIterator& it) const
02850   {
02851     return this->getRealImplementation(it);
02852   }
02853 
02854 private:
02855 
02856 #if HAVE_MPI
02857   CollectiveCommunication<MPI_Comm> ccobj; 
02858 #else
02859   CollectiveCommunication<YaspGrid> ccobj; 
02860 #endif
02861 
02862   std::vector<YaspLevelIndexSet<const YaspGrid<dim,dim> >*> indexsets;
02863   std::vector<YaspLeafIndexSet<const YaspGrid<dim,dim> >*> theleafindexset;
02864   std::vector<YaspGlobalIdSet<const YaspGrid<dim,dim> >*> theglobalidset;
02865 
02866   // Index classes need access to the real entity
02867   friend class Dune::YaspLevelIndexSet<const Dune::YaspGrid<dim,dim> >;
02868   friend class Dune::YaspLeafIndexSet<const Dune::YaspGrid<dim,dim> >;
02869   friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim,dim> >;
02870 
02871   template<class T>
02872   void deallocatePointers(T& container)
02873   {
02874     typedef typename T::iterator Iterator;
02875     
02876     for(Iterator entry=container.begin(); entry != container.end(); ++entry)
02877       delete (*entry);
02878   }
02879   
02880   template<int codim>
02881   YaspEntity<codim,dim,const YaspGrid<dim,dim> >& 
02882   getRealEntity(typename Traits::template Codim<codim>::Entity& e )
02883   {
02884     return this->getRealImplementation(e);
02885   }
02886   
02887   template<int codim>
02888   const YaspEntity<codim,dim,const YaspGrid<dim,dim> >& 
02889   getRealEntity(const typename Traits::template Codim<codim>::Entity& e ) const
02890   {
02891     return this->getRealImplementation(e);
02892   }
02893 
02894   template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
02895   friend class Entity;
02896 
02897   template<class DT>
02898   class MessageBuffer {
02899   public:
02900     // Constructor
02901     MessageBuffer (DT *p)
02902     {
02903       a=p;
02904       i=0;
02905       j=0;
02906     }
02907 
02908     // write data to message buffer, acts like a stream !
02909     template<class Y>
02910     void write (const Y& data)
02911     {
02912       dune_static_assert(( is_same<DT,Y>::value ), "DataType missmatch");
02913       a[i++] = data;
02914     }
02915     
02916     // read data from message buffer, acts like a stream !
02917     template<class Y>
02918     void read (Y& data) const 
02919     {
02920       dune_static_assert(( is_same<DT,Y>::value ), "DataType missmatch");
02921       data = a[j++];
02922     }
02923 
02924   private:
02925     DT *a;
02926     int i;
02927     mutable int j;
02928   };
02929 
02930   void setsizes ()
02931   {
02932     for (YGLI g=MultiYGrid<dim,ctype>::begin(); g!=MultiYGrid<dim,ctype>::end(); ++g)
02933       {
02934         // codim 0 (elements)
02935         sizes[g.level()][0] = 1;
02936         for (int i=0; i<dim; ++i) 
02937           sizes[g.level()][0] *= g.cell_overlap().size(i);
02938 
02939         // codim 1 (faces)
02940         if (dim>1)
02941           {
02942             sizes[g.level()][1] = 0;
02943             for (int i=0; i<dim; ++i) 
02944               {
02945                 int s=g.cell_overlap().size(i)+1;
02946                 for (int j=0; j<dim; ++j)
02947                   if (j!=i)
02948                     s *= g.cell_overlap().size(j);
02949                 sizes[g.level()][1] += s;
02950               }
02951           }
02952 
02953         // codim dim-1 (edges)
02954         if (dim>2)
02955           {
02956             sizes[g.level()][dim-1] = 0;
02957             for (int i=0; i<dim; ++i) 
02958               {
02959                 int s=g.cell_overlap().size(i);
02960                 for (int j=0; j<dim; ++j)
02961                   if (j!=i)
02962                     s *= g.cell_overlap().size(j)+1;
02963                 sizes[g.level()][dim-1] += s;
02964               }
02965           }
02966 
02967         // codim dim (vertices)
02968         sizes[g.level()][dim] = 1;
02969         for (int i=0; i<dim; ++i) 
02970           sizes[g.level()][dim] *= g.vertex_overlapfront().size(i);
02971       }
02972   }
02973 
02975   template<int cd, PartitionIteratorType pitype>
02976   YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
02977   {
02978         dune_static_assert( cd == dim || cd == 0 ,
02979           "YaspGrid only supports Entities with codim=dim and codim=0");
02980         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02981         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
02982         if (cd==0) // the elements
02983           {
02984                 if (pitype<=InteriorBorder_Partition) 
02985                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_interior().tsubbegin());
02986                 if (pitype<=All_Partition)            
02987                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_overlap().tsubbegin());
02988           }
02989         if (cd==dim) // the vertices
02990           {
02991                 if (pitype==Interior_Partition)       
02992                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interior().tsubbegin());
02993                 if (pitype==InteriorBorder_Partition) 
02994                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interiorborder().tsubbegin());
02995                 if (pitype==Overlap_Partition)        
02996                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlap().tsubbegin());
02997                 if (pitype<=All_Partition)            
02998                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlapfront().tsubbegin());
02999           }
03000         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
03001   }
03002 
03004   template<int cd, PartitionIteratorType pitype>
03005   YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
03006   {
03007         dune_static_assert( cd == dim || cd == 0 ,
03008           "YaspGrid only supports Entities with codim=dim and codim=0");
03009         YGLI g = MultiYGrid<dim,ctype>::begin(level);
03010         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
03011         if (cd==0) // the elements
03012           {
03013                 if (pitype<=InteriorBorder_Partition) 
03014                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_interior().tsubend());
03015                 if (pitype<=All_Partition)            
03016                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_overlap().tsubend());
03017           }
03018         if (cd==dim) // the vertices
03019           {
03020                 if (pitype==Interior_Partition)       
03021                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interior().tsubend());
03022                 if (pitype==InteriorBorder_Partition) 
03023                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interiorborder().tsubend());
03024                 if (pitype==Overlap_Partition)        
03025                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlap().tsubend());
03026                 if (pitype<=All_Partition)            
03027                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlapfront().tsubend());
03028           }
03029         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
03030   }
03031 
03032   int sizes[MAXL][dim+1]; // total number of entities per level and codim
03033   bool keep_ovlp;
03034 };
03035 
03036 namespace Capabilities
03037 {
03038 
03050   template<int dim, int dimw>
03051   struct hasEntity< YaspGrid<dim,dimw>, 0 >
03052   {
03053     static const bool v = true;
03054   };
03055   
03059   template<int dim, int dimw>
03060   struct hasEntity< YaspGrid<dim,dimw>, dim >
03061   {
03062     static const bool v = true;
03063   };
03064 
03068   template<int dim, int dimw>
03069   struct isParallel< YaspGrid<dim,dimw> >
03070   {
03071     static const bool v = true;
03072   };
03073 
03077   template<int dim, int dimw>
03078   struct isLevelwiseConforming< YaspGrid<dim,dimw> >
03079   {
03080     static const bool v = true;
03081   };
03082 
03086   template<int dim, int dimw>
03087   struct isLeafwiseConforming< YaspGrid<dim,dimw> >
03088   {
03089     static const bool v = true;
03090   };
03091 
03095   template<int dim, int dimw>
03096   struct IsUnstructured< YaspGrid<dim,dimw> >
03097   {
03098     static const bool v = false;
03099   };
03100 
03104   template<int dim, int dimw>
03105   struct hasHangingNodes< YaspGrid<dim,dimw> >
03106   {
03107     static const bool v = false;
03108   };
03109 
03110 }
03111 
03112 } // end namespace
03113 
03114 
03115 #endif
03116 

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