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

Generated on 12 Dec 2007 with Doxygen (ver 1.5.1)