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 
00009 #include <dune/grid/common/grid.hh>     // the grid base classes
00010 #include <dune/grid/yaspgrid/grids.hh>  // the yaspgrid base classes
00011 #include <dune/common/stack.hh> // the stack class
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   bool checkInside (const FieldVector<ctype, mydim>& local) const
00243   {
00244         for (int i=0; i<mydim; i++)
00245           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00246         return true;
00247   }
00248 
00250   YaspGeometry (const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, int& m)
00251         : midpoint(p), extension(h), missing(m)
00252   {
00253         if (cdim!=mydim+1)
00254           DUNE_THROW(GridError, "general YaspGeometry assumes cdim=mydim+1");
00255   }
00256 
00258   void print (std::ostream& s) const
00259   {
00260         s << "YaspGeometry<"<<mydim<<","<<cdim<< "> ";
00261         s << "midpoint";
00262         for (int i=0; i<cdim; i++)
00263                   s << " " << midpoint[i];
00264                 s << " extension";
00265         for (int i=0; i<cdim; i++)
00266                   s << " " << extension[i];
00267                 s << " missing is " << missing;
00268   }
00269 private:
00270   // the element is fully defined by its midpoint the extension
00271   // in each direction and the missing direction.
00272   // References are used because this information
00273   // is known outside the element in many cases.
00274   // Note cdim==cdim+1
00275 
00276   // IMPORTANT midpoint and extension are references,
00277   // YaspGeometry can't be copied
00278   const FieldVector<ctype, cdim> & midpoint;    // the midpoint
00279   const FieldVector<ctype, cdim> & extension;   // the extension
00280   int& missing;                         // the missing, i.e. constant direction
00281 
00282   // In addition we need memory in order to return references.
00283   // Possibly we should change this in the interface ...
00284   mutable FieldMatrix<ctype,  Power_m_p<2,mydim>::power, cdim> coord_;   // the coordinates 
00285 
00286   const YaspGeometry<mydim,cdim,GridImp>&
00287   operator = (const YaspGeometry<mydim,cdim,GridImp>& g);
00288   
00289 };
00290 
00291 
00292 
00294 template<int mydim, class GridImp>
00295 class YaspGeometry<mydim,mydim,GridImp> : public GeometryDefaultImplementation<mydim,mydim,GridImp,YaspGeometry>
00296 {
00297 public:
00298   typedef typename GridImp::ctype ctype;
00300   typedef Geometry<mydim,mydim,GridImp,Dune::YaspGeometry> ReferenceGeometry;
00301   
00303   GeometryType type () const
00304   {
00305       return GeometryType(GeometryType::cube,mydim);
00306   }
00307 
00309   int corners () const
00310   {
00311         return 1<<mydim;
00312   }
00313 
00315   const FieldVector<ctype, mydim>& operator[] (int i) const
00316   {
00317         assert( i >= 0 && i < (int) coord_.N() ); 
00318         FieldVector<ctype, mydim>& c = coord_[i];
00319         for (int k=0; k<mydim; k++)
00320           if (i&(1<<k))
00321                 c[k] = midpoint[k]+0.5*extension[k]; // kth bit is 1 in i
00322           else
00323                 c[k] = midpoint[k]-0.5*extension[k]; // kth bit is 0 in i
00324         return c;
00325   }
00326 
00328   FieldVector<ctype, mydim> global (const FieldVector<ctype, mydim>& local) const
00329   {
00330       FieldVector<ctype,mydim> g;
00331         for (int k=0; k<mydim; k++)
00332           g[k] = midpoint[k] + (local[k]-0.5)*extension[k];
00333         return g;
00334   }
00335 
00337   FieldVector<ctype, mydim> local (const FieldVector<ctype,mydim>& global) const
00338   {
00339       FieldVector<ctype, mydim> l; // result
00340         for (int k=0; k<mydim; k++)
00341           l[k] = (global[k]-midpoint[k])/extension[k] + 0.5;
00342         return l;
00343   }
00344 
00347   ctype integrationElement (const FieldVector<ctype, mydim>& local) const
00348   {
00349     return volume();
00350   }
00351 
00353   ctype volume () const
00354   {
00355     ctype vol=1.0;
00356     for (int k=0; k<mydim; k++) vol *= extension[k];
00357     return vol;
00358   }
00359 
00361   FieldMatrix<ctype,mydim,mydim>& jacobianInverseTransposed (const FieldVector<ctype, mydim>& local) const
00362   {
00363         for (int i=0; i<mydim; ++i)
00364           {
00365                 Jinv[i] = 0.0;                // set column to zero
00366                 Jinv[i][i] = 1.0/extension[i]; // set diagonal element
00367           }
00368         return Jinv;
00369   }
00370 
00372   bool checkInside (const FieldVector<ctype, mydim>& local) const
00373   {
00374         for (int i=0; i<mydim; i++)
00375           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00376         return true;
00377   }
00378 
00380   YaspGeometry (const FieldVector<ctype, mydim>& p, const FieldVector<ctype, mydim>& h)
00381         : midpoint(p), extension(h)
00382   {}
00383 
00385   void print (std::ostream& s) const
00386   {
00387         s << "YaspGeometry<"<<mydim<<","<<mydim<< "> ";
00388         s << "midpoint";
00389         for (int i=0; i<mydim; i++)
00390                   s << " " << midpoint[i];
00391                 s << " extension";
00392         for (int i=0; i<mydim; i++)
00393                   s << " " << extension[i];
00394   }
00395 
00396 private:
00397   // the element is fully defined by midpoint and the extension
00398   // in each direction. References are used because this information
00399   // is known outside the element in many cases.
00400   // Note mydim==cdim
00401 
00402   // IMPORTANT midpoint and extension are references,
00403   // YaspGeometry can't be copied
00404   const FieldVector<ctype, mydim> & midpoint;   // the midpoint
00405   const FieldVector<ctype, mydim> & extension;  // the extension
00406 
00407   // In addition we need memory in order to return references.
00408   // Possibly we should change this in the interface ...
00409   mutable FieldMatrix<ctype,mydim,mydim> Jinv;   // the jacobian inverse
00410   mutable FieldMatrix<ctype,  Power_m_p<2,mydim>::power, mydim> coord_;   // the coordinates 
00411 
00412   // disable copy
00413   const YaspGeometry<mydim,mydim,GridImp>&
00414   operator = (const YaspGeometry<mydim,mydim,GridImp>& g);
00415 };
00416 
00418 template<int cdim, class GridImp>
00419 class YaspGeometry<0,cdim,GridImp> : public GeometryDefaultImplementation<0,cdim,GridImp,YaspGeometry>
00420 {
00421 public:
00422   typedef typename GridImp::ctype ctype;
00423   
00425   GeometryType type () const
00426   {
00427       return GeometryType(GeometryType::cube,0);
00428   }
00429 
00431   int corners () const
00432   {
00433         return 1;
00434   }
00435 
00437   const FieldVector<ctype, cdim>& operator[] (int i) const
00438   {
00439         return position;
00440   }
00441 
00444   ctype integrationElement (const FieldVector<ctype, 0>& local) const
00445   {
00446         return 1.0;
00447   }
00448 
00450   YaspGeometry (const FieldVector<ctype, cdim>& p) : position(p)
00451   {}
00452 
00454   void print (std::ostream& s) const
00455   {
00456         s << "YaspGeometry<"<<0<<","<<cdim<< "> ";
00457         s << "position " << position;
00458   }
00459 
00460 private:
00461   // IMPORTANT position is a reference,
00462   // YaspGeometry can't be copied
00463   const FieldVector<ctype, cdim> & position; 
00464 
00465   const YaspGeometry<0,cdim,GridImp>&
00466   operator = (const YaspGeometry<0,cdim,GridImp>& g);
00467 };
00468 
00469 // operator<< for all YaspGeometrys
00470 template <int mydim, int cdim, class GridImp>
00471 inline
00472 std::ostream& operator<< (std::ostream& s, YaspGeometry<mydim,cdim,GridImp>& e)
00473 {
00474   e.print(s);
00475   return s;
00476 }
00477 
00478 //========================================================================
00486 //========================================================================
00487 
00488 template<int codim, int dim, class GridImp>
00489 class YaspSpecialEntity :
00490   public GridImp::template Codim<codim>::Entity
00491 {
00492 public:
00493   typedef typename GridImp::ctype ctype;
00494   
00495   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00496   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00497 
00498   YaspSpecialEntity(const YGLI& g, const TSI& it) :
00499     GridImp::template Codim<codim>::Entity (YaspEntity<codim, dim, GridImp>(g,it))
00500     {};
00501   YaspSpecialEntity(const YaspEntity<codim, dim, GridImp>& e) :
00502     GridImp::template Codim<codim>::Entity (e)
00503     {};
00504   const TSI& transformingsubiterator () const
00505   {
00506         return this->realEntity.transformingsubiterator();
00507   }
00508   const YGLI& gridlevel () const
00509   {
00510         return this->realEntity.gridlevel();
00511   }
00512 };
00513 
00514 template<int codim, int dim, class GridImp>
00515 class YaspEntity 
00516 :  public EntityDefaultImplementation <codim,dim,GridImp,YaspEntity>
00517 {
00518 public:
00519   typedef typename GridImp::ctype ctype;
00520   
00521   typedef typename GridImp::template Codim<codim>::Geometry Geometry;
00523   int level () const
00524   {
00525         DUNE_THROW(GridError, "YaspEntity not implemented");
00526   }
00527 
00529   int index () const
00530   {
00531         DUNE_THROW(GridError, "YaspEntity not implemented");
00532   }
00533 
00535   const Geometry& geometry () const
00536   {
00537         DUNE_THROW(GridError, "YaspEntity not implemented");
00538   }
00539 
00541   PartitionType partitionType () const
00542   {
00543         DUNE_THROW(GridError, "YaspEntity not implemented");
00544   }
00545 
00546   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00547   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00548   YaspEntity (const YGLI& g, const TSI& it)
00549   {
00550         DUNE_THROW(GridError, "YaspEntity not implemented");
00551   }
00552 
00553   // IndexSets needs access to the private index methods
00554   friend class Dune::YaspLevelIndexSet<GridImp>;
00555   friend class Dune::YaspLeafIndexSet<GridImp>;
00556   friend class Dune::YaspGlobalIdSet<GridImp>;
00557   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00558 
00560   PersistentIndexType persistentIndex () const 
00561   { 
00562     DUNE_THROW(GridError, "YaspEntity not implemented");
00563   }
00564 
00566   int compressedIndex () const 
00567   {
00568     DUNE_THROW(GridError, "YaspEntity not implemented");
00569   }
00570 
00572   int compressedLeafIndex () const 
00573   {
00574     DUNE_THROW(GridError, "YaspEntity not implemented");
00575   }
00576 };
00577 
00578 
00579 // specialization for codim=0
00580 template<int dim, class GridImp>
00581 class YaspEntity<0,dim,GridImp> 
00582 : public EntityDefaultImplementation <0,dim,GridImp,YaspEntity>
00583 {
00584   enum { dimworld = GridImp::dimensionworld };
00585 public:
00586   typedef typename GridImp::ctype ctype;
00587   
00588   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00589   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00590 
00591   typedef YaspSpecialGeometry<dim-0,dim,GridImp> SpecialGeometry;
00592   
00593   typedef typename GridImp::template Codim<0>::Geometry Geometry;
00594   template <int cd>
00595   struct Codim
00596   {
00597     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
00598   };
00599   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
00600   typedef typename GridImp::template Codim<0>::LevelIntersectionIterator IntersectionIterator;
00601   typedef typename GridImp::template Codim<0>::LevelIntersectionIterator LevelIntersectionIterator;
00602   typedef typename GridImp::template Codim<0>::LeafIntersectionIterator  LeafIntersectionIterator;
00603   typedef typename GridImp::template Codim<0>::HierarchicIterator HierarchicIterator;
00604 
00606   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00607 
00609   typedef typename YGrid<dim,ctype>::iTupel iTupel;
00610 
00611   // constructor
00612   YaspEntity (const YGLI& g, const TSI& it)
00613     : _it(it), _g(g), _geometry(it.position(),it.meshsize())
00614   {
00615   }
00616 
00618   int level () const {return _g.level();}
00619 
00621   int index () const { return _it.superindex();} // superindex works also for iteration over subgrids
00622 
00624   int globalIndex () const {
00625     return _g.cell_global().index(_it.coord());
00626   }
00627 
00629   PartitionType partitionType () const
00630   {
00631         if (_g.cell_interior().inside(_it.coord())) return InteriorEntity;
00632         if (_g.cell_overlap().inside(_it.coord())) return OverlapEntity;
00633         return GhostEntity;
00634   }
00635 
00637   const Geometry& geometry () const { return _geometry; }
00638 
00641   template<int cc> int count () const
00642   {
00643     if (cc==dim) return 1<<dim;
00644     if (cc==1) return 2*dim;
00645     if (cc==dim-1) return dim*(1<<(dim-1));
00646     if (cc==0) return 1;
00647     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00648   }
00649 
00652   template<int cc>
00653   typename Codim<cc>::EntityPointer entity (int i) const
00654   {
00655         IsTrue< ( cc == dim || cc == 0 ) >::yes();
00656         // coordinates of the cell == coordinates of lower left corner
00657         if (cc==dim)
00658           {
00659                 iTupel coord = _it.coord();
00660 
00661                 // get corner from there
00662                 for (int k=0; k<dim; k++)
00663                   if (i&(1<<k)) (coord[k])++;
00664 
00665                 return YaspLevelIterator<cc,All_Partition,GridImp>(_g,_g.vertex_overlapfront().tsubbegin(coord));
00666           }
00667         if (cc==0)
00668           {
00669                 return YaspLevelIterator<cc,All_Partition,GridImp>(_g,_it);
00670           }
00671         DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00672   }
00673 
00675   EntityPointer father () const
00676   {
00677         // check if coarse level exists
00678         if (_g.level()<=0)
00679           DUNE_THROW(GridError, "tried to call father on level 0");
00680 
00681         // yes, get iterator to it
00682         YGLI cg = _g.coarser();
00683 
00684         // coordinates of the cell
00685         iTupel coord = _it.coord();
00686 
00687         // get coordinates on next coarser level
00688         for (int k=0; k<dim; k++) coord[k] = coord[k]/2;
00689 
00690         return YaspLevelIterator<0,All_Partition,GridImp>(cg,cg.cell_overlap().tsubbegin(coord));
00691   }
00692 
00702   const Geometry& geometryInFather () const
00703   {
00704     // determine which son we are
00705     int son = 0;
00706     for (int k=0; k<dim; k++)
00707       if (_it.coord(k)%2)
00708         son += (1<<k);
00709     
00710     // configure one of the 2^dim transformations
00711     return YaspFatherRelativeLocalElement<dim,GridImp>::getson(son);
00712   }
00713 
00714   const TSI& transformingsubiterator () const
00715   {
00716         return _it;
00717   }
00718 
00719   const YGLI& gridlevel () const
00720   {
00721         return _g;
00722   }
00723 
00724   bool isLeaf() const
00725   {
00726     return (_g.level() == _g.mg()->maxlevel());
00727   }
00728   
00730   IntersectionIterator ibegin () const
00731   {
00732         return YaspIntersectionIterator<GridImp>(*this,false);
00733   }
00734 
00736   LeafIntersectionIterator ileafbegin () const
00737   {
00738     // only if entity is leaf this iterator delivers intersections
00739     if (isLeaf())
00740       return ibegin();
00741     else
00742       return iend();
00743   }
00744 
00746   LevelIntersectionIterator ilevelbegin () const
00747   {
00748     return ibegin(); 
00749   }
00750 
00752   IntersectionIterator iend () const
00753   {
00754         return YaspIntersectionIterator<GridImp>(*this,true);
00755   }
00756 
00758   LeafIntersectionIterator ileafend () const
00759   {
00760     return iend();
00761   }
00762 
00764   LevelIntersectionIterator ilevelend () const
00765   {
00766     return iend();
00767   }
00768 
00773   HierarchicIterator hbegin (int maxlevel) const
00774   {
00775         return YaspHierarchicIterator<GridImp>(_g,_it,maxlevel);
00776   }
00777 
00779   HierarchicIterator hend (int maxlevel) const
00780   {
00781         return YaspHierarchicIterator<GridImp>(_g,_it,_g.level());
00782   }
00783 
00784 private:
00785   // IndexSets needs access to the private index methods
00786   friend class Dune::YaspLevelIndexSet<GridImp>;
00787   friend class Dune::YaspLeafIndexSet<GridImp>;
00788   friend class Dune::YaspGlobalIdSet<GridImp>;
00789 
00791   PersistentIndexType persistentIndex () const 
00792   { 
00793     // get size of global grid
00794     const iTupel& size =  _g.cell_global().size();
00795 
00796     // get coordinate correction for periodic boundaries
00797     int coord[dim];
00798     for (int i=0; i<dim; i++)
00799       {
00800         coord[i] = _it.coord(i);
00801         if (coord[i]<0) coord[i] += size[i];
00802         if (coord[i]>=size[i]) coord[i] -= size[i];
00803       }
00804 
00805     // encode codim
00806     PersistentIndexType id(0);
00807 
00808     // encode level
00809     id = id << yaspgrid_level_bits;
00810     id = id+PersistentIndexType(_g.level());
00811     
00812 
00813     // encode coordinates
00814     for (int i=dim-1; i>=0; i--)
00815       {
00816         id = id << yaspgrid_dim_bits;
00817         id = id+PersistentIndexType(coord[i]);
00818       }
00819     
00820     return id;
00821   }
00822 
00824   int compressedIndex () const 
00825   {
00826     return _it.superindex();
00827   }
00828 
00830   int compressedLeafIndex () const 
00831   {
00832     return _it.superindex();
00833   }
00834 
00836   template<int cc>
00837   PersistentIndexType subPersistentIndex (int i) const
00838   {
00839     if (cc==0)
00840       return persistentIndex();
00841 
00842     // get position of cell, note that global origin is zero
00843     // adjust for periodic boundaries
00844     int coord[dim];
00845     for (int k=0; k<dim; k++)
00846       {
00847         coord[k] = _it.coord(k);
00848         if (coord[k]<0) coord[k] += _g.cell_global().size(k);
00849         if (coord[k]>=_g.cell_global().size(k)) coord[k] -= _g.cell_global().size(k);
00850       }
00851     
00852     if (cc==dim)
00853       {
00854         // transform to vertex coordinates
00855         for (int k=0; k<dim; k++)
00856           if (i&(1<<k)) (coord[k])++;
00857 
00858         // determine min number of trailing zeroes
00859         int trailing = 1000;
00860         for (int i=0; i<dim; i++)
00861           {
00862             // count trailing zeros
00863             int zeros = 0;
00864             for (int j=0; j<_g.level(); j++)
00865               if (coord[i]&(1<<j))
00866                 break;
00867               else
00868                 zeros++;
00869             trailing = std::min(trailing,zeros);
00870           } 
00871 
00872         // determine the level of this vertex
00873         int level = _g.level()-trailing;
00874 
00875         // encode codim
00876         PersistentIndexType id(dim);
00877         
00878         // encode level
00879         id = id << yaspgrid_level_bits;
00880         id = id+PersistentIndexType(level);
00881         
00882         // encode coordinates
00883         for (int i=dim-1; i>=0; i--)
00884           {
00885             id = id << yaspgrid_dim_bits;
00886             id = id+PersistentIndexType(coord[i]>>trailing);
00887           }
00888         
00889         return id;
00890       }
00891 
00892     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
00893       {
00894         // Idea: Use the doubled grid to assign coordinates to faces
00895 
00896         // ivar is the direction that varies
00897         int ivar=i/2; 
00898 
00899         // compute position from cell position
00900         for (int k=0; k<dim; k++)
00901           coord[k] = coord[k]*2 + 1; // the doubled grid
00902         if (i%2) 
00903           coord[ivar] += 1;
00904         else
00905           coord[ivar] -= 1;
00906 
00907         // encode codim
00908         PersistentIndexType id(1);
00909         
00910         // encode level
00911         id = id << yaspgrid_level_bits;
00912         id = id+PersistentIndexType(_g.level());
00913         
00914         // encode coordinates
00915         for (int i=dim-1; i>=0; i--)
00916           {
00917             id = id << yaspgrid_dim_bits;
00918             id = id+PersistentIndexType(coord[i]);
00919           }
00920         
00921         return id;
00922       }
00923 
00924     if (cc==dim-1) // edges, exist only for dim>2
00925       {
00926         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
00927 
00928         // number of entities per direction
00929         int m=1<<(dim-1);
00930 
00931         // ifix is the direction that is fixed
00932         int ifix=(dim-1)-(i/m); 
00933 
00934         // compute position from cell position
00935         int bit=1;
00936         for (int k=0; k<dim; k++)
00937           {
00938             coord[k] = coord[k]*2+1; // cell position in doubled grid
00939             if (k==ifix) continue;
00940             if ((i%m)&bit) coord[k] += 1; else coord[k] -= 1;
00941             bit *= 2;
00942           } 
00943 
00944         // encode codim
00945         PersistentIndexType id(dim-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     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00962   }
00963 
00965   template<int cc>
00966   int subCompressedIndex (int i) const
00967   {
00968     if (cc==0)
00969       return compressedIndex();
00970     
00971     // get cell position relative to origin of local cell grid
00972     iTupel coord;
00973     for (int k=0; k<dim; ++k) 
00974       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
00975 
00976     if (cc==dim) // vertices
00977       {
00978         // transform cell coordinate to corner coordinate
00979         for (int k=0; k<dim; k++)
00980           if (i&(1<<k)) (coord[k])++;
00981 
00982         // do lexicographic numbering
00983         int index = coord[dim-1];
00984         for (int k=dim-2; k>=0; --k)
00985           index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
00986         return index;
00987       }
00988 
00989     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
00990       {
00991         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
00992 
00993         // ivar is the direction that varies
00994         int ivar=i/2; 
00995 
00996         // compute position from cell position
00997         if (i%2) coord[ivar] += 1; 
00998 
00999         // do lexicographic numbering
01000         int index = coord[dim-1];
01001         for (int k=dim-2; k>=0; --k)
01002           if (k==ivar)
01003             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01004           else
01005             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01006 
01007         // add size of all subsets for smaller directions
01008         for (int j=0; j<ivar; j++)
01009           {
01010             int n=_g.cell_overlap().size(j)+1;
01011             for (int l=0; l<dim; l++)
01012               if (l!=j) n *= _g.cell_overlap().size(l);
01013             index += n;
01014           }
01015 
01016         return index;
01017       }
01018 
01019     if (cc==dim-1) // edges, exist only for dim>2
01020       {
01021         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01022 
01023         // number of entities per direction
01024         int m=1<<(dim-1);
01025 
01026         // ifix is the direction that is fixed
01027         int ifix=(dim-1)-(i/m); 
01028 
01029         // compute position from cell position
01030         int bit=1;
01031         for (int k=0; k<dim; k++)
01032           {
01033             if (k==ifix) continue;
01034             if ((i%m)&bit) coord[k] += 1;
01035             bit *= 2;
01036           } 
01037 
01038         // do lexicographic numbering
01039         int index = coord[dim-1];
01040         for (int k=dim-2; k>=0; --k)
01041           if (k!=ifix)
01042             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01043           else
01044             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01045 
01046         // add size of all subsets for smaller directions
01047         for (int j=dim-1; j>ifix; j--)
01048           {
01049             int n=_g.cell_overlap().size(j);
01050             for (int l=0; l<dim; l++)
01051               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01052             index += n;
01053           }
01054 
01055         return index;
01056       }
01057 
01058     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01059   }
01060 
01062   template<int cc>
01063   int subCompressedLeafIndex (int i) const
01064   {
01065     if (cc==0)
01066       return compressedIndex();
01067     
01068     // get cell position relative to origin of local cell grid
01069     iTupel coord;
01070     for (int k=0; k<dim; ++k) 
01071       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
01072 
01073     if (cc==dim) // vertices
01074       {
01075         // transform cell coordinate to corner coordinate
01076         for (int k=0; k<dim; k++)
01077           if (i&(1<<k)) (coord[k])++;
01078 
01079         // move coordinates up to maxlevel
01080         for (int k=0; k<dim; k++)
01081           coord[k] = coord[k]<<(_g.mg()->maxlevel()-_g.level());
01082 
01083         // do lexicographic numbering
01084         int index = coord[dim-1];
01085         for (int k=dim-2; k>=0; --k)
01086           index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01087         return index;
01088       }
01089 
01090     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01091       {
01092         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
01093 
01094         // ivar is the direction that varies
01095         int ivar=i/2; 
01096 
01097         // compute position from cell position
01098         if (i%2) coord[ivar] += 1; 
01099 
01100         // do lexicographic numbering
01101         int index = coord[dim-1];
01102         for (int k=dim-2; k>=0; --k)
01103           if (k==ivar)
01104             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01105           else
01106             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01107 
01108         // add size of all subsets for smaller directions
01109         for (int j=0; j<ivar; j++)
01110           {
01111             int n=_g.cell_overlap().size(j)+1;
01112             for (int l=0; l<dim; l++)
01113               if (l!=j) n *= _g.cell_overlap().size(l);
01114             index += n;
01115           }
01116 
01117         return index;
01118       }
01119 
01120     if (cc==dim-1) // edges, exist only for dim>2
01121       {
01122         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01123 
01124         // number of entities per direction
01125         int m=1<<(dim-1);
01126 
01127         // ifix is the direction that is fixed
01128         int ifix=(dim-1)-(i/m); 
01129 
01130         // compute position from cell position
01131         int bit=1;
01132         for (int k=0; k<dim; k++)
01133           {
01134             if (k==ifix) continue;
01135             if ((i%m)&bit) coord[k] += 1;
01136             bit *= 2;
01137           } 
01138 
01139         // do lexicographic numbering
01140         int index = coord[dim-1];
01141         for (int k=dim-2; k>=0; --k)
01142           if (k!=ifix)
01143             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01144           else
01145             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01146 
01147         // add size of all subsets for smaller directions
01148         for (int j=dim-1; j>ifix; j--)
01149           {
01150             int n=_g.cell_overlap().size(j);
01151             for (int l=0; l<dim; l++)
01152               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01153             index += n;
01154           }
01155 
01156         return index;
01157       }
01158 
01159     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01160   }
01161 
01162 
01163   const TSI& _it;           // position in the grid level
01164   const YGLI& _g;           // access to grid level
01165   SpecialGeometry _geometry; // the element geometry
01166 };
01167 
01168 
01169 // specialization for codim=dim (vertex)
01170 template<int dim, class GridImp>
01171 class YaspEntity<dim,dim,GridImp> 
01172 : public EntityDefaultImplementation <dim,dim,GridImp,YaspEntity>
01173 {
01174   enum { dimworld = GridImp::dimensionworld };
01175 public:
01176   typedef typename GridImp::ctype ctype;
01177   
01178   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01179   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01180 
01181   typedef YaspSpecialGeometry<dim-dim,dim,GridImp> SpecialGeometry;
01182   
01183   typedef typename GridImp::template Codim<dim>::Geometry Geometry;
01184   template <int cd>
01185   struct Codim
01186   {
01187     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
01188   };
01189   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01190 
01192   typedef typename GridImp::PersistentIndexType PersistentIndexType;
01193 
01195   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01196 
01197   // constructor
01198   YaspEntity (const YGLI& g, const TSI& it)
01199         : _it(it), _g(g), _geometry(it.position())
01200   {  }
01201 
01203   int level () const {return _g.level();}
01204 
01206   int index () const {return _it.superindex();}
01207 
01209   int globalIndex () const { return _g.cell_global().index(_it.coord()); }
01210 
01211 
01213   const Geometry& geometry () const { return _geometry; }
01214 
01216   PartitionType partitionType () const
01217   {
01218         if (_g.vertex_interior().inside(_it.coord())) return InteriorEntity;
01219         if (_g.vertex_interiorborder().inside(_it.coord())) return BorderEntity;
01220         if (_g.vertex_overlap().inside(_it.coord())) return OverlapEntity;
01221         if (_g.vertex_overlapfront().inside(_it.coord())) return FrontEntity;
01222         return GhostEntity;
01223   }
01224 
01225 private:
01226   // IndexSets needs access to the private index methods
01227   friend class Dune::YaspLevelIndexSet<GridImp>;
01228   friend class Dune::YaspLeafIndexSet<GridImp>;
01229   friend class Dune::YaspGlobalIdSet<GridImp>;
01230 
01232   PersistentIndexType persistentIndex () const 
01233   { 
01234     // get coordinate and size of global grid
01235     const iTupel& size =  _g.vertex_global().size();
01236     int coord[dim];
01237 
01238     // correction for periodic boundaries
01239     for (int i=0; i<dim; i++)
01240       {
01241         coord[i] = _it.coord(i);
01242         if (coord[i]<0) coord[i] += size[i];
01243         if (coord[i]>=size[i]) coord[i] -= size[i];
01244       }
01245 
01246     // determine min number of trailing zeroes
01247     int trailing = 1000;
01248     for (int i=0; i<dim; i++)
01249       {
01250         // count trailing zeros
01251         int zeros = 0;
01252         for (int j=0; j<_g.level(); j++)
01253           if (coord[i]&(1<<j))
01254             break;
01255           else
01256             zeros++;
01257         trailing = std::min(trailing,zeros);
01258       } 
01259 
01260     // determine the level of this vertex
01261     int level = _g.level()-trailing;
01262 
01263     // encode codim
01264     PersistentIndexType id(dim);
01265 
01266     // encode level
01267     id = id << yaspgrid_level_bits;
01268     id = id+PersistentIndexType(level);
01269     
01270     // encode coordinates
01271     for (int i=dim-1; i>=0; i--)
01272       {
01273         id = id << yaspgrid_dim_bits;
01274         id = id+PersistentIndexType(coord[i]>>trailing);
01275       }
01276     
01277     return id;
01278   }
01279 
01281   int compressedIndex () const { return _it.superindex();}
01282 
01284   int compressedLeafIndex () const 
01285   { 
01286     if (_g.level()==_g.mg()->maxlevel())
01287       return _it.superindex();
01288 
01289     // not on leaf level, interpolate to finest grid
01290     int coord[dim];
01291     for (int i=0; i<dim; i++) coord[i] = _it.coord(i)-(_g).vertex_overlap().origin(i);
01292 
01293     // move coordinates up to maxlevel (multiply by 2 for each level
01294     for (int k=0; k<dim; k++)
01295       coord[k] = coord[k]*(1<<(_g.mg()->maxlevel()-_g.level()));
01296 
01297     // do lexicographic numbering
01298     int index = coord[dim-1];
01299     for (int k=dim-2; k>=0; --k)
01300       index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01301     return index;
01302   }
01303 
01304   const TSI& _it;                // position in the grid level
01305   const YGLI& _g;                // access to grid level
01306   SpecialGeometry  _geometry;     // the element geometry
01307   // temporary object
01308   mutable FieldVector<ctype, dim> loc;   // always computed before being returned
01309 };
01310 
01311 
01312 //========================================================================
01317 //========================================================================
01318 
01319 template<class GridImp>
01320 class YaspIntersectionIterator :
01321   public YaspEntityPointer<0,GridImp>,
01322   public IntersectionIteratorDefaultImplementation<GridImp,YaspIntersectionIterator>
01323 {
01324   enum { dim=GridImp::dimension };
01325   enum { dimworld=GridImp::dimensionworld };  
01326   typedef typename GridImp::ctype ctype;
01327   YaspIntersectionIterator();
01328 public:
01329   // types used from grids
01330   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01331   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01332   typedef typename GridImp::template Codim<0>::Entity Entity;
01333   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01334   typedef typename GridImp::template Codim<1>::Geometry Geometry;
01335   typedef typename GridImp::template Codim<1>::LocalGeometry LocalGeometry;
01336   typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01337   typedef YaspSpecialGeometry<dim-1,dimworld,GridImp> SpecialGeometry;
01338   typedef YaspSpecialGeometry<dim-1,dim,GridImp> SpecialLocalGeometry;
01339   
01341   void increment()
01342   {
01343         // check end
01344         if (_count==2*dim) return; // end iterator reached, we are done
01345         // update count, check end
01346         _count++;
01347 
01348         // update intersection iterator from current position
01349         if (_face==0) // direction remains valid
01350           {
01351                 _face = 1; // 0->1, _dir remains
01352 
01353                 // move transforming iterator
01354                 this->_it.move(_dir,2); // move two cells in positive direction
01355 
01356                 // make up faces
01357                 _pos_self_local[_dir] = 1.0;
01358                 _pos_nb_local[_dir] = 0.0;
01359                 _pos_world[_dir] += _myself.transformingsubiterator().meshsize(_dir);
01360 
01361                 // make up unit outer normal direction
01362                 _normal[_dir] = 1.0;
01363           }
01364         else // change direction
01365           {
01366                 // move transforming iterator
01367                 this->_it.move(_dir,-1); // move one cell back
01368                 if (_count==2*dim) return;
01369                 
01370                 // make up faces
01371                 _pos_self_local[_dir] = 0.5;
01372                 _pos_nb_local[_dir] = 0.5;
01373                 _pos_world[_dir] = _myself.transformingsubiterator().position(_dir);
01374 
01375                 // make up unit outer normal direction
01376                 _normal[_dir] = 0.0;
01377 
01378                 _face = 0;
01379                 _dir += 1;
01380 
01381                 // move transforming iterator
01382                 this->_it.move(_dir,-1); // move one cell in negative direction
01383 
01384                 // make up faces
01385                 _pos_self_local[_dir] = 0.0;
01386                 _pos_nb_local[_dir] = 1.0;
01387                 _pos_world[_dir] -= 0.5*_myself.transformingsubiterator().meshsize(_dir);
01388 
01389                 // make up unit outer normal direction
01390                 _normal[_dir] = -1.0;
01391           }
01392   }
01393 
01398   bool boundary () const
01399   {
01400     // The transforming iterator can be safely moved beyond the boundary.
01401     // So we only have to compare against the cell_global grid
01402     return (this->_it.coord(_dir)<_myself.gridlevel().cell_global().min(_dir)
01403             ||
01404             this->_it.coord(_dir)>_myself.gridlevel().cell_global().max(_dir));
01405   }
01406 
01408   bool neighbor () const
01409   {
01410     return
01411       (this->_it.coord(_dir)>=_myself.gridlevel().cell_overlap().min(_dir)
01412        &&
01413        this->_it.coord(_dir)<=_myself.gridlevel().cell_overlap().max(_dir));
01414     for (int d = 0; d < dim; d++)
01415     {
01416       if (this->_it.coord(_dir)<_myself.gridlevel().cell_overlap().min(_dir)
01417           ||
01418           this->_it.coord(_dir)>_myself.gridlevel().cell_overlap().max(_dir))
01419         return false;
01420     }
01421     return true;
01422   }
01423 
01426   EntityPointer inside() const
01427     {
01428       return _selfp;
01429     }
01430   
01433   EntityPointer outside() const
01434     {
01435       return *this;
01436     }
01437 
01440   int boundaryId() const
01441     {
01442       /*
01443       // this method returns 0 for 0th face which is wrong 
01444       if (this->_it.coord(_dir)<_myself.gridlevel().cell_global().min(_dir))
01445         return 2 * _dir + 1;
01446       if (this->_it.coord(_dir)>_myself.gridlevel().cell_global().max(_dir))
01447         return 2 * _dir + 2;
01448         */
01449       if(boundary()) return numberInSelf()+1;
01450       return 0;
01451     }
01452   
01454   FieldVector<ctype, dimworld> outerNormal (const FieldVector<ctype, dim-1>& local) const
01455   {
01456         return _normal;
01457   }
01458 
01459 
01461   FieldVector<ctype, dimworld> unitOuterNormal (const FieldVector<ctype, dim-1>& local) const
01462   {
01463         return _normal;
01464   }
01465 
01466 
01470   const LocalGeometry& intersectionSelfLocal () const
01471   {
01472         return _is_self_local;
01473   }
01474 
01478   const LocalGeometry& intersectionNeighborLocal () const
01479   {
01480         return _is_nb_local;
01481   }
01482 
01486   const Geometry& intersectionGlobal () const
01487   {
01488         return _is_global;
01489   }
01490 
01492   int numberInSelf () const
01493   {
01494         return _count;
01495   }
01496 
01498   int numberInNeighbor () const
01499   {
01500         return _count + 1-2*_face;
01501   }
01502 
01504   YaspIntersectionIterator (const YaspEntity<0,dim,GridImp>& myself, bool toend)
01505         : YaspEntityPointer<0,GridImp>(myself.gridlevel(),
01506                                        myself.transformingsubiterator()),
01507           _selfp(myself.gridlevel(),
01508                  myself.transformingsubiterator()),
01509           _myself(myself),
01510           _pos_self_local(0.5),
01511           _pos_nb_local(0.5),
01512           _pos_world(myself.transformingsubiterator().position()),
01513           _ext_local(1.0),
01514           _is_self_local(_pos_self_local,_ext_local,_dir),
01515           _is_nb_local(_pos_nb_local,_ext_local,_dir),
01516           _is_global(_pos_world,_myself.transformingsubiterator().meshsize(),_dir),
01517           _normal(0.0)
01518   {
01519         // making an end iterator?
01520         if (toend)
01521           {
01522                 // initialize end iterator
01523                 _count = 2*dim;
01524                 return;
01525           }
01526         // initialize to first neighbor
01527         _count = 0;
01528         _dir = 0;
01529         _face = 0;
01530 
01531         // move transforming iterator
01532         this->_it.move(_dir,-1);
01533 
01534         // make up faces
01535         _pos_self_local[0] = 0.0;
01536         _pos_nb_local[0] = 1.0;
01537         _pos_world[0] -= 0.5*_myself.transformingsubiterator().meshsize(0);
01538 
01539         // make up unit outer normal direction
01540         _normal[0] = -1.0;
01541   }
01542 
01544   YaspIntersectionIterator (const YaspIntersectionIterator& it)
01545         : YaspEntityPointer<0,GridImp>(it._g, it._it),
01546           _count(it._count),
01547           _dir(it._dir),
01548           _face(it._face),
01549           _selfp(it._selfp),
01550           _myself(it._myself),
01551           _pos_self_local(it._pos_self_local),
01552           _pos_nb_local(it._pos_nb_local),
01553           _pos_world(it._pos_world),
01554           _ext_local(it._ext_local),
01555     // Important: _is_* must be recreated -- not copied!
01556           _is_self_local(_pos_self_local,_ext_local,_dir),
01557           _is_nb_local(_pos_nb_local,_ext_local,_dir),
01558           _is_global(_pos_world,_myself.transformingsubiterator().meshsize(),_dir),
01559           _normal(it._normal)
01560   {}
01561 
01563   YaspIntersectionIterator & operator = (const YaspIntersectionIterator& it)
01564     {
01565       /* Assert same Iterator Context */
01566       if (! _selfp.equals(it._selfp))
01567         DUNE_THROW(GridError, "assignment of YaspIntersectionIterator "
01568                    << "with different inside Entity");
01569 
01570       /* Copy baseclass */
01571       YaspEntityPointer<0,GridImp>::operator=(it);
01572       
01573       /* Assign current position */
01574       _count = it._count;
01575       _dir = it._dir;
01576       _face = it._face;
01577       _pos_self_local = it._pos_self_local;
01578       _pos_nb_local = it._pos_nb_local;
01579       _pos_world = it._pos_world;
01580       _ext_local = it._ext_local;
01581       _normal = it._normal;
01582 
01583       return *this;
01584   }
01585 
01586 private:
01587   /* current position */
01588   int _count;                             
01589   int _dir;                               
01590   int _face;                              
01591   /* neighbouring Entity/Poiner (get automatically updated) */
01592   const YaspEntityPointer<0,GridImp> _selfp; 
01593   const YaspEntity<0,dim,GridImp>&
01594                           _myself;        
01595   /* precalculated data for current position */
01596   FieldVector<ctype, dim> _pos_self_local;
01597   FieldVector<ctype, dim> _pos_nb_local;  
01598   FieldVector<ctype, dimworld>_pos_world;     
01599   FieldVector<ctype, dim> _ext_local;     
01600   /* geometry objects (get automatically updated) */
01601   SpecialLocalGeometry    _is_self_local; 
01602   SpecialLocalGeometry    _is_nb_local;   
01603   SpecialGeometry         _is_global;     
01604   /* normal vector */
01605   FieldVector<ctype, dimworld> _normal;   
01606 };
01607 
01608 
01609 //========================================================================
01613 //========================================================================
01614 
01615 template<class GridImp>
01616 class YaspHierarchicIterator :
01617   public YaspEntityPointer<0,GridImp>,
01618   public HierarchicIteratorDefaultImplementation <GridImp,YaspHierarchicIterator>
01619 {
01620   enum { dim=GridImp::dimension };
01621   enum { dimworld=GridImp::dimensionworld };  
01622   typedef typename GridImp::ctype ctype;
01623 public:
01624   // types used from grids
01625   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01626   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01627   typedef typename GridImp::template Codim<0>::Entity Entity;
01628   typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01629 
01631   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01632 
01634   YaspHierarchicIterator (const YGLI& g, const TSI& it, int maxlevel) :
01635     YaspEntityPointer<0,GridImp>(g,it)
01636   {
01637         // now iterator points to current cell
01638 
01639         // determine maximum level
01640         _maxlevel = std::min(maxlevel,this->_g.mg()->maxlevel());
01641 
01642         // if maxlevel not reached then push yourself and sons
01643         if (this->_g.level()<_maxlevel)
01644           {
01645                 StackElem se(this->_g);
01646                 se.coord = this->_it.coord();
01647                 stack.push(se);
01648                 push_sons();
01649           }
01650 
01651         // and make iterator point to first son if stack is not empty
01652         if (!stack.empty())
01653           pop_tos();
01654   }
01655 
01657   YaspHierarchicIterator (const YaspHierarchicIterator& it) :
01658     YaspEntityPointer<0,GridImp>(it._g, it._it),
01659     _maxlevel(it._maxlevel), stack(it.stack)
01660   {}
01661 
01663   void increment ()
01664   {
01665         // sanity check: do nothing when stack is empty
01666         if (stack.empty()) return;
01667 
01668         // if maxlevel not reached then push sons
01669         if (this->_g.level()<_maxlevel)
01670           push_sons();
01671 
01672         // in any case pop one element
01673         pop_tos();
01674   }
01675 
01676   void print (std::ostream& s) const
01677   {
01678         s << "HIER: " << "level=" << this->_g.level()
01679           << " position=" << this->_it.coord()
01680           << " superindex=" << this->_it.superindex()
01681           << " maxlevel=" << this->_maxlevel
01682           << " stacksize=" << stack.size()
01683           << std::endl;
01684   }
01685 
01686 private:
01687   int _maxlevel;         
01688 
01689   struct StackElem {
01690         YGLI g;       // grid level of the element
01691         iTupel coord; // and the coordinates
01692         StackElem(YGLI gg) : g(gg) {}
01693   };
01694   Stack<StackElem> stack;      
01695 
01696   // push sons of current element on the stack
01697   void push_sons ()
01698   {
01699         // yes, process all 1<<dim sons
01700         StackElem se(this->_g.finer());
01701         for (int i=0; i<(1<<dim); i++)
01702           {
01703                 for (int k=0; k<dim; k++)
01704                   if (i&(1<<k))
01705                         se.coord[k] = this->_it.coord(k)*2+1;
01706                   else
01707                         se.coord[k] = this->_it.coord(k)*2;
01708                 stack.push(se);
01709           }
01710   }
01711 
01712   // make TOS the current element
01713   void pop_tos ()
01714   {
01715         StackElem se = stack.pop();
01716         this->_g = se.g;
01717         this->_it.reinit(this->_g.cell_overlap(),se.coord);
01718   }
01719 };
01720 
01721 //========================================================================
01731 //========================================================================
01732 template<int codim, class GridImp>
01733 class YaspEntityPointer :
01734   public EntityPointerDefaultImplementation<codim,GridImp,
01735                               Dune::YaspEntityPointer<codim,GridImp> >
01736 {
01738   enum { dim=GridImp::dimension };
01740   enum { dimworld=GridImp::dimensionworld };
01741   typedef typename GridImp::ctype ctype;
01742 public:
01743   typedef typename GridImp::template Codim<codim>::Entity Entity;
01744   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01745   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01746   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
01747   typedef YaspEntityPointer<codim,GridImp> Base;
01748 
01750   YaspEntityPointer (const YGLI & g, const TSI & it) : _g(g), _it(it), _entity(_g,_it)
01751   {
01752         if (codim>0 && codim<dim)
01753           {
01754                 DUNE_THROW(GridError, "YaspLevelIterator: codim not implemented");
01755           }
01756   }
01757 
01759   YaspEntityPointer (const YaspEntityPointer& rhs) : _g(rhs._g), _it(rhs._it), _entity(_g,_it)
01760   {
01761         if (codim>0 && codim<dim)
01762           {
01763                 DUNE_THROW(GridError, "YaspLevelIterator: codim not implemented");
01764           }
01765   }
01766 
01768   bool equals (const YaspEntityPointer& rhs) const
01769   {
01770         return (_it==rhs._it && _g == rhs._g);
01771   }
01772 
01774   Entity& dereference() const
01775   {
01776         return _entity;
01777   }
01778 
01780   int level () const {return _g.level();}
01781 
01782   const YaspEntityPointer&
01783   operator = (const YaspEntityPointer& rhs)
01784     {
01785       _g = rhs._g;
01786       _it = rhs._it;
01787       /* _entity = i._entity
01788        * is done implicitely, as the entity is completely
01789        * defined via the interator it belongs to
01790        */
01791       return *this;
01792     }  
01793   
01794 protected:
01795   YGLI _g;               // access to grid level
01796   TSI _it;               // position in the grid level
01797   mutable SpecialEntity _entity; 
01798 };
01799 
01800 //========================================================================
01808 //========================================================================
01809 
01810 template<int codim, PartitionIteratorType pitype, class GridImp>
01811 class YaspLevelIterator :
01812   public YaspEntityPointer<codim,GridImp>,
01813   public LevelIteratorDefaultImplementation<codim,pitype,GridImp,YaspLevelIterator>
01814 {
01816   enum { dim=GridImp::dimension };
01818   enum { dimworld=GridImp::dimensionworld };
01819   typedef typename GridImp::ctype ctype;
01820 public:
01821   typedef typename GridImp::template Codim<codim>::Entity Entity;
01822   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01823   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01824   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
01825 
01827   YaspLevelIterator (const YGLI & g, const TSI & it) :
01828     YaspEntityPointer<codim,GridImp>(g,it) {}
01829 
01831   YaspLevelIterator (const YaspLevelIterator& i) :
01832     YaspEntityPointer<codim,GridImp>(i) {}
01833 
01835   void increment()
01836   {
01837         ++(this->_it);
01838   }
01839 };
01840 
01841 
01842 //========================================================================
01847 //========================================================================
01848 
01849 template <class GridImp> 
01850 struct YaspLevelIndexSetTypes
01851 {
01853   template<int cd>
01854   struct Codim
01855   {
01856     template<PartitionIteratorType pitype>
01857     struct Partition
01858     {
01859       typedef typename GridImp::Traits::template Codim<cd>::template Partition<pitype>::LevelIterator Iterator;
01860     };
01861   };
01862 };
01863 
01864 template<class GridImp>
01865 class YaspLevelIndexSet : public IndexSetDefaultImplementation<GridImp,YaspLevelIndexSet<GridImp>,YaspLevelIndexSetTypes<GridImp> >
01866 {
01867   typedef IndexSetDefaultImplementation<GridImp,YaspLevelIndexSet<GridImp>,YaspLevelIndexSetTypes<GridImp> > Base;
01868 public:
01869 
01871   YaspLevelIndexSet (const GridImp& g, int l) : grid(g), level(l)
01872   {
01873     // contains a single element type;
01874     for (int codim=0; codim<=GridImp::dimension; codim++)
01875       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
01876   }
01877 
01879   template<int cd>
01880   int index (const typename GridImp::Traits::template Codim<cd>::Entity& e) const 
01881   {
01882     return grid.template getRealEntity<cd>(e).compressedIndex(); 
01883   }
01884 
01886   template<int cc>
01887   int subIndex (const typename GridImp::Traits::template Codim<0>::Entity& e, int i) const
01888   {
01889     return grid.template getRealEntity<0>(e).template subCompressedIndex<cc>(i);
01890   }
01891 
01893   int size (GeometryType type) const
01894   {
01895       return grid.size(level,GridImp::dimension-type.dim());
01896   }
01897 
01899   int size (int codim) const
01900   {
01901     return Base::size(codim);
01902   }
01903 
01905   const std::vector<GeometryType>& geomTypes (int codim) const
01906   {
01907     return mytypes[codim];
01908   }
01909 
01911   template<int cd, PartitionIteratorType pitype>
01912   typename Base::template Codim<cd>::template Partition<pitype>::Iterator begin () const
01913   {
01914     return grid.template lbegin<cd,pitype>(level);
01915   }
01916   
01918   template<int cd, PartitionIteratorType pitype>
01919   typename Base::template Codim<cd>::template Partition<pitype>::Iterator end () const
01920   {
01921     return grid.template lend<cd,pitype>(level);
01922   }
01923 
01924 private:
01925   const GridImp& grid;
01926   int level;
01927   std::vector<GeometryType> mytypes[GridImp::dimension+1];
01928 };
01929 
01930 
01931 // Leaf Index Set
01932 
01933 template <class GridImp> 
01934 struct YaspLeafIndexSetTypes
01935 {
01937   template<int cd>
01938   struct Codim
01939   {
01940     template<PartitionIteratorType pitype>
01941     struct Partition
01942     {
01943       typedef typename GridImp::Traits::template Codim<cd>::template Partition<pitype>::LevelIterator Iterator;
01944     };
01945   };
01946 };
01947 
01948 template<class GridImp>
01949 class YaspLeafIndexSet : public IndexSetDefaultImplementation<GridImp,YaspLeafIndexSet<GridImp>,YaspLeafIndexSetTypes<GridImp> >
01950 {
01951   typedef IndexSetDefaultImplementation<GridImp,YaspLeafIndexSet<GridImp>,YaspLeafIndexSetTypes<GridImp> > Base;
01952 public:
01953 
01955   YaspLeafIndexSet (const GridImp& g) : grid(g)
01956   {
01957     // contains a single element type;
01958     for (int codim=0; codim<=GridImp::dimension; codim++)
01959       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
01960   }
01961 
01963   /*
01964     We use the remove_const to extract the Type from the mutable class,
01965     because the const class is not instantiated yet.
01966   */
01967   template<int cd>
01968   int index (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const 
01969   {
01970     assert(cd==0 || cd==GridImp::dimension);
01971     return grid.template getRealEntity<cd>(e).compressedLeafIndex(); 
01972   }
01973 
01975   /*
01976     We use the remove_const to extract the Type from the mutable class,
01977     because the const class is not instantiated yet.
01978   */
01979   template<int cc>
01980   int subIndex (const typename remove_const<GridImp>::type::Traits::template Codim<0>::Entity& e, int i) const
01981   {
01982     return grid.template getRealEntity<0>(e).template subCompressedLeafIndex<cc>(i);
01983   }
01984 
01986   int size (GeometryType type) const
01987   {
01988       return grid.size(grid.maxLevel(),GridImp::dimension-type.dim());
01989   }
01990 
01992   int size (int codim) const
01993   {
01994     return Base::size(codim);
01995   }
01996 
01998   const std::vector<GeometryType>& geomTypes (int codim) const
01999   {
02000     return mytypes[codim];
02001   }
02002 
02004   template<int cd, PartitionIteratorType pitype>
02005   typename Base::template Codim<cd>::template Partition<pitype>::Iterator begin () const
02006   {
02007     return grid.template lbegin<cd,pitype>(grid.maxLevel());
02008   }
02009   
02011   template<int cd, PartitionIteratorType pitype>
02012   typename Base::template Codim<cd>::template Partition<pitype>::Iterator end () const
02013   {
02014     return grid.template lend<cd,pitype>(grid.maxLevel());
02015   }
02016 
02017 private:
02018   const GridImp& grid;
02019   enum { ncodim = remove_const<GridImp>::type::dimension+1 };
02020   std::vector<GeometryType> mytypes[ncodim];
02021 };
02022 
02023 
02024 
02025 
02026 //========================================================================
02031 //========================================================================
02032 
02033 template<class GridImp>
02034 class YaspGlobalIdSet : public IdSetDefaultImplementation<GridImp,YaspGlobalIdSet<GridImp>,
02035                                      typename remove_const<GridImp>::type::PersistentIndexType >
02036 /*
02037   We used the remove_const to extract the Type from the mutable class,
02038   because the const class is not instantiated yet.
02039 */
02040 {
02041 public:
02043   typedef typename remove_const<GridImp>::type::PersistentIndexType IdType;
02044 
02046   YaspGlobalIdSet (const GridImp& g) : grid(g) {}
02047 
02049   /*
02050     We use the remove_const to extract the Type from the mutable class,
02051     because the const class is not instantiated yet.
02052   */
02053   template<int cd>
02054   IdType id (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const 
02055   {
02056     return grid.template getRealEntity<cd>(e).persistentIndex();
02057   }
02058 
02060   /*
02061     We use the remove_const to extract the Type from the mutable class,
02062     because the const class is not instantiated yet.
02063   */
02064   template<int cc>
02065   IdType subId (const typename remove_const<GridImp>::type::Traits::template Codim<0>::Entity& e, int i) const
02066   {
02067     return grid.template getRealEntity<0>(e).template subPersistentIndex<cc>(i);
02068   }
02069 
02070 private:
02071   const GridImp& grid;
02072 };
02073 
02074 
02075 template<int dim, int dimworld>
02076 struct YaspGridFamily
02077 {
02078 #if HAVE_MPI
02079   typedef CollectiveCommunication<MPI_Comm> CCType;
02080 #else
02081   typedef CollectiveCommunication<Dune::YaspGrid<dim,dimworld> > CCType;
02082 #endif
02083 
02084   typedef GridTraits<dim,dimworld,Dune::YaspGrid<dim,dimworld>,
02085                      YaspGeometry,YaspEntity,
02086                      YaspEntityPointer,YaspLevelIterator,
02087                      YaspIntersectionIterator, // leaf  intersection iter 
02088                      YaspIntersectionIterator, // level intersection iter 
02089                      YaspHierarchicIterator,
02090                      YaspLevelIterator,
02091                      YaspLevelIndexSet<const YaspGrid<dim,dimworld> >,
02092                      YaspLevelIndexSetTypes<const YaspGrid<dim,dimworld> >,
02093                      YaspLeafIndexSet<const YaspGrid<dim,dimworld> >,
02094                      YaspLeafIndexSetTypes<const YaspGrid<dim,dimworld> >,
02095                      YaspGlobalIdSet<const YaspGrid<dim,dimworld> >, 
02096                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02097                      YaspGlobalIdSet<const YaspGrid<dim,dimworld> >, 
02098                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02099                      CCType> 
02100   Traits;  
02101 };
02102 
02103 template<int dim, int codim>
02104 struct YaspCommunicateMeta {
02105   template<class G, class DataHandle>
02106   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02107   {
02108     if (data.contains(dim,codim))
02109       g.template communicateCodim<DataHandle,codim>(data,iftype,dir,level);
02110     YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
02111   }
02112 };
02113 
02114 template<int dim>
02115 struct YaspCommunicateMeta<dim,0> {
02116   template<class G, class DataHandle>
02117   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02118   {
02119     if (data.contains(dim,0))
02120       g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
02121   }
02122 };
02123 
02124 
02125 //************************************************************************
02139 template<int dim, int dimworld>
02140 class YaspGrid :
02141   public GridDefaultImplementation<dim,dimworld,yaspgrid_ctype,YaspGridFamily<dim,dimworld> >,
02142   public MultiYGrid<dim,yaspgrid_ctype>
02143 {
02144   typedef const YaspGrid<dim,dimworld> GridImp;
02145 public:
02147   typedef yaspgrid_ctype ctype;
02148 
02149   // define the persistent index type
02150   typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits> PersistentIndexType;
02151 
02153   typedef YaspGridFamily<dim,dimworld> GridFamily;
02155   typedef typename YaspGridFamily<dim,dimworld>::Traits Traits;
02156 
02157   // need for friend declarations in entity
02158   typedef YaspLevelIndexSet<YaspGrid<dim,dimworld> > LevelIndexSetType;
02159   typedef YaspLeafIndexSet<YaspGrid<dim,dimworld> > LeafIndexSetType;
02160   typedef YaspGlobalIdSet<YaspGrid<dim,dimworld> > GlobalIdSetType;
02161 
02163   enum { MAXL=64 };
02164 
02166   typedef MultiYGrid<dim,ctype> YMG;
02167   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
02168   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
02169   typedef typename MultiYGrid<dim,ctype>::Intersection IS;
02170   typedef typename std::deque<IS>::const_iterator ISIT;
02171 
02173   std::string name() const { return "YaspGrid"; };
02174 
02182 #if HAVE_MPI
02183   YaspGrid (MPI_Comm comm, Dune::FieldVector<ctype, dim> L, 
02184             Dune::FieldVector<int, dim> s, 
02185             Dune::FieldVector<bool, dim> periodic, int overlap) 
02186     : MultiYGrid<dim,ctype>(comm,L,s,periodic,overlap), ccobj(comm)
02187 #else
02188   YaspGrid (Dune::FieldVector<ctype, dim> L, 
02189             Dune::FieldVector<int, dim> s, 
02190             Dune::FieldVector<bool, dim> periodic, int overlap) 
02191     : MultiYGrid<dim,ctype>(L,s,periodic,overlap)
02192 #endif
02193   {
02194     setsizes();
02195     indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dimworld> >(*this,0) );
02196     theleafindexset.push_back( new YaspLeafIndexSet<const YaspGrid<dim,dimworld> >(*this) );
02197     theglobalidset.push_back( new YaspGlobalIdSet<const YaspGrid<dim,dimworld> >(*this) );
02198   }
02199 
02200   ~YaspGrid()
02201   {
02202     deallocatePointers(indexsets);
02203     deallocatePointers(theleafindexset);
02204     deallocatePointers(theglobalidset);
02205   }
02206 
02207   private:
02208   // do not copy this class 
02209   YaspGrid(const YaspGrid&);
02210 
02211   public:
02212   
02216   int maxLevel() const {return MultiYGrid<dim,ctype>::maxlevel();} // delegate
02217 
02219   void globalRefine (int refCount)
02220   {
02221     bool b=true;
02222     for (int k=0; k<refCount; k++)
02223       {
02224         MultiYGrid<dim,ctype>::refine(b);
02225         setsizes();
02226         indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dimworld> >(*this,maxLevel()) );
02227       }
02228   }
02229 
02230   // set options for refinement
02231   void refineOptions (bool b)
02232   {
02233     keep_ovlp=b;
02234   }
02235   
02237   void refine (bool b)
02238   {
02239     MultiYGrid<dim,ctype>::refine(b);
02240     setsizes();
02241     indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim,dimworld> >(*this,maxLevel()) );
02242   }
02243 
02245   bool adapt ()    
02246   {
02247     globalRefine(1);
02248     return true; 
02249   }
02250   
02252   template<int cd, PartitionIteratorType pitype>
02253   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
02254     {
02255       return levelbegin<cd,pitype>(level);
02256     }
02257 
02259   template<int cd, PartitionIteratorType pitype>
02260   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
02261     {
02262       return levelend<cd,pitype>(level);
02263     }
02264 
02266   template<int cd>
02267   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
02268     {
02269       return levelbegin<cd,All_Partition>(level);
02270     }
02271 
02273   template<int cd>
02274   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
02275     {
02276       return levelend<cd,All_Partition>(level);
02277     }
02278 
02280   template<int cd, PartitionIteratorType pitype>
02281   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const
02282     {
02283       return levelbegin<cd,pitype>(maxLevel());
02284     }
02285 
02287   template<int cd, PartitionIteratorType pitype>
02288   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const
02289     {
02290       return levelend<cd,pitype>(maxLevel());
02291     }
02292 
02294   template<int cd>
02295   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
02296     {
02297       return levelbegin<cd,All_Partition>(maxLevel());
02298     };
02299 
02301   template<int cd>
02302   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
02303     {
02304       return levelend<cd,All_Partition>(maxLevel());
02305     }
02306 
02308   int overlapSize (int level, int codim) const
02309   {
02310         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02311         return g.overlap();
02312   }
02313 
02315   int overlapSize (int codim) const
02316   {
02317         YGLI g = MultiYGrid<dim,ctype>::begin(maxLevel());
02318         return g.overlap();
02319   }
02320 
02322   int ghostSize (int level, int codim) const
02323   {
02324         return 0;
02325   }
02326 
02328   int ghostSize (int codim) const
02329   {
02330         return 0;
02331   }
02332 
02334   int size (int level, int codim) const
02335   {
02336     return sizes[level][codim];
02337   }
02338 
02340   int size (int codim) const
02341   {
02342     return sizes[maxLevel()][codim];
02343   }
02344 
02346   int size (int level, GeometryType type) const
02347   {
02348       return (type.isCube()) ? sizes[level][dim-type.dim()] : 0;
02349   }
02350 
02352   int size (GeometryType type) const
02353   {
02354     return size(maxLevel(),type);
02355   }
02356 
02368 #if HAVE_MPI
02369 //   template<class T, template<class> class P, int codim>
02370 //   void communicate (T& t, InterfaceType iftype, CommunicationDirection dir, int level) const
02371 //   {
02372 //         IsTrue< ( codim == dim || codim == 0 ) >::yes();
02373 //         // access to grid level
02374 //         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02375 
02376 //         // find send/recv lists or throw error
02377 //         const std::deque<IS>* sendlist;
02378 //         const std::deque<IS>* recvlist;
02379 //         if (codim==0) // the elements
02380 //           {
02381 //                 if (iftype==InteriorBorder_InteriorBorder_Interface)
02382 //                   return; // there is nothing to do in this case
02383 //                 if (iftype==InteriorBorder_All_Interface)
02384 //                   {
02385 //                         sendlist = &g.send_cell_interior_overlap();
02386 //                         recvlist = &g.recv_cell_overlap_interior();
02387 //                   }
02388 //                 if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface || iftype==All_All_Interface)
02389 //                   {
02390 //                         sendlist = &g.send_cell_overlap_overlap();
02391 //                         recvlist = &g.recv_cell_overlap_overlap();
02392 //                   }
02393 //           }
02394 //         if (codim==dim) // the vertices
02395 //           {
02396 //                 if (iftype==InteriorBorder_InteriorBorder_Interface)
02397 //                   {
02398 //                         sendlist = &g.send_vertex_interiorborder_interiorborder();
02399 //                         recvlist = &g.recv_vertex_interiorborder_interiorborder();
02400 //                   }
02401 
02402 //                 if (iftype==InteriorBorder_All_Interface)
02403 //                   {
02404 //                         sendlist = &g.send_vertex_interiorborder_overlapfront();
02405 //                         recvlist = &g.recv_vertex_overlapfront_interiorborder();
02406 //                   }
02407 //                 if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface)
02408 //                   {
02409 //                         sendlist = &g.send_vertex_overlap_overlapfront();
02410 //                         recvlist = &g.recv_vertex_overlapfront_overlap();
02411 //                   }
02412 //                 if (iftype==All_All_Interface)
02413 //                   {
02414 //                         sendlist = &g.send_vertex_overlapfront_overlapfront();
02415 //                         recvlist = &g.recv_vertex_overlapfront_overlapfront();
02416 //                   }
02417 //           }
02418 //         if (codim>0 && codim<dim)
02419 //           {
02420 //                 DUNE_THROW(GridError, "interface communication not implemented");
02421 //           }
02422 
02423 //         // change communication direction?
02424 //         if (dir==BackwardCommunication)
02425 //           std::swap(sendlist,recvlist);
02426 
02427 //         // allocate & fill the send buffers & store send request
02428 //         std::vector<P<T>*> sends; // store pointers to send buffers
02429 //         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02430 //           {
02431 //                 // allocate send buffer
02432 //                 P<T> *buf = new P<T>[is->grid.totalsize()];
02433 
02434 //                 // remember send buffer
02435 //                 sends.push_back(buf);
02436 
02437 //                 // fill send buffer; iterate over cells in intersection
02438 //                 typename SubYGrid<dim,ctype>::SubIterator subend = is->grid.subend();
02439 //                 for (typename SubYGrid<dim,ctype>::SubIterator i=is->grid.subbegin(); i!=subend; ++i)
02440 //                   buf[i.index()].gather(t,i.superindex());
02441 
02442 //                 // hand over send request to torus class
02443 //                 MultiYGrid<dim,ctype>::torus().send(is->rank,buf,is->grid.totalsize()*sizeof(P<T>));
02444 //           }
02445 
02446 //         // allocate recv buffers and store receive request
02447 //         std::vector<P<T>*> recvs; // pointers to receive buffers
02448 //         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02449 //           {
02450 //                 // allocate recv buffer
02451 //                 P<T> *buf = new P<T>[is->grid.totalsize()];
02452 
02453 //                 // remember recv buffer
02454 //                 recvs.push_back(buf);
02455 
02456 //                 // hand over recv request to torus class
02457 //                 MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(P<T>));
02458 //           }
02459 
02460 //         // exchange all buffers now
02461 //         MultiYGrid<dim,ctype>::torus().exchange();
02462 
02463 //         // release send buffers
02464 //         for (int i=0; i<sends.size(); i++)
02465 //           delete[] sends[i];
02466 
02467 //         // process receive buffers and delete them
02468 //         int i=0;
02469 //         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02470 //           {
02471 //                 // get recv buffer
02472 //                 P<T> *buf = recvs[i++];
02473 
02474 //                 // copy data from receive buffer; iterate over cells in intersection
02475 //                 typename SubYGrid<dim,ctype>::SubIterator subend = is->grid.subend();
02476 //                 for (typename SubYGrid<dim,ctype>::SubIterator i=is->grid.subbegin(); i!=subend; ++i)
02477 //                   buf[i.index()].scatter(t,i.superindex());
02478 
02479 //                 // delete buffer
02480 //                 delete[] buf;
02481 //           }
02482 //   }
02483 #endif
02484 
02489   template<class DataHandleImp, class DataType>
02490   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir, int level) const
02491   {
02492     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
02493   }
02494 
02499   template<class DataHandleImp, class DataType>
02500   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir) const
02501   {
02502     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
02503   }
02504 
02509   template<class DataHandle, int codim>
02510   void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
02511   {
02512     // check input
02513     if (!data.contains(dim,codim)) return; // should have been checked outside
02514 
02515     // data types
02516     typedef typename DataHandle::DataType DataType;
02517 
02518     // access to grid level
02519     YGLI g = MultiYGrid<dim,ctype>::begin(level);
02520 
02521     // find send/recv lists or throw error
02522     const std::deque<IS>* sendlist=0;
02523     const std::deque<IS>* recvlist=0;
02524     if (codim==0) // the elements
02525       {
02526         if (iftype==InteriorBorder_InteriorBorder_Interface)
02527           return; // there is nothing to do in this case
02528         if (iftype==InteriorBorder_All_Interface)
02529           {
02530             sendlist = &g.send_cell_interior_overlap();
02531             recvlist = &g.recv_cell_overlap_interior();
02532           }
02533         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface || iftype==All_All_Interface)
02534           {
02535             sendlist = &g.send_cell_overlap_overlap();
02536             recvlist = &g.recv_cell_overlap_overlap();
02537           }
02538       }
02539     if (codim==dim) // the vertices
02540       {
02541         if (iftype==InteriorBorder_InteriorBorder_Interface)
02542           {
02543             sendlist = &g.send_vertex_interiorborder_interiorborder();
02544             recvlist = &g.recv_vertex_interiorborder_interiorborder();
02545           }
02546 
02547         if (iftype==InteriorBorder_All_Interface)
02548           {
02549             sendlist = &g.send_vertex_interiorborder_overlapfront();
02550             recvlist = &g.recv_vertex_overlapfront_interiorborder();
02551           }
02552         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface)
02553           {
02554             sendlist = &g.send_vertex_overlap_overlapfront();
02555             recvlist = &g.recv_vertex_overlapfront_overlap();
02556           }
02557         if (iftype==All_All_Interface)
02558           {
02559             sendlist = &g.send_vertex_overlapfront_overlapfront();
02560             recvlist = &g.recv_vertex_overlapfront_overlapfront();
02561           }
02562       }
02563     if (codim>0 && codim<dim)
02564       {
02565         DUNE_THROW(GridError, "interface communication not implemented");
02566       }
02567 
02568     // change communication direction?
02569     if (dir==BackwardCommunication)
02570       std::swap(sendlist,recvlist);
02571 
02572     // Size computation (requires communication if variable size)
02573     std::map<int,int> send_size;      // map rank to total number of objects (of type DataType) to be sent
02574     std::map<int,int> recv_size;      // map rank to total number of objects (of type DataType) to be recvd
02575     std::map<int,size_t*> send_sizes; // map rank to array giving number of objects per entity to be sent
02576     std::map<int,size_t*> recv_sizes; // map rank to array giving number of objects per entity to be recvd
02577     if (data.fixedsize(dim,codim))
02578       {
02579         // fixed size: just take a dummy entity, size can be computed without communication
02580         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02581           {
02582             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02583               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02584             send_size[is->rank] = is->grid.totalsize() * data.size(*it);
02585           }
02586         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02587           {
02588             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02589               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02590             recv_size[is->rank] = is->grid.totalsize() * data.size(*it);
02591           }
02592       }
02593     else
02594       {
02595         // variable size case: sender side determines the size
02596         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02597           {
02598             // allocate send buffer for sizes per entitiy
02599             size_t *buf = new size_t[is->grid.totalsize()];
02600             send_sizes[is->rank] = buf;
02601 
02602             // loop over entities and ask for size
02603             int i=0; size_t n=0;
02604             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02605               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02606             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02607               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02608             for ( ; it!=tsubend; ++it)
02609               {
02610                 buf[i] = data.size(*it);
02611                 n += buf[i];
02612                 i++;
02613               }
02614 
02615             // now we know the size for this rank
02616             send_size[is->rank] = n;
02617 
02618             // hand over send request to torus class
02619             MultiYGrid<dim,ctype>::torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
02620           }
02621 
02622         // allocate recv buffers for sizes and store receive request
02623         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02624           {
02625             // allocate recv buffer
02626             size_t *buf = new size_t[is->grid.totalsize()];
02627             recv_sizes[is->rank] = buf;
02628             
02629             // hand over recv request to torus class
02630             MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
02631           }
02632 
02633         // exchange all size buffers now
02634         MultiYGrid<dim,ctype>::torus().exchange();
02635 
02636         // release send size buffers
02637         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02638           delete[] send_sizes[is->rank];
02639 
02640         // process receive size buffers
02641         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02642           {
02643             // get recv buffer
02644             size_t *buf = recv_sizes[is->rank];
02645 
02646             // compute total size
02647             size_t n=0;
02648             for (int i=0; i<is->grid.totalsize(); ++i)
02649               n += buf[i];
02650 
02651             // ... and store it
02652             recv_size[is->rank] = n;
02653           }
02654       }
02655 
02656 
02657     // allocate & fill the send buffers & store send request
02658     std::map<int,DataType*> sends; // store pointers to send buffers
02659     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02660       {
02661 //      std::cout << "[" << this->comm().rank() << "] "
02662 //                << " send " << " dest=" << is->rank
02663 //                << " size=" << send_size[is->rank] 
02664 //                << std::endl;
02665 
02666         // allocate send buffer
02667         DataType *buf = new DataType[send_size[is->rank]];
02668 
02669         // remember send buffer
02670         sends[is->rank] = buf;
02671 
02672         // make a message buffer
02673         MessageBuffer<DataType> mb(buf);
02674 
02675         // fill send buffer; iterate over cells in intersection
02676         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02677           it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02678         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02679           tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02680         for ( ; it!=tsubend; ++it)
02681           data.gather(mb,*it);
02682 
02683         // hand over send request to torus class
02684         MultiYGrid<dim,ctype>::torus().send(is->rank,buf,send_size[is->rank]*sizeof(DataType));
02685       }
02686 
02687     // allocate recv buffers and store receive request
02688     std::map<int,DataType*> recvs; // store pointers to send buffers
02689     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02690       {
02691 //      std::cout << "[" << this->comm().rank() << "] "
02692 //                << " recv " << "  src=" << is->rank
02693 //                << " size=" << recv_size[is->rank] 
02694 //                << std::endl;
02695 
02696         // allocate recv buffer
02697         DataType *buf = new DataType[recv_size[is->rank]];
02698 
02699         // remember recv buffer
02700         recvs[is->rank] = buf;
02701 
02702         // hand over recv request to torus class
02703         MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,recv_size[is->rank]*sizeof(DataType));
02704       }
02705 
02706     // exchange all buffers now
02707     MultiYGrid<dim,ctype>::torus().exchange();
02708 
02709     // release send buffers
02710     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02711       delete[] sends[is->rank];
02712 
02713     // process receive buffers and delete them
02714     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02715       {
02716         // get recv buffer
02717         DataType *buf = recvs[is->rank];
02718 
02719         // make a message buffer
02720         MessageBuffer<DataType> mb(buf);
02721 
02722         // copy data from receive buffer; iterate over cells in intersection
02723         if (data.fixedsize(dim,codim))
02724           {
02725             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02726               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02727             size_t n=data.size(*it);
02728             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02729               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02730             for ( ; it!=tsubend; ++it)
02731               data.scatter(mb,*it,n);
02732           }
02733         else
02734           {
02735             int i=0;
02736             size_t *sbuf = recv_sizes[is->rank];
02737             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02738               it(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubbegin()));
02739             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02740               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(g,is->grid.tsubend()));
02741             for ( ; it!=tsubend; ++it)
02742               data.scatter(mb,*it,sbuf[i++]);
02743             delete[] sbuf;
02744           }
02745 
02746         // delete buffer
02747         delete[] buf;
02748       }
02749   }
02750 
02751   // The new index sets from DDM 11.07.2005
02752   const typename Traits::GlobalIdSet& globalIdSet() const
02753   {
02754     return *(theglobalidset[0]);
02755   }
02756   
02757   const typename Traits::LocalIdSet& localIdSet() const
02758   {
02759     return *(theglobalidset[0]);
02760   }
02761 
02762   const typename Traits::LevelIndexSet& levelIndexSet(int level) const
02763   {
02764     if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
02765     return *(indexsets[level]);
02766   }
02767 
02768   const typename Traits::LeafIndexSet& leafIndexSet() const
02769   {
02770     return *(theleafindexset[0]);
02771   }
02772 
02773 #if HAVE_MPI
02774 
02776   const CollectiveCommunication<MPI_Comm>& comm () const
02777   {
02778     return ccobj;
02779   }
02780 #else
02781 
02783   const CollectiveCommunication<YaspGrid>& comm () const
02784   {
02785     return ccobj;
02786   }
02787 #endif
02788 
02789   YaspIntersectionIterator<const YaspGrid<dim, dimworld> >&
02790   getRealIntersectionIterator(typename Traits::LevelIntersectionIterator& it)
02791   {
02792     return this->getRealImplementation(it);
02793   }
02794 
02795   const YaspIntersectionIterator<const YaspGrid<dim, dimworld> >&
02796   getRealIntersectionIterator(const typename Traits::LevelIntersectionIterator& it) const
02797   {
02798     return this->getRealImplementation(it);
02799   }
02800 
02801 private:
02802 
02803 #if HAVE_MPI
02804   CollectiveCommunication<MPI_Comm> ccobj; 
02805 #else
02806   CollectiveCommunication<YaspGrid> ccobj; 
02807 #endif
02808 
02809   std::vector<YaspLevelIndexSet<const YaspGrid<dim,dimworld> >*> indexsets;
02810   std::vector<YaspLeafIndexSet<const YaspGrid<dim,dimworld> >*> theleafindexset;
02811   std::vector<YaspGlobalIdSet<const YaspGrid<dim,dimworld> >*> theglobalidset;
02812 
02813   // Index classes need access to the real entity
02814   friend class Dune::YaspLevelIndexSet<const Dune::YaspGrid<dim,dimworld> >;
02815   friend class Dune::YaspLeafIndexSet<const Dune::YaspGrid<dim,dimworld> >;
02816   friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim,dimworld> >;
02817 
02818   template<class T>
02819   void deallocatePointers(T& container)
02820   {
02821     typedef typename T::iterator Iterator;
02822     
02823     for(Iterator entry=container.begin(); entry != container.end(); ++entry)
02824       delete (*entry);
02825   }
02826   
02827   template<int codim>
02828   YaspEntity<codim,dim,const YaspGrid<dim,dimworld> >& 
02829   getRealEntity(typename Traits::template Codim<codim>::Entity& e )
02830   {
02831     return this->getRealImplementation(e);
02832   }
02833   
02834   template<int codim>
02835   const YaspEntity<codim,dim,const YaspGrid<dim,dimworld> >& 
02836   getRealEntity(const typename Traits::template Codim<codim>::Entity& e ) const
02837   {
02838     return this->getRealImplementation(e);
02839   }
02840 
02841   template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
02842   friend class Entity;
02843 
02844   template<class DT>
02845   class MessageBuffer {
02846   public:
02847     // Constructor
02848     MessageBuffer (DT *p)
02849     {
02850       a=p;
02851       i=0;
02852       j=0;
02853     }
02854 
02855     // write data to message buffer, acts like a stream !
02856     template<class Y>
02857     void write (const Y& data)
02858     {
02859       IsTrue< ( is_same<DT,Y>::value ) >::yes();
02860       a[i++] = data;
02861     }
02862     
02863     // read data from message buffer, acts like a stream !
02864     template<class Y>
02865     void read (Y& data) const 
02866     {
02867       IsTrue< ( is_same<DT,Y>::value ) >::yes();
02868       data = a[j++];
02869     }
02870 
02871   private:
02872     DT *a;
02873       int i;
02874     mutable int j;
02875   };
02876 
02877   void setsizes ()
02878   {
02879     for (YGLI g=MultiYGrid<dim,ctype>::begin(); g!=MultiYGrid<dim,ctype>::end(); ++g)
02880       {
02881         // codim 0 (elements)
02882         sizes[g.level()][0] = 1;
02883         for (int i=0; i<dim; ++i) 
02884           sizes[g.level()][0] *= g.cell_overlap().size(i);
02885 
02886         // codim 1 (faces)
02887         if (dim>1)
02888           {
02889             sizes[g.level()][1] = 0;
02890             for (int i=0; i<dim; ++i) 
02891               {
02892                 int s=g.cell_overlap().size(i)+1;
02893                 for (int j=0; j<dim; ++j)
02894                   if (j!=i)
02895                     s *= g.cell_overlap().size(j);
02896                 sizes[g.level()][1] += s;
02897               }
02898           }
02899 
02900         // codim dim-1 (edges)
02901         if (dim>2)
02902           {
02903             sizes[g.level()][dim-1] = 0;
02904             for (int i=0; i<dim; ++i) 
02905               {
02906                 int s=g.cell_overlap().size(i);
02907                 for (int j=0; j<dim; ++j)
02908                   if (j!=i)
02909                     s *= g.cell_overlap().size(j)+1;
02910                 sizes[g.level()][dim-1] += s;
02911               }
02912           }
02913 
02914         // codim dim (vertices)
02915         sizes[g.level()][dim] = 1;
02916         for (int i=0; i<dim; ++i) 
02917           sizes[g.level()][dim] *= g.vertex_overlapfront().size(i);
02918       }
02919   }
02920 
02922   template<int cd, PartitionIteratorType pitype>
02923   YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
02924   {
02925         IsTrue< ( cd == dim || cd == 0 ) >::yes();
02926         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02927         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
02928         if (cd==0) // the elements
02929           {
02930                 if (pitype<=InteriorBorder_Partition) 
02931                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_interior().tsubbegin());
02932                 if (pitype<=All_Partition)            
02933                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_overlap().tsubbegin());
02934           }
02935         if (cd==dim) // the vertices
02936           {
02937                 if (pitype==Interior_Partition)       
02938                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interior().tsubbegin());
02939                 if (pitype==InteriorBorder_Partition) 
02940                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interiorborder().tsubbegin());
02941                 if (pitype==Overlap_Partition)        
02942                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlap().tsubbegin());
02943                 if (pitype<=All_Partition)            
02944                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlapfront().tsubbegin());
02945           }
02946         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
02947   }
02948 
02950   template<int cd, PartitionIteratorType pitype>
02951   YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
02952   {
02953         IsTrue< ( cd == dim || cd == 0 ) >::yes();
02954         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02955         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
02956         if (cd==0) // the elements
02957           {
02958                 if (pitype<=InteriorBorder_Partition) 
02959                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_interior().tsubend());
02960                 if (pitype<=All_Partition)            
02961                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.cell_overlap().tsubend());
02962           }
02963         if (cd==dim) // the vertices
02964           {
02965                 if (pitype==Interior_Partition)       
02966                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interior().tsubend());
02967                 if (pitype==InteriorBorder_Partition) 
02968                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_interiorborder().tsubend());
02969                 if (pitype==Overlap_Partition)        
02970                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlap().tsubend());
02971                 if (pitype<=All_Partition)            
02972                   return YaspLevelIterator<cd,pitype,GridImp>(g,g.vertex_overlapfront().tsubend());
02973           }
02974         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
02975   }
02976 
02977   int sizes[MAXL][dim+1]; // total number of entities per level and codim
02978   bool keep_ovlp;
02979 };
02980 
02981 namespace Capabilities
02982 {
02983 
02995   template<int dim, int dimw>
02996   struct hasEntity< YaspGrid<dim,dimw>, 0 >
02997   {
02998     static const bool v = true;
02999   };
03000   
03004   template<int dim, int dimw>
03005   struct hasEntity< YaspGrid<dim,dimw>, dim >
03006   {
03007     static const bool v = true;
03008   };
03009 
03013   template<int dim, int dimw>
03014   struct isParallel< YaspGrid<dim,dimw> >
03015   {
03016     static const bool v = true;
03017   };
03018 
03022   template<int dim, int dimw>
03023   struct isLevelwiseConforming< YaspGrid<dim,dimw> >
03024   {
03025     static const bool v = true;
03026   };
03027 
03031   template<int dim, int dimw>
03032   struct isLeafwiseConforming< YaspGrid<dim,dimw> >
03033   {
03034     static const bool v = true;
03035   };
03036 
03040   template<int dim, int dimw>
03041   struct IsUnstructured< YaspGrid<dim,dimw> >
03042   {
03043     static const bool v = false;
03044   };
03045 
03049   template<int dim, int dimw>
03050   struct hasHangingNodes< YaspGrid<dim,dimw> >
03051   {
03052     static const bool v = false;
03053   };
03054 
03055 }
03056 
03057 } // end namespace
03058 
03059 
03060 #endif
03061 

Generated on 9 Apr 2008 with Doxygen (ver 1.5.2) [logfile].