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<algorithm>
00007 #include<stack>
00008 
00009 // either include stdint.h or provide fallback for uint8_t
00010 #if HAVE_STDINT_H
00011 #include<stdint.h>
00012 #else
00013 typedef unsigned char uint8_t;
00014 #endif
00015 
00016 #include <dune/grid/common/grid.hh>     // the grid base classes
00017 #include <dune/grid/yaspgrid/grids.hh>  // the yaspgrid base classes
00018 #include <dune/grid/common/capabilities.hh> // the capabilities
00019 #include <dune/common/misc.hh>
00020 #include <dune/common/bigunsignedint.hh>
00021 #include <dune/common/typetraits.hh>
00022 #include <dune/common/collectivecommunication.hh>
00023 #include <dune/common/mpihelper.hh>
00024 #include <dune/grid/common/indexidset.hh>
00025 #include <dune/grid/common/datahandleif.hh>
00026 #include <dune/grid/genericgeometry/topologytypes.hh>
00027 
00028 
00029 #if HAVE_MPI
00030 #include <dune/common/mpicollectivecommunication.hh>
00031 #endif
00032 
00040 namespace Dune {
00041 
00042 //************************************************************************
00046 typedef double yaspgrid_ctype;
00047 static const yaspgrid_ctype yasptolerance=1E-13; // tolerance in coordinate computations
00048 
00049 /* some sizes for building global ids
00050  */
00051   const int yaspgrid_dim_bits = 24;   // bits for encoding each dimension
00052   const int yaspgrid_level_bits = 6;  // bits for encoding level number
00053   const int yaspgrid_codim_bits = 4;  // bits for encoding codimension
00054 
00055 
00056 //************************************************************************
00057 // forward declaration of templates
00058 
00059 template<int dim>                             class YaspGrid;
00060 template<int mydim, int cdim, class GridImp>  class YaspGeometry;
00061 template<int codim, int dim, class GridImp>   class YaspEntity;
00062 template<int codim, class GridImp>            class YaspEntityPointer;
00063 template<int codim, PartitionIteratorType pitype, class GridImp> class YaspLevelIterator;
00064 template<class GridImp>            class YaspIntersectionIterator;
00065 template<class GridImp>            class YaspIntersection;
00066 template<class GridImp>            class YaspHierarchicIterator;
00067 template<class GridImp>            class YaspLevelIndexSet;
00068 template<class GridImp>            class YaspLeafIndexSet;
00069 template<class GridImp>            class YaspGlobalIdSet;
00070 
00071 //========================================================================
00079 //========================================================================
00080 
00081 template<int mydim, int cdim, class GridImp>
00082 class YaspSpecialGeometry : public Geometry<mydim, cdim, GridImp, YaspGeometry>
00083 {
00085   typedef typename GridImp::ctype ctype;
00086 public:
00087   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, uint8_t& m) :
00088     Geometry<mydim, cdim, GridImp, YaspGeometry>(YaspGeometry<mydim, cdim, GridImp>(p,h,m))
00089     {}
00090   YaspSpecialGeometry() :
00091     Geometry<mydim, cdim, GridImp, YaspGeometry>(YaspGeometry<mydim, cdim, GridImp>(false))
00092     {};
00093 };
00094 
00095 template<int mydim, class GridImp>
00096 class YaspSpecialGeometry<mydim,mydim,GridImp> : public Geometry<mydim, mydim, GridImp, YaspGeometry>
00097 {
00099   typedef typename GridImp::ctype ctype;
00100 public:
00101   YaspSpecialGeometry(const FieldVector<ctype, mydim>& p, const FieldVector<ctype, mydim>& h) :
00102     Geometry<mydim, mydim, GridImp, YaspGeometry>(YaspGeometry<mydim, mydim, GridImp>(p,h))
00103     {}
00104   YaspSpecialGeometry() :
00105     Geometry<mydim, mydim, GridImp, YaspGeometry>(YaspGeometry<mydim, mydim, GridImp>(false))
00106     {};
00107 };
00108 
00109 template<int cdim, class GridImp>
00110 class YaspSpecialGeometry<0,cdim,GridImp> : public Geometry<0, cdim, GridImp, YaspGeometry>
00111 {
00113   typedef typename GridImp::ctype ctype;
00114 public:
00115   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p) :
00116     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(p))
00117     {}
00118   YaspSpecialGeometry(const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, uint8_t& m) :
00119     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(p))
00120     {}
00121   YaspSpecialGeometry() :
00122     Geometry<0, cdim, GridImp, YaspGeometry>(YaspGeometry<0, cdim, GridImp>(false))
00123     {};
00124 };
00125 
00126 //========================================================================
00127 // The transformation describing the refinement rule
00128 
00129 template<int dim, class GridImp>
00130 class YaspFatherRelativeLocalElement {
00131 public:
00132   static FieldVector<yaspgrid_ctype, dim> midpoint;  // data neded for the refelem below
00133   static FieldVector<yaspgrid_ctype, dim> extension; // data needed for the refelem below
00134   static YaspSpecialGeometry<dim,dim,GridImp> geo;
00135   static YaspSpecialGeometry<dim,dim,GridImp>& getson (int i)
00136   {
00137     for (int k=0; k<dim; k++)
00138       if (i&(1<<k))
00139         midpoint[k] = 0.75;
00140       else
00141         midpoint[k] = 0.25;
00142     return geo;
00143   }
00144 };
00145 
00146 // initialize static variable with bool constructor (which makes reference elements)
00147 template<int dim, class GridImp>
00148 YaspSpecialGeometry<dim,dim,GridImp> 
00149 YaspFatherRelativeLocalElement<dim,GridImp>::geo(YaspFatherRelativeLocalElement<dim,GridImp>::midpoint,
00150                                                  YaspFatherRelativeLocalElement<dim,GridImp>::extension);
00151 template<int dim, class GridImp>
00152 FieldVector<yaspgrid_ctype,dim> YaspFatherRelativeLocalElement<dim,GridImp>::midpoint(0.25);
00153 
00154 template<int dim, class GridImp>
00155 FieldVector<yaspgrid_ctype,dim> YaspFatherRelativeLocalElement<dim,GridImp>::extension(0.5);
00156 
00158 template<int mydim,int cdim, class GridImp>
00159 class YaspGeometry : public GeometryDefaultImplementation<mydim,cdim,GridImp,YaspGeometry>
00160 {
00161 public:
00163   typedef typename GridImp::ctype ctype;
00165   typedef Geometry<mydim,mydim,GridImp,Dune::YaspGeometry> ReferenceGeometry;
00166   
00168   GeometryType type () const
00169   {
00170       return GeometryType(GeometryType::cube,mydim);
00171   }
00172 
00174   bool affine() const { return true; }
00175 
00177   int corners () const
00178   {
00179         return 1<<mydim;
00180   }
00181 
00183   const FieldVector<ctype, cdim>& operator[] (int i) const DUNE_DEPRECATED
00184   {
00185       return corner(i);
00186   }
00187 
00189   FieldVector< ctype, cdim > corner ( const int i ) const
00190   {
00191       assert( i >= 0 && i < (int) coord_.N() ); 
00192       FieldVector<ctype, cdim>& c = coord_[i];
00193       int bit=0;
00194       for (int k=0; k<cdim; k++) // run over all directions in world
00195       {
00196           if (k==missing)
00197           {
00198               c[k] = midpoint[k];
00199               continue;
00200           }
00201           //k is not the missing direction
00202           if (i&(1<<bit))  // check whether bit is set or not
00203               c[k] = midpoint[k]+0.5*extension[k]; // bit is 1 in i
00204           else
00205               c[k] = midpoint[k]-0.5*extension[k]; // bit is 0 in i
00206           bit++; // we have processed a direction
00207       }
00208       
00209       return c;
00210   }
00211 
00213   FieldVector< ctype, cdim > center ( ) const
00214   {
00215       return midpoint;
00216   }
00217 
00219   FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
00220   {
00221       FieldVector<ctype, cdim> g;
00222       int bit=0;
00223       for (int k=0; k<cdim; k++)
00224           if (k==missing)
00225               g[k] = midpoint[k];
00226           else
00227           {
00228               g[k] = midpoint[k] + (local[bit]-0.5)*extension[k];
00229               bit++;
00230           }
00231       return g;
00232   }
00233 
00235   FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& global) const
00236   {
00237       FieldVector<ctype, mydim> l; // result
00238       int bit=0;
00239       for (int k=0; k<cdim; k++)
00240           if (k!=missing)
00241           {
00242               l[bit] = (global[k]-midpoint[k])/extension[k] + 0.5;
00243               bit++;
00244           }
00245       return l;
00246   }
00247 
00249   ctype volume () const
00250   {
00251     ctype volume=1.0;
00252     for (int k=0; k<cdim; k++)
00253       if (k!=missing) volume *= extension[k];
00254     return volume;
00255   }
00256 
00259   ctype integrationElement (const FieldVector<ctype, mydim>& local) const
00260   {
00261     return volume();
00262   }
00263 
00265   FieldMatrix<ctype,mydim,cdim>& jacobianTransposed (const FieldVector<ctype, mydim>& local) const
00266   {
00267       JT = 0.0;
00268       int k=0;
00269       for (int i=0; i<cdim; ++i)
00270       {
00271           if (i != missing)
00272           {
00273               JT[k][i] = extension[i]; // set diagonal element
00274               k++;
00275           }
00276       }
00277       return JT;
00278   }
00280   FieldMatrix<ctype,cdim,mydim>& jacobianInverseTransposed (const FieldVector<ctype, mydim>& local) const
00281   {
00282       Jinv = 0.0;
00283       int k=0;
00284       for (int i=0; i<cdim; ++i)
00285       {
00286           if (i != missing)
00287           {
00288               Jinv[i][k] = 1.0/extension[i]; // set diagonal element
00289               k++;
00290           }
00291       }
00292       return Jinv;
00293   }
00294 
00295 #if 0
00296 
00297   bool checkInside (const FieldVector<ctype, mydim>& local) const
00298   {
00299         for (int i=0; i<mydim; i++)
00300           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00301         return true;
00302   }
00303 #endif
00304 
00306   YaspGeometry (const FieldVector<ctype, cdim>& p, const FieldVector<ctype, cdim>& h, uint8_t& m)
00307         : midpoint(p), extension(h), missing(m)
00308   {
00309         if (cdim!=mydim+1)
00310           DUNE_THROW(GridError, "general YaspGeometry assumes cdim=mydim+1");
00311   }
00312 
00314   YaspGeometry (const YaspGeometry& other) 
00315         : midpoint(other.midpoint), 
00316           extension(other.extension), 
00317           missing(other.missing)
00318   {
00319   }
00320 
00322   void print (std::ostream& s) const
00323   {
00324         s << "YaspGeometry<"<<mydim<<","<<cdim<< "> ";
00325         s << "midpoint";
00326         for (int i=0; i<cdim; i++)
00327                   s << " " << midpoint[i];
00328                 s << " extension";
00329         for (int i=0; i<cdim; i++)
00330                   s << " " << extension[i];
00331                 s << " missing is " << missing;
00332   }
00333 private:
00334   // the element is fully defined by its midpoint the extension
00335   // in each direction and the missing direction.
00336   // References are used because this information
00337   // is known outside the element in many cases.
00338   // Note cdim==cdim+1
00339 
00340   // IMPORTANT midpoint and extension are references,
00341   // YaspGeometry can't be copied
00342   const FieldVector<ctype, cdim> & midpoint;    // the midpoint
00343   const FieldVector<ctype, cdim> & extension;   // the extension
00344   uint8_t& missing;                       // the missing, i.e. constant direction
00345 
00346   // In addition we need memory in order to return references.
00347   // Possibly we should change this in the interface ...
00348   mutable FieldMatrix<ctype, mydim, cdim> JT;     // the transposed of the jacobian 
00349   mutable FieldMatrix<ctype, cdim, mydim> Jinv;   // the transposed of the jacobian inverse 
00350   mutable FieldMatrix<ctype, Power_m_p<2,mydim>::power, cdim> coord_;   // the coordinates 
00351 
00352   const YaspGeometry<mydim,cdim,GridImp>&
00353   operator = (const YaspGeometry<mydim,cdim,GridImp>& g);
00354   
00355 };
00356 
00357 
00358 
00360 template<int mydim, class GridImp>
00361 class YaspGeometry<mydim,mydim,GridImp> : public GeometryDefaultImplementation<mydim,mydim,GridImp,YaspGeometry>
00362 {
00363 public:
00364   typedef typename GridImp::ctype ctype;
00366   typedef Geometry<mydim,mydim,GridImp,Dune::YaspGeometry> ReferenceGeometry;
00367   
00369   GeometryType type () const
00370   {
00371       return GeometryType(GeometryType::cube,mydim);
00372   }
00373 
00375   bool affine() const { return true; }
00376 
00378   int corners () const
00379   {
00380         return 1<<mydim;
00381   }
00382 
00384   const FieldVector<ctype, mydim>& operator[] (int i) const
00385   {
00386       return corner(i);
00387   }
00388 
00390   FieldVector< ctype, mydim > corner ( const int i ) const
00391   {
00392       assert( i >= 0 && i < (int) coord_.N() ); 
00393       FieldVector<ctype, mydim>& c = coord_[i];
00394       for (int k=0; k<mydim; k++)
00395           if (i&(1<<k))
00396               c[k] = midpoint[k]+0.5*extension[k]; // kth bit is 1 in i
00397           else
00398               c[k] = midpoint[k]-0.5*extension[k]; // kth bit is 0 in i
00399       return c;
00400   }
00401 
00403   FieldVector< ctype, mydim > center ( ) const
00404   {
00405       return midpoint;
00406   }
00407 
00409   FieldVector<ctype, mydim> global (const FieldVector<ctype, mydim>& local) const
00410   {
00411       FieldVector<ctype,mydim> g;
00412         for (int k=0; k<mydim; k++)
00413           g[k] = midpoint[k] + (local[k]-0.5)*extension[k];
00414         return g;
00415   }
00416 
00418   FieldVector<ctype, mydim> local (const FieldVector<ctype,mydim>& global) const
00419   {
00420       FieldVector<ctype, mydim> l; // result
00421         for (int k=0; k<mydim; k++)
00422           l[k] = (global[k]-midpoint[k])/extension[k] + 0.5;
00423         return l;
00424   }
00425 
00428   ctype integrationElement (const FieldVector<ctype, mydim>& local) const
00429   {
00430     return volume();
00431   }
00432 
00434   ctype volume () const
00435   {
00436     ctype vol=1.0;
00437     for (int k=0; k<mydim; k++) vol *= extension[k];
00438     return vol;
00439   }
00440 
00442   FieldMatrix<ctype,mydim,mydim>& jacobianTransposed (const FieldVector<ctype, mydim>& local) const
00443   {
00444         for (int i=0; i<mydim; ++i)
00445           {
00446                 JT[i] = 0.0;                // set column to zero
00447                 JT[i][i] = extension[i]; // set diagonal element
00448           }
00449         return JT;
00450   }
00452   FieldMatrix<ctype,mydim,mydim>& jacobianInverseTransposed (const FieldVector<ctype, mydim>& local) const
00453   {
00454         for (int i=0; i<mydim; ++i)
00455           {
00456                 Jinv[i] = 0.0;                // set column to zero
00457                 Jinv[i][i] = 1.0/extension[i]; // set diagonal element
00458           }
00459         return Jinv;
00460   }
00461 
00462 #if 0
00463 
00464   bool checkInside (const FieldVector<ctype, mydim>& local) const
00465   {
00466         for (int i=0; i<mydim; i++)
00467           if (local[i]<-yasptolerance || local[i]>1+yasptolerance) return false;
00468         return true;
00469   }
00470 #endif
00471 
00473   YaspGeometry (const FieldVector<ctype, mydim>& p, const FieldVector<ctype, mydim>& h)
00474         : midpoint(p), extension(h)
00475   {}
00476 
00478   YaspGeometry (const YaspGeometry& other) 
00479         : midpoint(other.midpoint), 
00480           extension(other.extension) 
00481   {
00482   }
00483 
00485   void print (std::ostream& s) const
00486   {
00487         s << "YaspGeometry<"<<mydim<<","<<mydim<< "> ";
00488         s << "midpoint";
00489         for (int i=0; i<mydim; i++)
00490                   s << " " << midpoint[i];
00491                 s << " extension";
00492         for (int i=0; i<mydim; i++)
00493                   s << " " << extension[i];
00494   }
00495 
00496 private:
00497   // the element is fully defined by midpoint and the extension
00498   // in each direction. References are used because this information
00499   // is known outside the element in many cases.
00500   // Note mydim==cdim
00501 
00502   // IMPORTANT midpoint and extension are references,
00503   // YaspGeometry can't be copied
00504   const FieldVector<ctype, mydim> & midpoint;   // the midpoint
00505   const FieldVector<ctype, mydim> & extension;  // the extension
00506 
00507   // In addition we need memory in order to return references.
00508   // Possibly we should change this in the interface ...
00509   mutable FieldMatrix<ctype, mydim, mydim> Jinv,JT;   // the transpose of the jacobian and its inverse inverse
00510   mutable FieldMatrix<ctype, Power_m_p<2,mydim>::power, mydim> coord_;   // the coordinates 
00511 
00512   // disable copy
00513   const YaspGeometry<mydim,mydim,GridImp>&
00514   operator = (const YaspGeometry<mydim,mydim,GridImp>& g);
00515 };
00516 
00518 template<int cdim, class GridImp>
00519 class YaspGeometry<0,cdim,GridImp> : public GeometryDefaultImplementation<0,cdim,GridImp,YaspGeometry>
00520 {
00521 public:
00522   typedef typename GridImp::ctype ctype;
00523   
00525   GeometryType type () const
00526   {
00527       return GeometryType(GeometryType::cube,0);
00528   }
00529 
00531   bool affine() const { return true; }
00532 
00534   int corners () const
00535   {
00536         return 1;
00537   }
00538 
00540   const FieldVector<ctype, cdim>& operator[] (int i) const
00541   {
00542       return position;
00543   }
00544 
00546   FieldVector< ctype, cdim > corner ( const int i ) const
00547   {
00548       return position;
00549   }
00550 
00552   FieldVector< ctype, cdim > center ( ) const
00553   {
00554       return position;
00555   }
00556 
00559   ctype integrationElement (const FieldVector<ctype, 0>& local) const
00560   {
00561         return 1.0;
00562   }
00563 
00565   FieldMatrix<ctype,0,cdim>& jacobianTransposed (const FieldVector<ctype, 0>& local) const
00566   {
00567       static FieldMatrix<ctype,0,cdim> JT(0.0);
00568       return JT;
00569   }
00571   FieldMatrix<ctype,cdim,0>& jacobianInverseTransposed (const FieldVector<ctype, 0>& local) const
00572   {
00573       static FieldMatrix<ctype,cdim,0> Jinv(0.0);
00574       return Jinv;
00575   }
00576 
00578   YaspGeometry (const FieldVector<ctype, cdim>& p) : position(p)
00579   {}
00580 
00582   void print (std::ostream& s) const
00583   {
00584         s << "YaspGeometry<"<<0<<","<<cdim<< "> ";
00585         s << "position " << position;
00586   }
00587 
00588 private:
00589   // IMPORTANT position is a reference,
00590   // YaspGeometry can't be copied
00591   const FieldVector<ctype, cdim> & position; 
00592 
00593   const YaspGeometry<0,cdim,GridImp>&
00594   operator = (const YaspGeometry<0,cdim,GridImp>& g);
00595 };
00596 
00597 // operator<< for all YaspGeometrys
00598 template <int mydim, int cdim, class GridImp>
00599 inline
00600 std::ostream& operator<< (std::ostream& s, YaspGeometry<mydim,cdim,GridImp>& e)
00601 {
00602   e.print(s);
00603   return s;
00604 }
00605 
00606 //========================================================================
00614 //========================================================================
00615 
00616 template<int codim, int dim, class GridImp>
00617 class YaspSpecialEntity :
00618   public GridImp::template Codim<codim>::Entity
00619 {
00620 public:
00621   typedef typename GridImp::ctype ctype;
00622   
00623   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00624   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00625 
00626   YaspSpecialEntity(const GridImp* yg, const YGLI& g, const TSI& it) :
00627     GridImp::template Codim<codim>::Entity (YaspEntity<codim, dim, GridImp>(yg,g,it))
00628     {};
00629   YaspSpecialEntity(const YaspEntity<codim, dim, GridImp>& e) :
00630     GridImp::template Codim<codim>::Entity (e)
00631     {};
00632   const TSI& transformingsubiterator () const
00633   {
00634         return this->realEntity.transformingsubiterator();
00635   }
00636   const YGLI& gridlevel () const
00637   {
00638         return this->realEntity.gridlevel();
00639   }
00640   const GridImp * yaspgrid() const
00641   {
00642       return this->realEntity.yaspgrid();
00643   }
00644 };
00645 
00646 template<int codim, int dim, class GridImp>
00647 class YaspEntity 
00648 :  public EntityDefaultImplementation <codim,dim,GridImp,YaspEntity>
00649 {
00650 public:
00651   typedef typename GridImp::ctype ctype;
00652   
00653   typedef typename GridImp::template Codim<codim>::Geometry Geometry;
00655   int level () const
00656   {
00657         DUNE_THROW(GridError, "YaspEntity not implemented");
00658   }
00659 
00661   int index () const
00662   {
00663         DUNE_THROW(GridError, "YaspEntity not implemented");
00664   }
00665 
00667   const Geometry& geometry () const
00668   {
00669         DUNE_THROW(GridError, "YaspEntity not implemented");
00670   }
00671 
00673   PartitionType partitionType () const
00674   {
00675         DUNE_THROW(GridError, "YaspEntity not implemented");
00676   }
00677 
00678   const GridImp * yaspgrid() const
00679   {
00680         DUNE_THROW(GridError, "YaspEntity not implemented");
00681   }
00682   
00683   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00684   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00685   YaspEntity (const GridImp* yg, const YGLI& g, const TSI& it)
00686   {
00687         DUNE_THROW(GridError, "YaspEntity not implemented");
00688   }
00689 
00690   // IndexSets needs access to the private index methods
00691   friend class Dune::YaspLevelIndexSet<GridImp>;
00692   friend class Dune::YaspLeafIndexSet<GridImp>;
00693   friend class Dune::YaspGlobalIdSet<GridImp>;
00694   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00695 
00697   PersistentIndexType persistentIndex () const 
00698   { 
00699     DUNE_THROW(GridError, "YaspEntity not implemented");
00700   }
00701 
00703   int compressedIndex () const 
00704   {
00705     DUNE_THROW(GridError, "YaspEntity not implemented");
00706   }
00707 
00709   int compressedLeafIndex () const 
00710   {
00711     DUNE_THROW(GridError, "YaspEntity not implemented");
00712   }
00713 };
00714 
00715 
00716 // specialization for codim=0
00717 template<int dim, class GridImp>
00718 class YaspEntity<0,dim,GridImp> 
00719 : public EntityDefaultImplementation <0,dim,GridImp,YaspEntity>
00720 {
00721   enum { dimworld = GridImp::dimensionworld };
00722 public:
00723   typedef typename GridImp::ctype ctype;
00724   
00725   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
00726   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
00727 
00728   typedef YaspSpecialGeometry<dim-0,dim,GridImp> SpecialGeometry;
00729   
00730   typedef typename GridImp::template Codim<0>::Geometry Geometry;
00731   template <int cd>
00732   struct Codim
00733   {
00734     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
00735   };
00736   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
00737   typedef typename GridImp::LevelIntersectionIterator IntersectionIterator;
00738   typedef typename GridImp::LevelIntersectionIterator LevelIntersectionIterator;
00739   typedef typename GridImp::LeafIntersectionIterator  LeafIntersectionIterator;
00740   typedef typename GridImp::HierarchicIterator HierarchicIterator;
00741 
00743   typedef typename GridImp::PersistentIndexType PersistentIndexType;
00744 
00746   typedef typename YGrid<dim,ctype>::iTupel iTupel;
00747 
00748   // constructor
00749   YaspEntity (const GridImp * yg, const YGLI& g, const TSI& it)
00750     : _yg(yg), _it(it), _g(g), _geometry(it.position(),it.meshsize())
00751   {
00752   }
00753 
00755   int level () const { return _g.level(); }
00756 
00758   int index () const { return _it.superindex(); } // superindex works also for iteration over subgrids
00759 
00761   int globalIndex () const {
00762     return _g.cell_global().index(_it.coord());
00763   }
00764 
00766   PartitionType partitionType () const
00767   {
00768         if (_g.cell_interior().inside(_it.coord())) return InteriorEntity;
00769         if (_g.cell_overlap().inside(_it.coord())) return OverlapEntity;
00770         DUNE_THROW(GridError, "Impossible GhostEntity " << _it.coord() << "\t"
00771                   << _g.cell_interior().origin() << "/" << _g.cell_interior().size());
00772         return GhostEntity;
00773   }
00774 
00776   const Geometry& geometry () const { return _geometry; }
00777 
00780   template<int cc> int count () const
00781   {
00782     if (cc==dim) return 1<<dim;
00783     if (cc==1) return 2*dim;
00784     if (cc==dim-1) return dim*(1<<(dim-1));
00785     if (cc==0) return 1;
00786     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00787   }
00788 
00791   template<int cc>
00792   typename Codim<cc>::EntityPointer subEntity (int i) const
00793   {
00794       dune_static_assert( cc == dim || cc == 0 ,
00795           "YaspGrid only supports Entities with codim=dim and codim=0");
00796       // coordinates of the cell == coordinates of lower left corner
00797       if (cc==dim)
00798           {
00799                 iTupel coord = _it.coord();
00800 
00801                 // get corner from there
00802                 for (int k=0; k<dim; k++)
00803                   if (i&(1<<k)) (coord[k])++;
00804 
00805                 return YaspEntityPointer<cc,GridImp>(_yg,_g,_g.vertex_overlapfront().tsubbegin(coord));
00806           }
00807       if (cc==0)
00808           {
00809                 return YaspEntityPointer<cc,GridImp>(_yg,_g,_it);
00810           }
00811       DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
00812   }
00813 
00815   EntityPointer father () const
00816   {
00817         // check if coarse level exists
00818         if (_g.level()<=0)
00819           DUNE_THROW(GridError, "tried to call father on level 0");
00820 
00821         // yes, get iterator to it
00822         YGLI cg = _g.coarser();
00823 
00824         // coordinates of the cell
00825         iTupel coord = _it.coord();
00826 
00827         // get coordinates on next coarser level
00828         for (int k=0; k<dim; k++) coord[k] = coord[k]/2;
00829 
00830         return YaspEntityPointer<0,GridImp>(_yg,cg,cg.cell_overlap().tsubbegin(coord));
00831   }
00832   
00834   bool hasFather () const
00835   {
00836       return (_g.level()>0);
00837   }
00838 
00848   const Geometry& geometryInFather () const
00849   {
00850     // determine which son we are
00851     int son = 0;
00852     for (int k=0; k<dim; k++)
00853       if (_it.coord(k)%2)
00854         son += (1<<k);
00855     
00856     // configure one of the 2^dim transformations
00857     return YaspFatherRelativeLocalElement<dim,GridImp>::getson(son);
00858   }
00859 
00860   const TSI& transformingsubiterator () const
00861   {
00862         return _it;
00863   }
00864 
00865   const YGLI& gridlevel () const
00866   {
00867         return _g;
00868   }
00869 
00870   const GridImp* yaspgrid () const
00871   {
00872         return _yg;
00873   }
00874 
00875   bool isLeaf() const
00876   {
00877     return (_g.level() == _g.mg()->maxlevel());
00878   }
00879   
00882   bool isNew () const { return _yg->adaptRefCount > 0 && _g.mg()->maxlevel() < _g.level() + _yg->adaptRefCount; }
00883   
00886   bool mightVanish () const { return false; }
00887   // { return _yg->adaptRefCount < 0 && _g.mg()->maxlevel() < _g.level() - _yg->adaptRefCount; }
00888 
00890   IntersectionIterator ibegin () const
00891   {
00892     return YaspIntersectionIterator<GridImp>(*this,false);
00893   }
00894 
00896   LeafIntersectionIterator ileafbegin () const
00897   {
00898     // only if entity is leaf this iterator delivers intersections
00899     return YaspIntersectionIterator<GridImp>(*this, ! isLeaf() );
00900   }
00901 
00903   LevelIntersectionIterator ilevelbegin () const
00904   {
00905     return ibegin(); 
00906   }
00907 
00909   IntersectionIterator iend () const
00910   {
00911     return YaspIntersectionIterator<GridImp>(*this,true);
00912   }
00913 
00915   LeafIntersectionIterator ileafend () const
00916   {
00917     return iend();
00918   }
00919 
00921   LevelIntersectionIterator ilevelend () const
00922   {
00923     return iend();
00924   }
00925 
00930   HierarchicIterator hbegin (int maxlevel) const
00931   {
00932     return YaspHierarchicIterator<GridImp>(_yg,_g,_it,maxlevel);
00933   }
00934 
00936   HierarchicIterator hend (int maxlevel) const
00937   {
00938     return YaspHierarchicIterator<GridImp>(_yg,_g,_it,_g.level());
00939   }
00940 
00941 private:
00942   // IndexSets needs access to the private index methods
00943   friend class Dune::YaspLevelIndexSet<GridImp>;
00944   friend class Dune::YaspLeafIndexSet<GridImp>;
00945   friend class Dune::YaspGlobalIdSet<GridImp>;
00946 
00948   PersistentIndexType persistentIndex () const 
00949   { 
00950     // get size of global grid
00951     const iTupel& size =  _g.cell_global().size();
00952 
00953     // get coordinate correction for periodic boundaries
00954     int coord[dim];
00955     for (int i=0; i<dim; i++)
00956       {
00957         coord[i] = _it.coord(i);
00958         if (coord[i]<0) coord[i] += size[i];
00959         if (coord[i]>=size[i]) coord[i] -= size[i];
00960       }
00961 
00962     // encode codim
00963     PersistentIndexType id(0);
00964 
00965     // encode level
00966     id = id << yaspgrid_level_bits;
00967     id = id+PersistentIndexType(_g.level());
00968     
00969 
00970     // encode coordinates
00971     for (int i=dim-1; i>=0; i--)
00972       {
00973         id = id << yaspgrid_dim_bits;
00974         id = id+PersistentIndexType(coord[i]);
00975       }
00976     
00977     return id;
00978   }
00979 
00981   int compressedIndex () const 
00982   {
00983     return _it.superindex();
00984   }
00985 
00987   int compressedLeafIndex () const 
00988   {
00989     return _it.superindex();
00990   }
00991 
00993   PersistentIndexType subPersistentIndex (int i, int cc) const
00994   {
00995     if (cc==0)
00996       return persistentIndex();
00997 
00998     // get position of cell, note that global origin is zero
00999     // adjust for periodic boundaries
01000     int coord[dim];
01001     for (int k=0; k<dim; k++)
01002       {
01003         coord[k] = _it.coord(k);
01004         if (coord[k]<0) coord[k] += _g.cell_global().size(k);
01005         if (coord[k]>=_g.cell_global().size(k)) coord[k] -= _g.cell_global().size(k);
01006       }
01007     
01008     if (cc==dim)
01009       {
01010         // transform to vertex coordinates
01011         for (int k=0; k<dim; k++)
01012           if (i&(1<<k)) (coord[k])++;
01013 
01014         // determine min number of trailing zeroes
01015         int trailing = 1000;
01016         for (int i=0; i<dim; i++)
01017           {
01018             // count trailing zeros
01019             int zeros = 0;
01020             for (int j=0; j<_g.level(); j++)
01021               if (coord[i]&(1<<j))
01022                 break;
01023               else
01024                 zeros++;
01025             trailing = std::min(trailing,zeros);
01026           } 
01027 
01028         // determine the level of this vertex
01029         int level = _g.level()-trailing;
01030 
01031         // encode codim
01032         PersistentIndexType id(dim);
01033         
01034         // encode level
01035         id = id << yaspgrid_level_bits;
01036         id = id+PersistentIndexType(level);
01037         
01038         // encode coordinates
01039         for (int i=dim-1; i>=0; i--)
01040           {
01041             id = id << yaspgrid_dim_bits;
01042             id = id+PersistentIndexType(coord[i]>>trailing);
01043           }
01044         
01045         return id;
01046       }
01047 
01048     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01049       {
01050         // Idea: Use the doubled grid to assign coordinates to faces
01051 
01052         // ivar is the direction that varies
01053         int ivar=i/2; 
01054 
01055         // compute position from cell position
01056         for (int k=0; k<dim; k++)
01057           coord[k] = coord[k]*2 + 1; // the doubled grid
01058         if (i%2) 
01059           coord[ivar] += 1;
01060         else
01061           coord[ivar] -= 1;
01062 
01063         // encode codim
01064         PersistentIndexType id(1);
01065         
01066         // encode level
01067         id = id << yaspgrid_level_bits;
01068         id = id+PersistentIndexType(_g.level());
01069         
01070         // encode coordinates
01071         for (int i=dim-1; i>=0; i--)
01072           {
01073             id = id << yaspgrid_dim_bits;
01074             id = id+PersistentIndexType(coord[i]);
01075           }
01076         
01077         return id;
01078       }
01079 
01080     // map to old numbering
01081     static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
01082     i = edge[i];
01083 
01084     if (cc==dim-1) // edges, exist only for dim>2
01085       {
01086         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01087 
01088         // number of entities per direction
01089         int m=1<<(dim-1);
01090 
01091         // ifix is the direction that is fixed
01092         int ifix=(dim-1)-(i/m); 
01093 
01094         // compute position from cell position
01095         int bit=1;
01096         for (int k=0; k<dim; k++)
01097           {
01098             coord[k] = coord[k]*2+1; // cell position in doubled grid
01099             if (k==ifix) continue;
01100             if ((i%m)&bit) coord[k] += 1; else coord[k] -= 1;
01101             bit *= 2;
01102           } 
01103 
01104         // encode codim
01105         PersistentIndexType id(dim-1);
01106         
01107         // encode level
01108         id = id << yaspgrid_level_bits;
01109         id = id+PersistentIndexType(_g.level());
01110         
01111         // encode coordinates
01112         for (int i=dim-1; i>=0; i--)
01113           {
01114             id = id << yaspgrid_dim_bits;
01115             id = id+PersistentIndexType(coord[i]);
01116           }
01117         
01118         return id;
01119       }
01120 
01121     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01122   }
01123 
01125   int subCompressedIndex (int i, int cc) const
01126   {
01127     if (cc==0)
01128       return compressedIndex();
01129     
01130     // get cell position relative to origin of local cell grid
01131     iTupel coord;
01132     for (int k=0; k<dim; ++k) 
01133       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
01134 
01135     if (cc==dim) // vertices
01136       {
01137         // transform cell coordinate to corner coordinate
01138         for (int k=0; k<dim; k++)
01139           if (i&(1<<k)) (coord[k])++;
01140 
01141         // do lexicographic numbering
01142         int index = coord[dim-1];
01143         for (int k=dim-2; k>=0; --k)
01144           index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
01145         return index;
01146       }
01147 
01148     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01149       {
01150         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
01151 
01152         // ivar is the direction that varies
01153         int ivar=i/2; 
01154 
01155         // compute position from cell position
01156         if (i%2) coord[ivar] += 1; 
01157 
01158         // do lexicographic numbering
01159         int index = coord[dim-1];
01160         for (int k=dim-2; k>=0; --k)
01161           if (k==ivar)
01162             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01163           else
01164             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01165 
01166         // add size of all subsets for smaller directions
01167         for (int j=0; j<ivar; j++)
01168           {
01169             int n=_g.cell_overlap().size(j)+1;
01170             for (int l=0; l<dim; l++)
01171               if (l!=j) n *= _g.cell_overlap().size(l);
01172             index += n;
01173           }
01174 
01175         return index;
01176       }
01177 
01178     // map to old numbering
01179     static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
01180     i = edge[i];
01181 
01182     if (cc==dim-1) // edges, exist only for dim>2
01183       {
01184         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01185 
01186         // number of entities per direction
01187         int m=1<<(dim-1);
01188 
01189         // ifix is the direction that is fixed
01190         int ifix=(dim-1)-(i/m); 
01191 
01192         // compute position from cell position
01193         int bit=1;
01194         for (int k=0; k<dim; k++)
01195           {
01196             if (k==ifix) continue;
01197             if ((i%m)&bit) coord[k] += 1;
01198             bit *= 2;
01199           } 
01200 
01201         // do lexicographic numbering
01202         int index = coord[dim-1];
01203         for (int k=dim-2; k>=0; --k)
01204           if (k!=ifix)
01205             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01206           else
01207             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01208 
01209         // add size of all subsets for smaller directions
01210         for (int j=dim-1; j>ifix; j--)
01211           {
01212             int n=_g.cell_overlap().size(j);
01213             for (int l=0; l<dim; l++)
01214               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01215             index += n;
01216           }
01217 
01218         return index;
01219       }
01220 
01221     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01222   }
01223 
01225   int subCompressedLeafIndex (int i, int cc) const
01226   {
01227     if (cc==0)
01228       return compressedIndex();
01229     
01230     // get cell position relative to origin of local cell grid
01231     iTupel coord;
01232     for (int k=0; k<dim; ++k) 
01233       coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
01234 
01235     if (cc==dim) // vertices
01236       {
01237         // transform cell coordinate to corner coordinate
01238         for (int k=0; k<dim; k++)
01239           if (i&(1<<k)) (coord[k])++;
01240 
01241         // move coordinates up to maxlevel
01242         for (int k=0; k<dim; k++)
01243           coord[k] = coord[k]<<(_g.mg()->maxlevel()-_g.level());
01244 
01245         // do lexicographic numbering
01246         int index = coord[dim-1];
01247         for (int k=dim-2; k>=0; --k)
01248           index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01249         return index;
01250       }
01251 
01252     if (cc==1) // faces, i.e. for dim=2 codim=1 is treated as a face
01253       {
01254         // Idea: direction ivar varies, all others are fixed, i.e. 2 possibilities per direction
01255 
01256         // ivar is the direction that varies
01257         int ivar=i/2; 
01258 
01259         // compute position from cell position
01260         if (i%2) coord[ivar] += 1; 
01261 
01262         // do lexicographic numbering
01263         int index = coord[dim-1];
01264         for (int k=dim-2; k>=0; --k)
01265           if (k==ivar)
01266             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01267           else
01268             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01269 
01270         // add size of all subsets for smaller directions
01271         for (int j=0; j<ivar; j++)
01272           {
01273             int n=_g.cell_overlap().size(j)+1;
01274             for (int l=0; l<dim; l++)
01275               if (l!=j) n *= _g.cell_overlap().size(l);
01276             index += n;
01277           }
01278 
01279         return index;
01280       }
01281 
01282     // map to old numbering
01283     static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
01284     i = edge[i];
01285 
01286     if (cc==dim-1) // edges, exist only for dim>2
01287       {
01288         // Idea: direction i is fixed, all others are vary, i.e. 2^(dim-1) possibilities per direction
01289 
01290         // number of entities per direction
01291         int m=1<<(dim-1);
01292 
01293         // ifix is the direction that is fixed
01294         int ifix=(dim-1)-(i/m); 
01295 
01296         // compute position from cell position
01297         int bit=1;
01298         for (int k=0; k<dim; k++)
01299           {
01300             if (k==ifix) continue;
01301             if ((i%m)&bit) coord[k] += 1;
01302             bit *= 2;
01303           } 
01304 
01305         // do lexicographic numbering
01306         int index = coord[dim-1];
01307         for (int k=dim-2; k>=0; --k)
01308           if (k!=ifix)
01309             index = (index*(_g.cell_overlap().size(k)+1))+coord[k]; // one more
01310           else
01311             index = (index*(_g.cell_overlap().size(k)))+coord[k];
01312 
01313         // add size of all subsets for smaller directions
01314         for (int j=dim-1; j>ifix; j--)
01315           {
01316             int n=_g.cell_overlap().size(j);
01317             for (int l=0; l<dim; l++)
01318               if (l!=j) n *= _g.cell_overlap().size(l)+1;
01319             index += n;
01320           }
01321 
01322         return index;
01323       }
01324 
01325     DUNE_THROW(GridError, "codim " << cc << " (dim=" << dim << ") not (yet) implemented");
01326   }
01327 
01328   const GridImp * _yg;      // access to YaspGrid
01329   const TSI& _it;           // position in the grid level
01330   const YGLI& _g;           // access to grid level
01331   SpecialGeometry _geometry; // the element geometry
01332 };
01333 
01334 
01335 // specialization for codim=dim (vertex)
01336 template<int dim, class GridImp>
01337 class YaspEntity<dim,dim,GridImp> 
01338 : public EntityDefaultImplementation <dim,dim,GridImp,YaspEntity>
01339 {
01340   enum { dimworld = GridImp::dimensionworld };
01341 public:
01342   typedef typename GridImp::ctype ctype;
01343   
01344   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01345   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01346 
01347   typedef YaspSpecialGeometry<dim-dim,dim,GridImp> SpecialGeometry;
01348   
01349   typedef typename GridImp::template Codim<dim>::Geometry Geometry;
01350   template <int cd>
01351   struct Codim
01352   {
01353     typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
01354   };
01355   typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01356 
01358   typedef typename GridImp::PersistentIndexType PersistentIndexType;
01359 
01361   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01362 
01363   // constructor
01364   YaspEntity (const GridImp* yg, const YGLI& g, const TSI& it)
01365     : _yg(yg), _it(it), _g(g), _geometry(it.position())
01366   {}
01367 
01369   int level () const {return _g.level();}
01370 
01372   int index () const {return _it.superindex();}
01373 
01375   int globalIndex () const { return _g.cell_global().index(_it.coord()); }
01376 
01377 
01379   const Geometry& geometry () const { return _geometry; }
01380 
01382   PartitionType partitionType () const
01383   {
01384         if (_g.vertex_interior().inside(_it.coord())) return InteriorEntity;
01385         if (_g.vertex_interiorborder().inside(_it.coord())) return BorderEntity;
01386         if (_g.vertex_overlap().inside(_it.coord())) return OverlapEntity;
01387         if (_g.vertex_overlapfront().inside(_it.coord())) return FrontEntity;
01388         return GhostEntity;
01389   }
01390 
01391 private:
01392   // IndexSets needs access to the private index methods
01393   friend class Dune::YaspLevelIndexSet<GridImp>;
01394   friend class Dune::YaspLeafIndexSet<GridImp>;
01395   friend class Dune::YaspGlobalIdSet<GridImp>;
01396 
01398   PersistentIndexType persistentIndex () const 
01399   { 
01400     // get coordinate and size of global grid
01401     const iTupel& size =  _g.vertex_global().size();
01402     int coord[dim];
01403 
01404     // correction for periodic boundaries
01405     for (int i=0; i<dim; i++)
01406       {
01407         coord[i] = _it.coord(i);
01408         if (coord[i]<0) coord[i] += size[i];
01409         if (coord[i]>=size[i]) coord[i] -= size[i];
01410       }
01411 
01412     // determine min number of trailing zeroes
01413     int trailing = 1000;
01414     for (int i=0; i<dim; i++)
01415       {
01416         // count trailing zeros
01417         int zeros = 0;
01418         for (int j=0; j<_g.level(); j++)
01419           if (coord[i]&(1<<j))
01420             break;
01421           else
01422             zeros++;
01423         trailing = std::min(trailing,zeros);
01424       } 
01425 
01426     // determine the level of this vertex
01427     int level = _g.level()-trailing;
01428 
01429     // encode codim
01430     PersistentIndexType id(dim);
01431 
01432     // encode level
01433     id = id << yaspgrid_level_bits;
01434     id = id+PersistentIndexType(level);
01435     
01436     // encode coordinates
01437     for (int i=dim-1; i>=0; i--)
01438       {
01439         id = id << yaspgrid_dim_bits;
01440         id = id+PersistentIndexType(coord[i]>>trailing);
01441       }
01442     
01443     return id;
01444   }
01445 
01447   int compressedIndex () const { return _it.superindex();}
01448 
01450   int compressedLeafIndex () const 
01451   { 
01452     if (_g.level()==_g.mg()->maxlevel())
01453       return _it.superindex();
01454 
01455     // not on leaf level, interpolate to finest grid
01456     int coord[dim];
01457     for (int i=0; i<dim; i++) coord[i] = _it.coord(i)-(_g).vertex_overlap().origin(i);
01458 
01459     // move coordinates up to maxlevel (multiply by 2 for each level
01460     for (int k=0; k<dim; k++)
01461       coord[k] = coord[k]*(1<<(_g.mg()->maxlevel()-_g.level()));
01462 
01463     // do lexicographic numbering
01464     int index = coord[dim-1];
01465     for (int k=dim-2; k>=0; --k)
01466       index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
01467     return index;
01468   }
01469 
01470 public:  
01471   const TSI& transformingsubiterator() const { return _it; }
01472   const YGLI& gridlevel() const { return _g; }
01473   const GridImp * yaspgrid() const { return _yg; }
01474 protected:
01475   const GridImp * _yg;            // access to YaspGrid
01476   const TSI& _it;                 // position in the grid level
01477   const YGLI& _g;                 // access to grid level
01478   SpecialGeometry  _geometry;     // the element geometry
01479   // temporary object
01480   mutable FieldVector<ctype, dim> loc;   // always computed before being returned
01481 };
01482 
01483 //========================================================================
01488 //========================================================================
01489 
01490 template<class GridImp>
01491 class YaspIntersection
01492 {
01493     enum { dim=GridImp::dimension };
01494     enum { dimworld=GridImp::dimensionworld };  
01495     typedef typename GridImp::ctype ctype;
01496     YaspIntersection();
01497     YaspIntersection& operator = (const YaspIntersection&);
01498 public:
01499     // types used from grids
01500     typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01501     typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01502     typedef typename GridImp::template Codim<0>::Entity Entity;
01503     typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
01504     typedef typename GridImp::template Codim<1>::Geometry Geometry;
01505     typedef typename GridImp::template Codim<1>::LocalGeometry LocalGeometry;
01506     typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01507     typedef YaspSpecialGeometry<dim-1,dimworld,GridImp> SpecialGeometry;
01508     typedef YaspSpecialGeometry<dim-1,dim,GridImp> SpecialLocalGeometry;
01509     typedef Dune::Intersection<const GridImp, Dune::YaspIntersectionIterator> Intersection;
01510   
01511     void update() const {
01512         const_cast<YaspIntersection*>(this)->update();
01513     }
01514     void update() {
01515         if (_count == 2*_dir + _face || _count >= 2*dim)
01516             return;
01517 
01518 #if 0
01519         // update intersection iterator from current position
01520         if (_count == 2*_dir + 1 && _face==0) // direction remains valid
01521         {
01522             _face = 1; // 0->1, _dir remains
01523             
01524             // move transforming iterator
01525             _outside.transformingsubiterator().move(_dir,2); // move two cells in positive direction
01526                     
01527             // make up faces
01528             _pos_world[_dir] += _inside.transformingsubiterator().meshsize(_dir);
01529 
01530             return;
01531         }
01532         if (_count == 2*_dir + 1 && _face==1) // change direction
01533         {
01534             // move transforming iterator
01535             _outside.transformingsubiterator().move(_dir,-1); // move one cell back
01536                     
01537             // make up faces
01538             _pos_world[_dir] = _inside.transformingsubiterator().position(_dir);
01539                     
01540             _face = 0;
01541             _dir += 1;
01542                     
01543             // move transforming iterator
01544             _outside.transformingsubiterator().move(_dir,-1); // move one cell in negative direction
01545                     
01546             // make up faces
01547             _pos_world[_dir] -= 0.5*_inside.transformingsubiterator().meshsize(_dir);
01548 
01549             return;
01550         }
01551 #endif
01552                 
01553         // cleanup old stuff
01554         _outside.transformingsubiterator().move(_dir,1-2*_face); // move home
01555         _pos_world[_dir] = _inside.transformingsubiterator().position(_dir);
01556                 
01557         // update face info
01558         _dir = _count / 2;
01559         _face = _count % 2;
01560                 
01561         // move transforming iterator
01562         _outside.transformingsubiterator().move(_dir,-1+2*_face);
01563                 
01564         // make up faces
01565         _pos_world[_dir] += (-0.5+_face)*_inside.transformingsubiterator().meshsize(_dir);
01566     }
01567 
01569     void increment() {
01570         _count += (_count < 2*dim);
01571     }
01572 
01574     bool equals (const YaspIntersection& other) const
01575     {
01576         return _inside.equals(other._inside) && _count == other._count;
01577     }
01578 
01583     bool boundary () const
01584     {
01585 #if 1
01586         return (_inside.transformingsubiterator().coord(_count/2) + 2*(_count%2) - 1 < _inside.gridlevel().cell_global().min(_count/2)
01587             ||
01588             _inside.transformingsubiterator().coord(_count/2) + 2*(_count%2) - 1 > _inside.gridlevel().cell_global().max(_count/2));
01589 #else
01590         update();
01591         // The transforming iterator can be safely moved beyond the boundary.
01592         // So we only have to compare against the cell_global grid
01593         return (_outside.transformingsubiterator().coord(_dir) < _inside.gridlevel().cell_global().min(_dir)
01594             ||
01595             _outside.transformingsubiterator().coord(_dir) > _inside.gridlevel().cell_global().max(_dir));
01596 #endif
01597     }
01598 
01600     bool neighbor () const
01601     {
01602 #if 1
01603         return (_inside.transformingsubiterator().coord(_count/2) + 2*(_count%2) - 1 >= _inside.gridlevel().cell_overlap().min(_count/2)
01604             &&
01605             _inside.transformingsubiterator().coord(_count/2) + 2*(_count%2) - 1 <= _inside.gridlevel().cell_overlap().max(_count/2));
01606 #else
01607         update();
01608         return (_outside.transformingsubiterator().coord(_dir) >= _inside.gridlevel().cell_overlap().min(_dir)
01609             &&
01610             _outside.transformingsubiterator().coord(_dir) <= _inside.gridlevel().cell_overlap().max(_dir));
01611 #endif
01612     }
01613 
01615     bool conforming () const
01616     {
01617         return true;
01618     }
01619 
01622     EntityPointer inside() const
01623     {
01624         return _inside;
01625     }
01626   
01629     EntityPointer outside() const
01630     {
01631         update();
01632         return _outside;
01633     }
01634 
01637     int boundaryId() const
01638     {
01639         if(boundary()) return indexInInside()+1;
01640         return 0;
01641     }
01642   
01645     int boundarySegmentIndex() const
01646     {
01647         if(! boundary())
01648             DUNE_THROW(GridError, "called boundarySegmentIndex while boundary() == false");
01649         update();
01650         // size of local macro grid
01651         const FieldVector<int, dim> & size = _inside.gridlevel().mg()->begin().cell_overlap().size();
01652         const FieldVector<int, dim> & origin = _inside.gridlevel().mg()->begin().cell_overlap().origin();
01653         FieldVector<int, dim> sides;
01654         {
01655             for (int i=0; i<dim; i++)
01656             {
01657                 sides[i] =
01658                     ((_inside.gridlevel().mg()->begin().cell_overlap().origin(i)
01659                         == _inside.gridlevel().mg()->begin().cell_global().origin(i))+
01660                         (_inside.gridlevel().mg()->begin().cell_overlap().origin(i) +
01661                             _inside.gridlevel().mg()->begin().cell_overlap().size(i)
01662                             == _inside.gridlevel().mg()->begin().cell_global().origin(i) +
01663                             _inside.gridlevel().mg()->begin().cell_global().size(i)));
01664             }
01665         }
01666         // gobal position of the cell on macro grid
01667         FieldVector<int, dim> pos = _inside.transformingsubiterator().coord();
01668         pos /= (1<<_inside.level());
01669         pos -= origin;
01670         // compute unit-cube-face-sizes
01671         FieldVector<int, dim> fsize;
01672         {
01673             int vol = 1;
01674             for (int k=0; k<dim; k++)
01675                 vol *= size[k];
01676             for (int k=0; k<dim; k++)
01677                 fsize[k] = vol/size[k];
01678         }
01679         // compute index in the unit-cube-face
01680         int index = 0;
01681         {
01682             int localoffset = 1;
01683             for (int k=dim-1; k>=0; k--)
01684             {
01685                 if (k == _dir) continue;
01686                 index += (pos[k]) * localoffset;
01687                 localoffset *= size[k];
01688             }
01689         }
01690         // add unit-cube-face-offsets
01691         {
01692             for (int k=0; k<_dir; k++)
01693                 index += sides[k] * fsize[k];
01694             // add fsize if we are on the right face and there is a left-face-boundary on this processor
01695             index += _face * (sides[_dir]>1) * fsize[_dir];
01696         }
01697 
01698         // int rank = 0;
01699         // MPI_Comm_rank(MPI_COMM_WORLD, &rank);
01700         // std::cout << rank << "... size: " << size << " sides: " << sides
01701         //           << " fsize: " << fsize
01702         //           << " pos: " << pos << " face: " << int(_dir) << "/" << int(_face)
01703         //           << " index: " << index << std::endl;
01704         
01705         return index;
01706     }
01707   
01709     FieldVector<ctype, dimworld> outerNormal (const FieldVector<ctype, dim-1>& local) const
01710     {
01711         return _faceInfo[_count].normal;
01712     }
01713 
01715     FieldVector<ctype, dimworld> unitOuterNormal (const FieldVector<ctype, dim-1>& local) const
01716     {
01717         return _faceInfo[_count].normal;
01718     }
01719 
01721     FieldVector<ctype, dimworld> centerUnitOuterNormal () const
01722     {
01723         return _faceInfo[_count].normal;
01724     }
01725 
01729     FieldVector<ctype, dimworld> integrationOuterNormal (const FieldVector<ctype, dim-1>& local) const
01730     {
01731         FieldVector<ctype, dimworld> n = _faceInfo[_count].normal;
01732         n *= geometry().volume();
01733         return n;
01734     }
01735  
01739     const LocalGeometry &geometryInInside () const
01740     {
01741         return _faceInfo[_count].geom_inside;
01742     }
01743 
01747     const LocalGeometry &geometryInOutside () const
01748     {
01749         return _faceInfo[_count].geom_outside;
01750     }
01751 
01755     const Geometry &geometry () const
01756     {
01757         update();
01758         return _is_global;
01759     }
01760 
01762     GeometryType type () const
01763     {
01764         static const GeometryType cube(GeometryType::cube, dim-1);
01765         return cube;
01766     }
01767 
01769     int indexInInside () const
01770     {
01771         return _count;
01772     }
01773 
01775     int indexInOutside () const
01776     {
01777         // flip the last bit
01778         return _count^1;
01779     }
01780 
01782     YaspIntersection (const YaspEntity<0,dim,GridImp>& myself, bool toend) :
01783         _inside(myself.yaspgrid(), myself.gridlevel(),
01784             myself.transformingsubiterator()),
01785         _outside(myself.yaspgrid(), myself.gridlevel(),
01786             myself.transformingsubiterator()),
01787         // initialize to first neighbor
01788         _count(0),
01789         _dir(0),
01790         _face(0),
01791         _pos_world(myself.transformingsubiterator().position()),
01792         _is_global(_pos_world,_inside.transformingsubiterator().meshsize(),_dir)
01793     {
01794         if (toend)
01795         {
01796             // initialize end iterator
01797             _count = 2*dim;
01798             return;
01799         }
01800         _count = 0;
01801         
01802         // move transforming iterator
01803         _outside.transformingsubiterator().move(_dir,-1);
01804 
01805         // make up faces
01806         _pos_world[0] -= 0.5*_inside.transformingsubiterator().meshsize(0);
01807     }
01808 
01810     YaspIntersection (const YaspIntersection& it) :
01811         _inside(it._inside),
01812         _outside(it._outside),
01813         _count(it._count),
01814         _dir(it._dir),
01815         _face(it._face),
01816         _pos_world(it._pos_world),
01817         // Important: geometry must be recreated -- not copied!
01818         _is_global(_pos_world,_inside.transformingsubiterator().meshsize(),_dir)
01819     {}
01820 
01822     void assign (const YaspIntersection& it)
01823     {
01824         _inside = it._inside;
01825         _outside = it._outside;
01826         _count = it._count;
01827         _dir = it._dir;
01828         _face = it._face;
01829         _pos_world = it._pos_world;
01830         // Important: geometry is automatically updated
01831     }
01832 
01833 private:
01834     /* EntityPointers (get automatically updated) */
01835     mutable YaspEntityPointer<0,GridImp> _inside;  
01836     mutable YaspEntityPointer<0,GridImp> _outside; 
01837     /* current position */
01838     uint8_t _count;                                
01839     mutable uint8_t _dir;                          
01840     mutable uint8_t _face;                         
01841     /* current position */
01842     FieldVector<ctype, dimworld> _pos_world;       
01843     /* geometry object (get automatically updated) */
01844     SpecialGeometry         _is_global;            
01845 
01846     /* static data */
01847     static const FieldVector<typename GridImp::ctype, GridImp::dimension> _ext_local;
01848     struct faceInfo
01849     {
01850         FieldVector<ctype, dim> pos_inside;    
01851         FieldVector<ctype, dim> pos_outside;   
01852         uint8_t                 dir;
01853         FieldVector<ctype, dimworld> normal;
01854         SpecialLocalGeometry    geom_inside;   
01855         SpecialLocalGeometry    geom_outside;  
01856         faceInfo() :
01857             geom_inside (pos_inside, _ext_local,dir),
01858             geom_outside(pos_outside,_ext_local,dir) {}
01859     };
01860 
01861     /* static face info */
01862     static const array<faceInfo, 2*GridImp::dimension> _faceInfo;
01863     
01864     static array<faceInfo, 2*dim> initFaceInfo()
01865     {
01866         array<faceInfo, 2*dim> I;
01867         for (int i=0; i<dim; i++)
01868         {
01869             // center of face
01870             FieldVector<ctype, dim> a(0.5); a[i] = 0.0;
01871             FieldVector<ctype, dim> b(0.5); b[i] = 1.0;
01872             I[2*i].pos_inside = a;
01873             I[2*i].pos_outside = b;
01874             I[2*i+1].pos_inside = b;
01875             I[2*i+1].pos_outside = a;
01876             // direction
01877             I[2*i].dir = i;
01878             I[2*i+1].dir = i;
01879             // normal vectors
01880             I[2*i].normal = 0.0;
01881             I[2*i+1].normal = 0.0;
01882             I[2*i].normal[i] = -1.0;
01883             I[2*i+1].normal[i] = +1.0;
01884         }
01885         return I;
01886     }
01887 };
01888 
01889 template<class GridImp>
01890 const FieldVector<typename GridImp::ctype, GridImp::dimension>
01891 YaspIntersection<GridImp>::_ext_local = FieldVector<typename GridImp::ctype, GridImp::dimension>(1.0);
01892 
01893 template<class GridImp>
01894 const array<typename YaspIntersection<GridImp>::faceInfo, 2*GridImp::dimension>
01895 YaspIntersection<GridImp>::_faceInfo =
01896     YaspIntersection<GridImp>::initFaceInfo();
01897 
01898 //========================================================================
01903 //========================================================================
01904 
01905 template<class GridImp>
01906 class YaspIntersectionIterator : public MakeableInterfaceObject< Dune::Intersection<const GridImp, Dune::YaspIntersection > >
01907 {
01908     enum { dim=GridImp::dimension };
01909     enum { dimworld=GridImp::dimensionworld };
01910     typedef typename GridImp::ctype ctype;
01911     YaspIntersectionIterator();
01912 public:
01913     // types used from grids
01914     typedef Dune::Intersection<const GridImp, Dune::YaspIntersection> Intersection;
01915     typedef MakeableInterfaceObject<Intersection> MakeableIntersection;
01916   
01918     void increment()
01919     {
01920         GridImp::getRealImplementation(*this).increment();
01921     }
01922 
01924     bool equals (const YaspIntersectionIterator& other) const
01925     {
01926         return GridImp::getRealImplementation(*this).equals(
01927             GridImp::getRealImplementation(other));
01928     }
01929 
01931     const Intersection & dereference() const
01932     {
01933         return *this;
01934     }
01935   
01937     YaspIntersectionIterator (const YaspEntity<0,dim,GridImp>& myself, bool toend) :
01938         MakeableIntersection(YaspIntersection<GridImp>(myself, toend))
01939     {}
01940 
01942     YaspIntersectionIterator (const YaspIntersectionIterator& it) :
01943         MakeableIntersection(it)
01944     {}
01945 
01947     YaspIntersectionIterator & operator = (const YaspIntersectionIterator& it)
01948     {
01949         GridImp::getRealImplementation(*this).assign(
01950             GridImp::getRealImplementation(it));
01951         return *this;
01952     }
01953 };
01954 
01955 //========================================================================
01959 //========================================================================
01960 
01961 template<class GridImp>
01962 class YaspHierarchicIterator :
01963   public YaspEntityPointer<0,GridImp>
01964 {
01965   enum { dim=GridImp::dimension };
01966   enum { dimworld=GridImp::dimensionworld };  
01967   typedef typename GridImp::ctype ctype;
01968 public:
01969   // types used from grids
01970   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
01971   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
01972   typedef typename GridImp::template Codim<0>::Entity Entity;
01973   typedef YaspSpecialEntity<0,dim,GridImp> SpecialEntity;
01974 
01976   typedef typename YGrid<dim,ctype>::iTupel iTupel;
01977 
01979   YaspHierarchicIterator (const GridImp* yg, const YGLI& g, const TSI& it, int maxlevel) :
01980     YaspEntityPointer<0,GridImp>(yg,g,it)
01981   {
01982       // now iterator points to current cell
01983       StackElem se(this->_g);
01984       se.coord = this->_it.coord();
01985       stack.push(se);
01986 
01987       // determine maximum level
01988       _maxlevel = std::min(maxlevel,this->_g.mg()->maxlevel());
01989 
01990       // if maxlevel not reached then push yourself and sons
01991       if (this->_g.level()<_maxlevel)
01992       {
01993         push_sons();
01994       }
01995 
01996       // and make iterator point to first son if stack is not empty
01997       if (!stack.empty())
01998         pop_tos();
01999   }
02000 
02002   YaspHierarchicIterator (const YaspHierarchicIterator& it) :
02003     YaspEntityPointer<0,GridImp>(it),
02004     _maxlevel(it._maxlevel), stack(it.stack)
02005   {}
02006 
02008   void increment ()
02009   {
02010         // sanity check: do nothing when stack is empty
02011         if (stack.empty()) return;
02012 
02013         // if maxlevel not reached then push sons
02014         if (this->_g.level()<_maxlevel)
02015           push_sons();
02016 
02017         // in any case pop one element
02018         pop_tos();
02019   }
02020 
02021   void print (std::ostream& s) const
02022   {
02023         s << "HIER: " << "level=" << this->_g.level()
02024           << " position=" << this->_it.coord()
02025           << " superindex=" << this->_it.superindex()
02026           << " maxlevel=" << this->_maxlevel
02027           << " stacksize=" << stack.size()
02028           << std::endl;
02029   }
02030 
02031 private:
02032   int _maxlevel;         
02033 
02034   struct StackElem {
02035         YGLI g;       // grid level of the element
02036         iTupel coord; // and the coordinates
02037         StackElem(YGLI gg) : g(gg) {}
02038   };
02039   std::stack<StackElem> stack;      
02040 
02041   // push sons of current element on the stack
02042   void push_sons ()
02043   {
02044         // yes, process all 1<<dim sons
02045         StackElem se(this->_g.finer());
02046         for (int i=0; i<(1<<dim); i++)
02047           {
02048                 for (int k=0; k<dim; k++)
02049                   if (i&(1<<k))
02050                         se.coord[k] = this->_it.coord(k)*2+1;
02051                   else
02052                         se.coord[k] = this->_it.coord(k)*2;
02053                 stack.push(se);
02054           }
02055   }
02056 
02057   // make TOS the current element
02058   void pop_tos ()
02059   {
02060         StackElem se = stack.top();
02061         stack.pop();
02062         this->_g = se.g;
02063         this->_it.reinit(this->_g.cell_overlap(),se.coord);
02064   }
02065 };
02066 
02067 //========================================================================
02077 //========================================================================
02078 template<int codim, class GridImp>
02079 class YaspEntityPointer
02080 {
02082   enum { dim=GridImp::dimension };
02084   enum { dimworld=GridImp::dimensionworld };
02085   typedef typename GridImp::ctype ctype;
02086 public:
02087   typedef typename GridImp::template Codim<codim>::Entity Entity;
02088   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
02089   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
02090   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
02091   typedef YaspEntityPointer<codim,GridImp> EntityPointerImp;
02092 protected:  
02093   typedef YaspEntity<codim, dim, GridImp> YaspEntityImp; 
02094 
02095 public:
02097   enum { codimension = codim };
02098   
02100   YaspEntityPointer (const GridImp * yg, const YGLI & g, const TSI & it)
02101     : _g(g), _it(it), _entity(yg, _g,_it)
02102   {
02103         if (codim>0 && codim<dim)
02104           {
02105                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
02106           }
02107   }
02108 
02110   YaspEntityPointer (const YaspEntityImp& entity) 
02111     : _g(entity.gridlevel()), 
02112       _it(entity.transformingsubiterator()), 
02113       _entity(entity.yaspgrid(),_g,_it)
02114   {
02115         if (codim>0 && codim<dim)
02116           {
02117                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
02118           }
02119   }
02120 
02122   YaspEntityPointer (const YaspEntityPointer& rhs)
02123     : _g(rhs._g), _it(rhs._it), _entity(rhs._entity.yaspgrid(),_g,_it)
02124   {
02125         if (codim>0 && codim<dim)
02126           {
02127                 DUNE_THROW(GridError, "YaspEntityPointer: codim not implemented");
02128           }
02129   }
02130 
02132   bool equals (const YaspEntityPointer& rhs) const
02133   {
02134         return (_it==rhs._it && _g == rhs._g);
02135   }
02136 
02138   Entity& dereference() const
02139   {
02140         return _entity;
02141   }
02142 
02144   void compactify() {}
02145 
02147   int level () const {return _g.level();}
02148 
02149   const YaspEntityPointer&
02150   operator = (const YaspEntityPointer& rhs)
02151     {
02152       _g = rhs._g;
02153       _it = rhs._it;
02154       /* _entity = i._entity
02155        * is done implicitely, as the entity is completely
02156        * defined via the interator it belongs to
02157        */
02158       return *this;
02159     }  
02160   
02161   const TSI& transformingsubiterator () const
02162   {
02163         return _it;
02164   }
02165 
02166   const YGLI& gridlevel () const
02167   {
02168         return _g;
02169   }
02170 
02171   TSI& transformingsubiterator ()
02172   {
02173         return _it;
02174   }
02175 
02176   YGLI& gridlevel ()
02177   {
02178         return _g;
02179   }
02180 
02181 protected:
02182   YGLI _g;               // access to grid level
02183   TSI _it;               // position in the grid level
02184   mutable SpecialEntity _entity; 
02185 };
02186 
02187 //========================================================================
02195 //========================================================================
02196 
02197 template<int codim, PartitionIteratorType pitype, class GridImp>
02198 class YaspLevelIterator :
02199   public YaspEntityPointer<codim,GridImp>
02200 {
02202   enum { dim=GridImp::dimension };
02204   enum { dimworld=GridImp::dimensionworld };
02205   typedef typename GridImp::ctype ctype;
02206 public:
02207   typedef typename GridImp::template Codim<codim>::Entity Entity;
02208   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
02209   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
02210   typedef YaspSpecialEntity<codim,dim,GridImp> SpecialEntity;
02211 
02213   YaspLevelIterator (const GridImp * yg, const YGLI & g, const TSI & it) :
02214     YaspEntityPointer<codim,GridImp>(yg,g,it) {}
02215 
02217   YaspLevelIterator (const YaspLevelIterator& i) :
02218     YaspEntityPointer<codim,GridImp>(i) {}
02219 
02221   void increment()
02222   {
02223         ++(this->_it);
02224   }
02225 };
02226 
02227 
02228 //========================================================================
02233 //========================================================================
02234 
02235 template<class GridImp>
02236 class YaspLevelIndexSet
02237   : public IndexSet< GridImp, YaspLevelIndexSet< GridImp >, unsigned int >
02238 {
02239   typedef YaspLevelIndexSet< GridImp > This;
02240   typedef IndexSet< GridImp, This, unsigned int > Base;
02241 
02242 public:
02243   typedef typename Base::IndexType IndexType;
02244 
02245   using Base::subIndex;
02246   
02248   YaspLevelIndexSet ( const GridImp &g, int l )
02249   : grid( g ),
02250     level( l )
02251   {
02252     // contains a single element type;
02253     for (int codim=0; codim<=GridImp::dimension; codim++)
02254       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
02255   }
02256 
02258   template<int cc>
02259   IndexType index (const typename GridImp::Traits::template Codim<cc>::Entity& e) const 
02260   {
02261     assert( cc == 0 || cc == GridImp::dimension );
02262     return grid.getRealImplementation(e).compressedIndex(); 
02263   }
02264 
02266   template< int cc >
02267   IndexType subIndex ( const typename remove_const< GridImp >::type::Traits::template Codim< cc >::Entity &e,
02268                        int i, unsigned int codim ) const
02269   {
02270     assert( cc == 0 || cc == GridImp::dimension );
02271     if( cc == GridImp::dimension )
02272       return grid.getRealImplementation(e).compressedIndex(); 
02273     else
02274       return grid.getRealImplementation(e).subCompressedIndex(i,codim);
02275   }
02276 
02278   int size (GeometryType type) const
02279   {
02280     return grid.size( level, type );
02281   }
02282 
02284   int size (int codim) const
02285   {
02286     return grid.size( level, codim );
02287   }
02288   
02290   template<class EntityType>
02291   bool contains (const EntityType& e) const
02292   {
02293     return e.level() == level;
02294   }
02295   
02297   const std::vector<GeometryType>& geomTypes (int codim) const
02298   {
02299     return mytypes[codim];
02300   }
02301 
02302 private:
02303   const GridImp& grid;
02304   int level;
02305   std::vector<GeometryType> mytypes[GridImp::dimension+1];
02306 };
02307 
02308 
02309 // Leaf Index Set
02310 
02311 template<class GridImp>
02312 class YaspLeafIndexSet
02313   : public IndexSet< GridImp, YaspLeafIndexSet< GridImp >, unsigned int >
02314 {
02315   typedef YaspLeafIndexSet< GridImp > This;
02316   typedef IndexSet< GridImp, This, unsigned int > Base;
02317 
02318 public:
02319   typedef typename Base::IndexType IndexType;
02320 
02321   using Base::subIndex;
02322   
02324   explicit YaspLeafIndexSet ( const GridImp &g )
02325   : grid( g )
02326   {
02327     // contains a single element type;
02328     for (int codim=0; codim<=GridImp::dimension; codim++)
02329       mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
02330   }
02331 
02333   template<int cc>
02334   IndexType index (const typename GridImp::Traits::template Codim<cc>::Entity& e) const 
02335   {
02336     assert( cc == 0 || cc == GridImp::dimension );
02337     return grid.getRealImplementation(e).compressedIndex(); 
02338   }
02339 
02341   template< int cc >
02342   IndexType subIndex ( const typename remove_const< GridImp >::type::Traits::template Codim< cc >::Entity &e,
02343                        int i, unsigned int codim ) const
02344   {
02345     assert( cc == 0 || cc == GridImp::dimension );
02346     if( cc == GridImp::dimension )
02347       return grid.getRealImplementation(e).compressedIndex(); 
02348     else
02349       return grid.getRealImplementation(e).subCompressedIndex(i,codim);
02350   }
02351 
02353   int size (GeometryType type) const
02354   {
02355     return grid.size( grid.maxLevel(), type );
02356   }
02357 
02359   int size (int codim) const
02360   {
02361     return grid.size( grid.maxLevel(), codim );
02362   }
02363 
02365   template<class EntityType>
02366   bool contains (const EntityType& e) const
02367   {
02368     //return e.isLeaf();
02369     return (e.level() == grid.maxLevel());
02370   }
02371   
02373   const std::vector<GeometryType>& geomTypes (int codim) const
02374   {
02375     return mytypes[codim];
02376   }
02377 
02378 private:
02379   const GridImp& grid;
02380   enum { ncodim = remove_const<GridImp>::type::dimension+1 };
02381   std::vector<GeometryType> mytypes[ncodim];
02382 };
02383 
02384 
02385 
02386 
02387 //========================================================================
02392 //========================================================================
02393 
02394 template<class GridImp>
02395 class YaspGlobalIdSet : public IdSet<GridImp,YaspGlobalIdSet<GridImp>,
02396                                      typename remove_const<GridImp>::type::PersistentIndexType >
02397 /*
02398   We used the remove_const to extract the Type from the mutable class,
02399   because the const class is not instantiated yet.
02400 */
02401 {
02402   typedef YaspGlobalIdSet< GridImp > This;
02403 
02404 public:
02406   typedef typename remove_const<GridImp>::type::PersistentIndexType IdType;
02407 
02408   using IdSet<GridImp, This, IdType>::subId;
02409 
02411   explicit YaspGlobalIdSet ( const GridImp &g )
02412   : grid( g )
02413   {}
02414 
02416   /*
02417     We use the remove_const to extract the Type from the mutable class,
02418     because the const class is not instantiated yet.
02419   */
02420   template<int cd>
02421   IdType id (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const 
02422   {
02423     return grid.getRealImplementation(e).persistentIndex();
02424   }
02425 
02427   /*
02428     We use the remove_const to extract the Type from the mutable class,
02429     because the const class is not instantiated yet.
02430   */
02431   IdType subId (const typename remove_const<GridImp>::type::Traits::template Codim< 0 >::Entity &e,
02432                  int i, unsigned int codim ) const
02433   {
02434     return grid.getRealImplementation(e).subPersistentIndex(i,codim);
02435   }
02436 
02437 private:
02438   const GridImp& grid;
02439 };
02440 
02441 
02442 template<int dim, int dimworld>
02443 struct YaspGridFamily
02444 {
02445 #if HAVE_MPI
02446   typedef CollectiveCommunication<MPI_Comm> CCType;
02447 #else
02448   typedef CollectiveCommunication<Dune::YaspGrid<dim> > CCType;
02449 #endif
02450 
02451   typedef GridTraits<dim,dimworld,Dune::YaspGrid<dim>,
02452                      YaspGeometry,YaspEntity,
02453                      YaspEntityPointer,YaspLevelIterator,
02454                      YaspIntersection, // leaf  intersection 
02455                      YaspIntersection, // level intersection 
02456                      YaspIntersectionIterator, // leaf  intersection iter 
02457                      YaspIntersectionIterator, // level intersection iter 
02458                      YaspHierarchicIterator,
02459                      YaspLevelIterator,
02460                      YaspLevelIndexSet< const YaspGrid< dim > >,
02461                      YaspLeafIndexSet< const YaspGrid< dim > >,
02462                      YaspGlobalIdSet<const YaspGrid<dim> >, 
02463                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02464                      YaspGlobalIdSet<const YaspGrid<dim> >, 
02465                      bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
02466                      CCType> 
02467   Traits;  
02468 };
02469 
02470 template<int dim, int codim>
02471 struct YaspCommunicateMeta {
02472   template<class G, class DataHandle>
02473   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02474   {
02475     if (data.contains(dim,codim))
02476       {
02477         DUNE_THROW(GridError, "interface communication not implemented");
02478       }
02479     YaspCommunicateMeta<dim,codim-1>::comm(g,data,iftype,dir,level);
02480   }
02481 };
02482 
02483 template<int dim>
02484 struct YaspCommunicateMeta<dim,dim> {
02485   template<class G, class DataHandle>
02486   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02487   {
02488     if (data.contains(dim,dim))
02489       g.template communicateCodim<DataHandle,dim>(data,iftype,dir,level);
02490     YaspCommunicateMeta<dim,dim-1>::comm(g,data,iftype,dir,level);
02491   }
02492 };
02493 
02494 template<int dim>
02495 struct YaspCommunicateMeta<dim,0> {
02496   template<class G, class DataHandle>
02497   static void comm (const G& g, DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level)
02498   {
02499     if (data.contains(dim,0))
02500       g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
02501   }
02502 };
02503 
02504 
02505 //************************************************************************
02521 template<int dim>
02522 class YaspGrid :
02523   public GridDefaultImplementation<dim,dim,yaspgrid_ctype,YaspGridFamily<dim,dim> >,
02524   public MultiYGrid<dim,yaspgrid_ctype>
02525 {
02526   typedef const YaspGrid<dim> GridImp;
02527 
02528   void init()
02529   {
02530     setsizes();
02531     indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim> >(*this,0) );
02532     theleafindexset.push_back( new YaspLeafIndexSet<const YaspGrid<dim> >(*this) );
02533     theglobalidset.push_back( new YaspGlobalIdSet<const YaspGrid<dim> >(*this) );
02534     boundarysegmentssize();
02535   }
02536 
02537   void boundarysegmentssize()
02538   {
02539       // sizes of local macro grid
02540       const FieldVector<int, dim> & size = MultiYGrid<dim,ctype>::begin().cell_overlap().size();
02541       FieldVector<int, dim> sides;
02542       {
02543           for (int i=0; i<dim; i++)
02544           {
02545               sides[i] =
02546                   ((MultiYGrid<dim,ctype>::begin().cell_overlap().origin(i)
02547                       == MultiYGrid<dim,ctype>::begin().cell_global().origin(i))+
02548                       (MultiYGrid<dim,ctype>::begin().cell_overlap().origin(i) +
02549                           MultiYGrid<dim,ctype>::begin().cell_overlap().size(i)
02550                           == MultiYGrid<dim,ctype>::begin().cell_global().origin(i) +
02551                           MultiYGrid<dim,ctype>::begin().cell_global().size(i)));
02552           }
02553       }
02554       nBSegments = 0;
02555       for (int k=0; k<dim; k++)
02556       {
02557           int offset = 1;
02558           for (int l=0; l<dim; l++)
02559           {
02560               if (l==k) continue;
02561               offset *= size[l];
02562           }
02563           nBSegments += sides[k]*offset;
02564       }
02565   }
02566 
02567 public:
02568 
02569   using MultiYGrid<dim,yaspgrid_ctype>::defaultLoadbalancer;
02570 
02572   typedef yaspgrid_ctype ctype;
02573 
02574   // define the persistent index type
02575   typedef bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits> PersistentIndexType;
02576 
02578   typedef YaspGridFamily<dim,dim> GridFamily;
02579   // the Traits 
02580   typedef typename YaspGridFamily<dim,dim>::Traits Traits;
02581 
02582   // need for friend declarations in entity
02583   typedef YaspLevelIndexSet<YaspGrid<dim> > LevelIndexSetType;
02584   typedef YaspLeafIndexSet<YaspGrid<dim> > LeafIndexSetType;
02585   typedef YaspGlobalIdSet<YaspGrid<dim> > GlobalIdSetType;
02586 
02588   enum { MAXL=64 };
02589 
02591   typedef MultiYGrid<dim,ctype> YMG;
02592   typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator YGLI;
02593   typedef typename SubYGrid<dim,ctype>::TransformingSubIterator TSI;
02594   typedef typename MultiYGrid<dim,ctype>::Intersection IS;
02595   typedef typename std::deque<IS>::const_iterator ISIT;
02596 
02605 #if HAVE_MPI
02606   YaspGrid (MPI_Comm comm, Dune::FieldVector<ctype, dim> L, 
02607             Dune::FieldVector<int, dim> s, 
02608             Dune::FieldVector<bool, dim> periodic, int overlap,
02609             const YLoadBalance<dim>* lb = defaultLoadbalancer())
02610       : YMG(comm,L,s,periodic,overlap,lb), ccobj(comm),
02611         keep_ovlp(true), adaptRefCount(0), adaptActive(false)
02612   {
02613     init();
02614   }
02615 
02616 #else
02617   YaspGrid (Dune::FakeMPIHelper::MPICommunicator, Dune::FieldVector<ctype, dim> L, 
02618             Dune::FieldVector<int, dim> s, 
02619             Dune::FieldVector<bool, dim> periodic, int overlap,
02620             const YLoadBalance<dim>* lb = defaultLoadbalancer())
02621       : YMG(L,s,periodic,overlap,lb),
02622         keep_ovlp(true), adaptRefCount(0), adaptActive(false)
02623   {
02624     init();
02625   }
02626 #endif
02627 
02628   
02641   YaspGrid (Dune::FieldVector<ctype, dim> L, 
02642             Dune::FieldVector<int, dim> s, 
02643             Dune::FieldVector<bool, dim> periodic, int overlap,
02644             const YLoadBalance<dim>* lb = YMG::defaultLoadbalancer())
02645 #if HAVE_MPI
02646       : MultiYGrid<dim,ctype>(MPI_COMM_SELF,L,s,periodic,overlap,lb), ccobj(MPI_COMM_SELF),
02647         keep_ovlp(true), adaptRefCount(0), adaptActive(false)
02648 #else
02649       : MultiYGrid<dim,ctype>(L,s,periodic,overlap,lb), 
02650         keep_ovlp(true), adaptRefCount(0), adaptActive(false)
02651 #endif
02652   {
02653     init();
02654   }
02655 
02656   ~YaspGrid()
02657   {
02658     deallocatePointers(indexsets);
02659     deallocatePointers(theleafindexset);
02660     deallocatePointers(theglobalidset);
02661   }
02662 
02663   private:
02664   // do not copy this class 
02665   YaspGrid(const YaspGrid&);
02666 
02667   public:
02668   
02672   int maxLevel() const {return MultiYGrid<dim,ctype>::maxlevel();} // delegate
02673 
02675   void globalRefine (int refCount)
02676   {
02677     if (refCount < -maxLevel())
02678       DUNE_THROW(GridError, "Only " << maxLevel() << " levels left. " <<
02679                  "Coarsening " << -refCount << " levels requested!");
02680     for (int k=refCount; k<0; k++)
02681       {
02682         MultiYGrid<dim,ctype>::coarsen();
02683         setsizes();
02684         indexsets.pop_back();
02685       }
02686     for (int k=0; k<refCount; k++)
02687       {
02688         MultiYGrid<dim,ctype>::refine(keep_ovlp);
02689         setsizes();
02690         indexsets.push_back( new YaspLevelIndexSet<const YaspGrid<dim> >(*this,maxLevel()) );
02691       }
02692   }
02693 
02698   void refineOptions (bool keepPhysicalOverlap)
02699   {
02700     keep_ovlp = keepPhysicalOverlap;
02701   }
02702 
02714   bool mark( int refCount, const typename Traits::template Codim<0>::Entity & e )
02715   {
02716     assert(adaptActive == false);
02717     if (e.level() != maxLevel()) return false;
02718     adaptRefCount = std::max(adaptRefCount, refCount);
02719     return true;
02720   }
02721 
02728   int getMark ( const typename Traits::template Codim<0>::Entity &e ) const
02729   {
02730     return ( e.level() == maxLevel() ) ? adaptRefCount : 0;
02731   }
02732   
02734   bool adapt ()   
02735   {
02736     globalRefine(adaptRefCount);
02737     return (adaptRefCount > 0);
02738   }
02739   
02741   bool preAdapt ()
02742   {
02743     adaptActive = true;
02744     adaptRefCount = comm().max(adaptRefCount);
02745     return (adaptRefCount < 0);
02746   }
02747 
02749   void postAdapt()
02750   {
02751     adaptActive = false;
02752     adaptRefCount = 0;
02753   }
02754 
02756   template<int cd, PartitionIteratorType pitype>
02757   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const
02758     {
02759       return levelbegin<cd,pitype>(level);
02760     }
02761 
02763   template<int cd, PartitionIteratorType pitype>
02764   typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const
02765     {
02766       return levelend<cd,pitype>(level);
02767     }
02768 
02770   template<int cd>
02771   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
02772     {
02773       return levelbegin<cd,All_Partition>(level);
02774     }
02775 
02777   template<int cd>
02778   typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
02779     {
02780       return levelend<cd,All_Partition>(level);
02781     }
02782 
02784   template<int cd, PartitionIteratorType pitype>
02785   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const
02786     {
02787       return levelbegin<cd,pitype>(maxLevel());
02788     }
02789 
02791   template<int cd, PartitionIteratorType pitype>
02792   typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const
02793     {
02794       return levelend<cd,pitype>(maxLevel());
02795     }
02796 
02798   template<int cd>
02799   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
02800     {
02801       return levelbegin<cd,All_Partition>(maxLevel());
02802     }
02803 
02805   template<int cd>
02806   typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
02807     {
02808       return levelend<cd,All_Partition>(maxLevel());
02809     }
02810 
02812   int overlapSize (int level, int codim) const
02813   {
02814         YGLI g = MultiYGrid<dim,ctype>::begin(level);
02815         return g.overlap();
02816   }
02817 
02819   int overlapSize (int codim) const
02820   {
02821         YGLI g = MultiYGrid<dim,ctype>::begin(maxLevel());
02822         return g.overlap();
02823   }
02824 
02826   int ghostSize (int level, int codim) const
02827   {
02828         return 0;
02829   }
02830 
02832   int ghostSize (int codim) const
02833   {
02834       return 0;
02835   }
02836 
02838   int size (int level, int codim) const
02839   {
02840       return sizes[level][codim];
02841   }
02842 
02844   int size (int codim) const
02845   {
02846       return sizes[maxLevel()][codim];
02847   }
02848 
02850   int size (int level, GeometryType type) const
02851   {
02852       return (type.isCube()) ? sizes[level][dim-type.dim()] : 0;
02853   }
02854 
02856   int size (GeometryType type) const
02857   {
02858       return size(maxLevel(),type);
02859   }
02860 
02862   size_t numBoundarySegments () const
02863   {
02864     return nBSegments;
02865   }
02866 
02871   template<class DataHandleImp, class DataType>
02872   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir, int level) const
02873   {
02874     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,level);
02875   }
02876 
02881   template<class DataHandleImp, class DataType>
02882   void communicate (CommDataHandleIF<DataHandleImp,DataType> & data, InterfaceType iftype, CommunicationDirection dir) const
02883   {
02884     YaspCommunicateMeta<dim,dim>::comm(*this,data,iftype,dir,this->maxLevel());
02885   }
02886 
02891   template<class DataHandle, int codim>
02892   void communicateCodim (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
02893   {
02894     // check input
02895     if (!data.contains(dim,codim)) return; // should have been checked outside
02896 
02897     // data types
02898     typedef typename DataHandle::DataType DataType;
02899 
02900     // access to grid level
02901     YGLI g = MultiYGrid<dim,ctype>::begin(level);
02902 
02903     // find send/recv lists or throw error
02904     const std::deque<IS>* sendlist=0;
02905     const std::deque<IS>* recvlist=0;
02906     if (codim==0) // the elements
02907       {
02908         if (iftype==InteriorBorder_InteriorBorder_Interface)
02909           return; // there is nothing to do in this case
02910         if (iftype==InteriorBorder_All_Interface)
02911           {
02912             sendlist = &g.send_cell_interior_overlap();
02913             recvlist = &g.recv_cell_overlap_interior();
02914           }
02915         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface || iftype==All_All_Interface)
02916           {
02917             sendlist = &g.send_cell_overlap_overlap();
02918             recvlist = &g.recv_cell_overlap_overlap();
02919           }
02920       }
02921     if (codim==dim) // the vertices
02922       {
02923         if (iftype==InteriorBorder_InteriorBorder_Interface)
02924           {
02925             sendlist = &g.send_vertex_interiorborder_interiorborder();
02926             recvlist = &g.recv_vertex_interiorborder_interiorborder();
02927           }
02928 
02929         if (iftype==InteriorBorder_All_Interface)
02930           {
02931             sendlist = &g.send_vertex_interiorborder_overlapfront();
02932             recvlist = &g.recv_vertex_overlapfront_interiorborder();
02933           }
02934         if (iftype==Overlap_OverlapFront_Interface || iftype==Overlap_All_Interface)
02935           {
02936             sendlist = &g.send_vertex_overlap_overlapfront();
02937             recvlist = &g.recv_vertex_overlapfront_overlap();
02938           }
02939         if (iftype==All_All_Interface)
02940           {
02941             sendlist = &g.send_vertex_overlapfront_overlapfront();
02942             recvlist = &g.recv_vertex_overlapfront_overlapfront();
02943           }
02944       }
02945 
02946     // change communication direction?
02947     if (dir==BackwardCommunication)
02948       std::swap(sendlist,recvlist);
02949 
02950     int cnt;
02951 
02952     // Size computation (requires communication if variable size)
02953     std::vector<int> send_size(sendlist->size(),-1);      // map rank to total number of objects (of type DataType) to be sent
02954     std::vector<int> recv_size(recvlist->size(),-1);      // map rank to total number of objects (of type DataType) to be recvd
02955     std::vector<size_t*> send_sizes(sendlist->size(),static_cast<size_t*>(0)); // map rank to array giving number of objects per entity to be sent
02956     std::vector<size_t*> recv_sizes(recvlist->size(),static_cast<size_t*>(0)); // map rank to array giving number of objects per entity to be recvd
02957     if (data.fixedsize(dim,codim))
02958       {
02959         // fixed size: just take a dummy entity, size can be computed without communication
02960         cnt=0;
02961         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02962           {
02963             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02964               it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
02965             send_size[cnt] = is->grid.totalsize() * data.size(*it);
02966             cnt++;
02967           }
02968         cnt=0;
02969         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
02970           {
02971             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02972               it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
02973             recv_size[cnt] = is->grid.totalsize() * data.size(*it);
02974             cnt++;
02975           }
02976       }
02977     else
02978       {
02979         // variable size case: sender side determines the size
02980         cnt=0;
02981         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
02982           {
02983             // allocate send buffer for sizes per entitiy
02984             size_t *buf = new size_t[is->grid.totalsize()];
02985             send_sizes[cnt] = buf;
02986 
02987             // loop over entities and ask for size
02988             int i=0; size_t n=0;
02989             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02990               it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
02991             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
02992               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubend()));
02993             for ( ; it!=tsubend; ++it)
02994               {
02995                 buf[i] = data.size(*it);
02996                 n += buf[i];
02997                 i++;
02998               }
02999 
03000             // now we know the size for this rank
03001             send_size[cnt] = n;
03002 
03003             // hand over send request to torus class
03004             MultiYGrid<dim,ctype>::torus().send(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
03005             cnt++;
03006           }
03007 
03008         // allocate recv buffers for sizes and store receive request
03009         cnt=0;
03010         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
03011           {
03012             // allocate recv buffer
03013             size_t *buf = new size_t[is->grid.totalsize()];
03014             recv_sizes[cnt] = buf;
03015             
03016             // hand over recv request to torus class
03017             MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,is->grid.totalsize()*sizeof(size_t));
03018             cnt++;
03019           }
03020 
03021         // exchange all size buffers now
03022         MultiYGrid<dim,ctype>::torus().exchange();
03023 
03024         // release send size buffers
03025         cnt=0;
03026         for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
03027           {
03028             delete[] send_sizes[cnt];
03029             send_sizes[cnt] = 0;
03030             cnt++;
03031           }
03032 
03033         // process receive size buffers
03034         cnt=0;
03035         for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
03036           {
03037             // get recv buffer
03038             size_t *buf = recv_sizes[cnt];
03039 
03040             // compute total size
03041             size_t n=0;
03042             for (int i=0; i<is->grid.totalsize(); ++i)
03043               n += buf[i];
03044 
03045             // ... and store it
03046             recv_size[cnt] = n;
03047             ++cnt;
03048           }
03049       }
03050 
03051 
03052     // allocate & fill the send buffers & store send request
03053     std::vector<DataType*> sends(sendlist->size(), static_cast<DataType*>(0)); // store pointers to send buffers
03054     cnt=0;
03055     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
03056       {
03057 //      std::cout << "[" << this->comm().rank() << "] "
03058 //                << " send " << " dest=" << is->rank
03059 //                << " size=" << send_size[cnt] 
03060 //                << std::endl;
03061 
03062         // allocate send buffer
03063         DataType *buf = new DataType[send_size[cnt]];
03064 
03065         // remember send buffer
03066         sends[cnt] = buf;
03067 
03068         // make a message buffer
03069         MessageBuffer<DataType> mb(buf);
03070 
03071         // fill send buffer; iterate over cells in intersection
03072         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03073           it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
03074         typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03075           tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubend()));
03076         for ( ; it!=tsubend; ++it)
03077           data.gather(mb,*it);
03078 
03079         // hand over send request to torus class
03080         MultiYGrid<dim,ctype>::torus().send(is->rank,buf,send_size[cnt]*sizeof(DataType));
03081         cnt++;
03082       }
03083 
03084     // allocate recv buffers and store receive request
03085     std::vector<DataType*> recvs(recvlist->size(),static_cast<DataType*>(0)); // store pointers to send buffers
03086     cnt=0;
03087     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
03088       {
03089 //      std::cout << "[" << this->comm().rank() << "] "
03090 //                << " recv " << "  src=" << is->rank
03091 //                << " size=" << recv_size[cnt] 
03092 //                << std::endl;
03093 
03094         // allocate recv buffer
03095         DataType *buf = new DataType[recv_size[cnt]];
03096 
03097         // remember recv buffer
03098         recvs[cnt] = buf;
03099 
03100         // hand over recv request to torus class
03101         MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,recv_size[cnt]*sizeof(DataType));
03102         cnt++;
03103       }
03104 
03105     // exchange all buffers now
03106     MultiYGrid<dim,ctype>::torus().exchange();
03107 
03108     // release send buffers
03109     cnt=0;
03110     for (ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
03111     {
03112       delete[] sends[cnt];
03113       sends[cnt] = 0;
03114       cnt++;
03115     }
03116 
03117     // process receive buffers and delete them
03118     cnt=0;
03119     for (ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
03120       {
03121         // get recv buffer
03122         DataType *buf = recvs[cnt];
03123 
03124         // make a message buffer
03125         MessageBuffer<DataType> mb(buf);
03126 
03127         // copy data from receive buffer; iterate over cells in intersection
03128         if (data.fixedsize(dim,codim))
03129           {
03130             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03131               it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
03132             size_t n=data.size(*it);
03133             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03134               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubend()));
03135             for ( ; it!=tsubend; ++it)
03136               data.scatter(mb,*it,n);
03137           }
03138         else
03139           {
03140             int i=0;
03141             size_t *sbuf = recv_sizes[cnt];
03142             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03143               it(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubbegin()));
03144             typename Traits::template Codim<codim>::template Partition<All_Partition>::LevelIterator
03145               tsubend(YaspLevelIterator<codim,All_Partition,GridImp>(this,g,is->grid.tsubend()));
03146             for ( ; it!=tsubend; ++it)
03147               data.scatter(mb,*it,sbuf[i++]);
03148             delete[] sbuf;
03149           }
03150 
03151         // delete buffer
03152         delete[] buf; // hier krachts !
03153         cnt++;
03154       }
03155   }
03156 
03157   // The new index sets from DDM 11.07.2005
03158   const typename Traits::GlobalIdSet& globalIdSet() const
03159   {
03160     return *(theglobalidset[0]);
03161   }
03162   
03163   const typename Traits::LocalIdSet& localIdSet() const
03164   {
03165     return *(theglobalidset[0]);
03166   }
03167 
03168   const typename Traits::LevelIndexSet& levelIndexSet(int level) const
03169   {
03170     if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
03171     return *(indexsets[level]);
03172   }
03173 
03174   const typename Traits::LeafIndexSet& leafIndexSet() const
03175   {
03176     return *(theleafindexset[0]);
03177   }
03178 
03179 #if HAVE_MPI
03180 
03182   const CollectiveCommunication<MPI_Comm>& comm () const
03183   {
03184     return ccobj;
03185   }
03186 #else
03187 
03189   const CollectiveCommunication<YaspGrid>& comm () const
03190   {
03191     return ccobj;
03192   }
03193 #endif
03194 
03195   YaspIntersectionIterator<const YaspGrid<dim> >&
03196   getRealIntersectionIterator(typename Traits::LevelIntersectionIterator& it)
03197   {
03198     return this->getRealImplementation(it);
03199   }
03200 
03201   const YaspIntersectionIterator<const YaspGrid<dim> >&
03202   getRealIntersectionIterator(const typename Traits::LevelIntersectionIterator& it) const
03203   {
03204     return this->getRealImplementation(it);
03205   }
03206 
03207 private:
03208 
03209 #if HAVE_MPI
03210   CollectiveCommunication<MPI_Comm> ccobj; 
03211 #else
03212   CollectiveCommunication<YaspGrid> ccobj; 
03213 #endif
03214 
03215   std::vector<YaspLevelIndexSet<const YaspGrid<dim> >*> indexsets;
03216   std::vector<YaspLeafIndexSet<const YaspGrid<dim> >*> theleafindexset;
03217   std::vector<YaspGlobalIdSet<const YaspGrid<dim> >*> theglobalidset;
03218   int nBSegments;
03219 
03220   // Index classes need access to the real entity
03221   friend class Dune::YaspLevelIndexSet<const Dune::YaspGrid<dim> >;
03222   friend class Dune::YaspLeafIndexSet<const Dune::YaspGrid<dim> >;
03223   friend class Dune::YaspGlobalIdSet<const Dune::YaspGrid<dim> >;
03224 
03225   friend class Dune::YaspIntersectionIterator<const Dune::YaspGrid<dim> >;
03226   friend class Dune::YaspIntersection<const Dune::YaspGrid<dim> >;
03227   friend class Dune::YaspEntity<0, dim, const Dune::YaspGrid<dim> >;
03228     
03229   template<class T>
03230   void deallocatePointers(T& container)
03231   {
03232     typedef typename T::iterator Iterator;
03233     
03234     for(Iterator entry=container.begin(); entry != container.end(); ++entry)
03235       delete (*entry);
03236   }
03237   
03238   template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
03239   friend class Entity;
03240 
03241   template<class DT>
03242   class MessageBuffer {
03243   public:
03244     // Constructor
03245     MessageBuffer (DT *p)
03246     {
03247       a=p;
03248       i=0;
03249       j=0;
03250     }
03251 
03252     // write data to message buffer, acts like a stream !
03253     template<class Y>
03254     void write (const Y& data)
03255     {
03256       dune_static_assert(( is_same<DT,Y>::value ), "DataType missmatch");
03257       a[i++] = data;
03258     }
03259     
03260     // read data from message buffer, acts like a stream !
03261     template<class Y>
03262     void read (Y& data) const 
03263     {
03264       dune_static_assert(( is_same<DT,Y>::value ), "DataType missmatch");
03265       data = a[j++];
03266     }
03267 
03268   private:
03269     DT *a;
03270     int i;
03271     mutable int j;
03272   };
03273 
03274   void setsizes ()
03275   {
03276     for (YGLI g=MultiYGrid<dim,ctype>::begin(); g!=MultiYGrid<dim,ctype>::end(); ++g)
03277       {
03278         // codim 0 (elements)
03279         sizes[g.level()][0] = 1;
03280         for (int i=0; i<dim; ++i) 
03281           sizes[g.level()][0] *= g.cell_overlap().size(i);
03282 
03283         // codim 1 (faces)
03284         if (dim>1)
03285           {
03286             sizes[g.level()][1] = 0;
03287             for (int i=0; i<dim; ++i) 
03288               {
03289                 int s=g.cell_overlap().size(i)+1;
03290                 for (int j=0; j<dim; ++j)
03291                   if (j!=i)
03292                     s *= g.cell_overlap().size(j);
03293                 sizes[g.level()][1] += s;
03294               }
03295           }
03296 
03297         // codim dim-1 (edges)
03298         if (dim>2)
03299           {
03300             sizes[g.level()][dim-1] = 0;
03301             for (int i=0; i<dim; ++i) 
03302               {
03303                 int s=g.cell_overlap().size(i);
03304                 for (int j=0; j<dim; ++j)
03305                   if (j!=i)
03306                     s *= g.cell_overlap().size(j)+1;
03307                 sizes[g.level()][dim-1] += s;
03308               }
03309           }
03310 
03311         // codim dim (vertices)
03312         sizes[g.level()][dim] = 1;
03313         for (int i=0; i<dim; ++i) 
03314           sizes[g.level()][dim] *= g.vertex_overlapfront().size(i);
03315       }
03316   }
03317 
03319   template<int cd, PartitionIteratorType pitype>
03320   YaspLevelIterator<cd,pitype,GridImp> levelbegin (int level) const
03321   {
03322         dune_static_assert( cd == dim || cd == 0 ,
03323           "YaspGrid only supports Entities with codim=dim and codim=0");
03324         YGLI g = MultiYGrid<dim,ctype>::begin(level);
03325         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
03326         if (pitype==Ghost_Partition)            
03327             return levelend <cd, pitype> (level);
03328         if (cd==0) // the elements
03329           {
03330                 if (pitype<=InteriorBorder_Partition) 
03331                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.cell_interior().tsubbegin());
03332                 if (pitype<=All_Partition)            
03333                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.cell_overlap().tsubbegin());
03334           }
03335         if (cd==dim) // the vertices
03336           {
03337                 if (pitype==Interior_Partition)       
03338                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_interior().tsubbegin());
03339                 if (pitype==InteriorBorder_Partition) 
03340                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_interiorborder().tsubbegin());
03341                 if (pitype==Overlap_Partition)        
03342                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_overlap().tsubbegin());
03343                 if (pitype<=All_Partition)            
03344                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_overlapfront().tsubbegin());
03345           }
03346         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
03347   }
03348 
03350   template<int cd, PartitionIteratorType pitype>
03351   YaspLevelIterator<cd,pitype,GridImp> levelend (int level) const
03352   {
03353         dune_static_assert( cd == dim || cd == 0 ,
03354           "YaspGrid only supports Entities with codim=dim and codim=0");
03355         YGLI g = MultiYGrid<dim,ctype>::begin(level);
03356         if (level<0 || level>maxLevel()) DUNE_THROW(RangeError, "level out of range");
03357         if (cd==0) // the elements
03358           {
03359                 if (pitype<=InteriorBorder_Partition) 
03360                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.cell_interior().tsubend());
03361                 if (pitype<=All_Partition || pitype == Ghost_Partition)
03362                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.cell_overlap().tsubend());
03363           }
03364         if (cd==dim) // the vertices
03365           {
03366                 if (pitype==Interior_Partition)       
03367                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_interior().tsubend());
03368                 if (pitype==InteriorBorder_Partition) 
03369                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_interiorborder().tsubend());
03370                 if (pitype==Overlap_Partition)        
03371                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_overlap().tsubend());
03372                 if (pitype<=All_Partition || pitype == Ghost_Partition)
03373                   return YaspLevelIterator<cd,pitype,GridImp>(this,g,g.vertex_overlapfront().tsubend());
03374           }
03375         DUNE_THROW(GridError, "YaspLevelIterator with this codim or partition type not implemented");
03376   }
03377 
03378   int sizes[MAXL][dim+1]; // total number of entities per level and codim
03379   bool keep_ovlp;
03380   int adaptRefCount;
03381   bool adaptActive;
03382 };
03383 
03384 namespace Capabilities
03385 {
03386 
03398   template<int dim>
03399   struct hasSingleGeometryType< YaspGrid<dim> >
03400   {
03401     static const bool v = true;
03402     static const unsigned int topologyId = GenericGeometry :: CubeTopology< dim > :: type :: id ;
03403   };
03404 
03408   template<int dim>
03409   struct isCartesian< YaspGrid<dim> >
03410   {
03411     static const bool v = true;
03412   };
03413 
03417   template<int dim>
03418   struct hasEntity< YaspGrid<dim>, 0 >
03419   {
03420     static const bool v = true;
03421   };
03422   
03426   template<int dim>
03427   struct hasEntity< YaspGrid<dim>, dim >
03428   {
03429     static const bool v = true;
03430   };
03431 
03432   template< int dim >
03433   struct canCommunicate< YaspGrid< dim >, 0 >
03434   {
03435     static const bool v = true;
03436   };
03437 
03438   template< int dim >
03439   struct canCommunicate< YaspGrid< dim >, dim >
03440   {
03441     static const bool v = true;
03442   };
03443 
03447   template<int dim>
03448   struct isParallel< YaspGrid<dim> >
03449   {
03450     static const bool v = true;
03451   };
03452 
03456   template<int dim>
03457   struct isLevelwiseConforming< YaspGrid<dim> >
03458   {
03459     static const bool v = true;
03460   };
03461 
03465   template<int dim>
03466   struct isLeafwiseConforming< YaspGrid<dim> >
03467   {
03468     static const bool v = true;
03469   };
03470 
03471 }
03472 
03473 } // end namespace
03474 
03475 
03476 #endif
03477 

Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].