grids.hh

Go to the documentation of this file.
00001 #ifndef DUNE_YGRIDS_HH
00002 #define DUNE_YGRIDS_HH
00003 
00004 // C++ includes
00005 #include<iostream>
00006 #include<cstdlib>
00007 #include<algorithm>
00008 #include<vector>
00009 #include<deque>
00010 
00011 // C includes
00012 #if HAVE_MPI
00013 #include<mpi.h>
00014 #endif
00015 #include<string.h>
00016 
00017 // local includes
00018 #include <dune/common/fvector.hh>
00019 #include <dune/common/stdstreams.hh>
00020 
00025 namespace Dune {
00026 
00027   // forward declarations
00028   template<int d, typename ct> class YGrid;
00029   template<int d, typename ct> class SubYGrid;
00030 
00031   static const double Ytolerance=1E-13;
00032  
00066   template<int d, typename ct>
00067   class YGrid {
00068   public:
00070         typedef FieldVector<int, d>  iTupel;
00071         typedef FieldVector<ct, d> fTupel;
00072         typedef FieldVector<bool, d> bTupel;
00073     
00075     virtual ~YGrid()
00076     {}
00077     
00079         YGrid ()
00080         {
00081           _origin = 0;
00082           _size = 0;
00083           _h = 0.0;
00084           _r = 0.0;
00085         }
00086 
00088         YGrid (iTupel o, iTupel s, fTupel h, fTupel r)
00089         {
00090           for (int i=0; i<d; ++i) 
00091                 {
00092                   _origin[i] = o[i];
00093                   _size[i] = s[i];
00094                   if (_size[i]<0) {
00095                         _size[i] = 0;
00096                   }
00097                   _h[i] = h[i];
00098                   _r[i] = r[i];
00099                 }
00100         }
00101 
00103         int origin (int i) const
00104         {
00105           return _origin[i];
00106         }
00107         
00109         void origin (int i, int oi) const
00110         {
00111           _origin[i] = oi;
00112         }
00113 
00115         const iTupel& origin () const
00116         {
00117           return _origin;
00118         }
00119 
00121         int size (int i) const
00122         {
00123           return _size[i];
00124         }
00125 
00127         void size (int i, int si) const
00128         {
00129           _size[i] = si;
00130           if (_size[i]<0) _size[i] = 0;
00131         }
00132         
00134         const iTupel& size () const
00135         {
00136           return _size;
00137         }
00138 
00140         int totalsize () const
00141         {
00142           int s=1;
00143           for (int i=0; i<d; ++i) s=s*_size[i];
00144           return s;
00145         }
00146 
00148         int min (int i) const
00149         {
00150           return _origin[i];
00151         }
00152 
00154         void min (int i, int mi) const
00155         {
00156           _size[i] = max(i)-mi+1;
00157           _origin[i] = mi;
00158           if (_size[i]<0) _size[i] = 0;
00159         }
00160 
00162         int max (int i) const
00163         {
00164           return _origin[i]+_size[i]-1;
00165         }
00166 
00168         void max (int i, int mi) const
00169         {
00170           _size[i] = mi-min(i)+1;
00171           if (_size[i]<0) _size[i] = 0;
00172         }
00173 
00175         const fTupel& meshsize () const
00176         {
00177           return _h;
00178         }
00179         
00181         ct meshsize (int i) const
00182         {
00183           return _h[i];
00184         }
00185         
00187         void meshsize (int i, int hi) const
00188         {
00189           _h[i] = hi;
00190         }
00191 
00193         const fTupel& shift () const
00194         {
00195           return _r;
00196         }
00197         
00199         ct shift (int i) const
00200         {
00201           return _r[i];
00202         }
00203         
00205         void shift (int i, int ri) const
00206         {
00207           _r[i] = ri;
00208         }
00209 
00211         bool empty () const
00212         {
00213           for (int i=0; i<d; ++i) if (_size[i]<=0) return true;
00214           return false;
00215         }
00216 
00218         int index (const iTupel& coord) const
00219         {
00220           int index = (coord[d-1]-_origin[d-1]);
00221           
00222           for (int i=d-2; i>=0; i--) 
00223                 index = index*_size[i] + (coord[i]-_origin[i]);
00224 
00225           return index;
00226         }
00227 
00229         bool inside (const iTupel& coord) const
00230         {
00231           for (int i=0; i<d; i++)
00232                 {
00233                   if (coord[i]<_origin[i] || coord[i]>=_origin[i]+_size[i]) return false;
00234                 }
00235           return true;
00236         }
00237 
00239         virtual SubYGrid<d,ct> intersection ( YGrid<d,ct>& r) const
00240         {
00241           // check if the two grids can be intersected, must have same mesh size and shift
00242           for (int i=0; i<d; i++)
00243                 if (fabs(meshsize(i)-r.meshsize(i))>Ytolerance) return SubYGrid<d,ct>();
00244           for (int i=0; i<d; i++)
00245                 if (fabs(shift(i)-r.shift(i))>Ytolerance) return SubYGrid<d,ct>();
00246 
00247           iTupel neworigin;
00248           iTupel newsize;
00249           iTupel offset;
00250 
00251           for (int i=0; i<d; ++i)
00252                 {
00253                   // intersect
00254                   neworigin[i] = std::max(min(i),r.min(i));
00255                   newsize[i] = std::min(max(i),r.max(i))-neworigin[i]+1;
00256                   if (newsize[i]<0) {
00257                         newsize[i] = 0;
00258                         neworigin[i] = min(i);
00259                   }
00260 
00261                   // offset to own origin
00262                   offset[i] = neworigin[i]-_origin[i];
00263                 }       
00264           return SubYGrid<d,ct>(neworigin,newsize,offset,_size,_h,_r);
00265         }
00266 
00268         YGrid<d,ct> move (iTupel v) const
00269         {
00270           for (int i=0; i<d; i++) v[i] += _origin[i];
00271           return YGrid<d,ct>(v,_size,_h,_r);
00272         }
00273 
00280         class Iterator {
00281         public:
00283           Iterator (const YGrid<d,ct>& r)
00284           {
00285                 // copy data coming from grid to iterate over 
00286                 for (int i=0; i<d; ++i) _origin[i] = r.origin(i);
00287                 for (int i=0; i<d; ++i) _end[i] = r.origin(i)+r.size(i)-1;
00288                 
00289                 // initialize to first position in index set
00290                 for (int i=0; i<d; ++i) _coord[i] = _origin[i];
00291                 _index = 0;
00292                 
00293                 // compute increments;
00294                 int inc = 1;
00295                 for (int i=0; i<d; ++i)
00296                   {
00297                         _increment[i] = inc;
00298                         inc *= r.size(i);
00299                   }
00300           }
00301 
00303           Iterator (const YGrid<d,ct>& r, const iTupel& coord)
00304           {
00305                 // copy data coming from grid to iterate over 
00306                 for (int i=0; i<d; ++i) _origin[i] = r.origin(i);
00307                 for (int i=0; i<d; ++i) _end[i] = r.origin(i)+r.size(i)-1;
00308                 
00309                 // compute increments;
00310                 int inc = 1;
00311                 for (int i=0; i<d; ++i)
00312                   {
00313                         _increment[i] = inc;
00314                         inc *= r.size(i);
00315                   }
00316 
00317                 // initialize to given position in index set
00318                 for (int i=0; i<d; ++i) _coord[i] = coord[i];
00319                 _index = r.index(coord);
00320           }
00321 
00323           void reinit (const YGrid<d,ct>& r, const iTupel& coord)
00324           {
00325                 // copy data coming from grid to iterate over 
00326                 for (int i=0; i<d; ++i) _origin[i] = r.origin(i);
00327                 for (int i=0; i<d; ++i) _end[i] = r.origin(i)+r.size(i)-1;
00328                 
00329                 // compute increments;
00330                 int inc = 1;
00331                 for (int i=0; i<d; ++i)
00332                   {
00333                         _increment[i] = inc;
00334                         inc *= r.size(i);
00335                   }
00336 
00337                 // initialize to given position in index set
00338                 for (int i=0; i<d; ++i) _coord[i] = coord[i];
00339                 _index = r.index(coord);
00340           }
00341 
00343           bool operator== (const Iterator& i) const
00344           {
00345                 return _index == i._index;
00346           }
00347 
00349           bool operator!= (const Iterator& i) const
00350           {
00351                 return _index != i._index;
00352           }
00353         
00355           int index () const
00356           {
00357                 return _index;
00358           }
00359           
00361           int coord (int i) const
00362           {
00363                 return _coord[i];
00364           }
00365 
00367           const iTupel& coord () const
00368           {
00369                 return _coord;
00370           }
00371 
00373           int neighbor (int i, int dist) const
00374           {
00375                 return _index+dist*_increment[i];
00376           }
00377 
00379           int down (int i) const
00380           {
00381                 return _index-_increment[i];
00382           }
00383 
00385           int up (int i) const
00386           {
00387                 return _index+_increment[i];
00388           }
00389 
00391           void move (int i, int dist)
00392           {
00393                 _coord[i] += dist;
00394                 _index += dist*_increment[i];
00395           }
00396 
00398           Iterator& operator++ ()
00399           {
00400                 ++_index;
00401                 for (int i=0; i<d; i++)
00402                   if (++(_coord[i])<=_end[i])
00403                         return *this;
00404                   else { _coord[i]=_origin[i]; }
00405                 return *this;
00406           }
00407           
00409           void print (std::ostream& s) const
00410           {
00411                 s << index() << " : [";
00412                 for (int i=0; i<d-1; i++) s << coord(i) << ",";
00413                 s << coord(d-1) << "]";
00414           }
00415 
00416         protected:
00417           int _index;          //< current lexicographic position in index set
00418           iTupel _coord;       //< current position in index set 
00419           iTupel _increment;   //< increment for next neighbor in direction i 
00420           iTupel _origin;      //< origin and
00421           iTupel _end;         //< last index in direction i
00422         };
00423 
00425         Iterator begin () const { return Iterator(*this); }
00426         
00428         Iterator end () const {
00429           iTupel last;
00430           for (int i=0; i<d; i++) last[i] = max(i);
00431           last[0] += 1;
00432           return Iterator(*this,last);
00433         }
00434 
00439         class TransformingIterator : public Iterator {
00440         public:
00442           TransformingIterator (const YGrid<d,ct>& r) : Iterator(r)
00443           {
00444                 for (int i=0; i<d; ++i) _h[i] = r.meshsize(i);
00445                 for (int i=0; i<d; ++i) _begin[i] = r.origin(i)*r.meshsize(i)+r.shift(i);
00446                 for (int i=0; i<d; ++i) _position[i] = _begin[i];
00447           }
00448 
00450           TransformingIterator (const YGrid<d,ct>& r, iTupel& coord) : Iterator(r,coord)
00451           {
00452                 for (int i=0; i<d; ++i) _h[i] = r.meshsize(i);
00453                 for (int i=0; i<d; ++i) _begin[i] = r.origin(i)*r.meshsize(i)+r.shift(i);
00454                 for (int i=0; i<d; ++i) _position[i] = coord[i]*r.meshsize(i)+r.shift(i);
00455           }
00456 
00458           TransformingIterator (Iterator i) : Iterator(i)
00459           {       }
00460 
00462           TransformingIterator& operator++ ()
00463           {
00464                 ++(this->_index);
00465                 for (int i=0; i<d; i++)
00466                   if (++(this->_coord[i])<=this->_end[i])
00467                         {
00468                           _position[i] += _h[i];
00469                           return *this;
00470                         }
00471                   else 
00472                         { 
00473                           this->_coord[i]=this->_origin[i]; 
00474                           _position[i] = _begin[i];
00475                         }
00476                 return *this;
00477           }
00478           
00480           ct position (int i) const
00481           {
00482                 return _position[i];
00483           }
00484 
00486           const fTupel& position () const
00487           {
00488                 return _position;
00489           }
00490 
00492           ct meshsize (int i) const
00493           {
00494                 return _h[i];
00495           }
00496 
00498           const fTupel& meshsize () const
00499           {
00500                 return _h;
00501           }
00502 
00504           void move (int i, int dist)
00505           {
00506                 Iterator::move(i,dist);
00507                 _position[i] += dist*_h[i];
00508           }
00509 
00511           void print (std::ostream& s) const
00512           {
00513                 Iterator::print(s);
00514                 s << " " << _position;
00515           }
00516 
00517         private:
00518           fTupel _h;        
00519           fTupel _begin;    
00520           fTupel _position; 
00521         };
00522 
00524         TransformingIterator tbegin () const
00525         {
00526           return TransformingIterator(*this);
00527         }
00528         
00530         TransformingIterator tend () const
00531         {
00532           iTupel last;
00533           for (int i=0; i<d; i++) last = max(i);
00534           last[0] += 1;
00535           return TransformingIterator(*this,last);
00536         }
00537         
00538   protected:
00540         iTupel _origin; 
00541         iTupel _size;
00542         fTupel _h;        
00543         fTupel _r;        
00544   };
00545 
00547   template <int d, typename ct>
00548   inline std::ostream& operator<< (std::ostream& s, YGrid<d,ct> e)
00549   {
00550         s << "{";
00551         for (int i=0; i<d-1; i++)
00552           s << "[" << e.min(i) << "," << e.max(i) << "]x";
00553         s << "[" << e.min(d-1) << "," << e.max(d-1) << "]";
00554         s << " = [";
00555         for (int i=0; i<d-1; i++) s << e.origin(i) << ",";
00556         s << e.origin(d-1) << "]x[";
00557         for (int i=0; i<d-1; i++) s << e.size(i) << ",";
00558         s << e.size(d-1) << "]";
00559         s << " h=[";
00560         for (int i=0; i<d-1; i++) s << e.meshsize(i) << ",";
00561         s << e.meshsize(d-1) << "]";
00562         s << " r=[";
00563         for (int i=0; i<d-1; i++) s << e.shift(i) << ",";
00564         s << e.shift(d-1) << "]";
00565         s << "}";
00566         return s;
00567   }
00568 
00570   template <int d, typename ct>
00571   inline std::ostream& operator<< (std::ostream& s, typename YGrid<d,ct>::Iterator& e)
00572   {
00573         e.print(s);
00574         return s;
00575   }
00576 
00577 
00588   template<int d, typename ct>
00589   class SubYGrid : public YGrid<d,ct> {
00590   public:
00591         typedef typename YGrid<d,ct>::iTupel iTupel;
00592         typedef typename YGrid<d,ct>::fTupel fTupel;
00593         typedef typename YGrid<d,ct>::bTupel bTupel;
00594 
00596     virtual ~SubYGrid()
00597     {}
00598     
00600         SubYGrid () {}
00601 
00603         SubYGrid (iTupel origin, iTupel size, iTupel offset, iTupel supersize, fTupel h, fTupel r) 
00604           : YGrid<d,ct>::YGrid(origin,size,h,r)
00605         {
00606           for (int i=0; i<d; ++i) 
00607                 {
00608                   _offset[i] = offset[i];
00609                   _supersize[i] = supersize[i];
00610                   if (offset[i]<0) 
00611                         std::cout << "warning: offset["
00612                                           << i <<"] negative in SubYGrid"
00613                                           << std::endl; 
00614                   if (-offset[i]+supersize[i]<size[i]) 
00615                         std::cout << "warning: subgrid larger than enclosing grid in direction "
00616                                           << i <<" in SubYGrid"
00617                                           << std::endl; 
00618                 }
00619         }
00620 
00622         SubYGrid (YGrid<d,ct> base) : YGrid<d,ct>(base)
00623         {
00624           for (int i=0; i<d; ++i) 
00625                 {
00626                   _offset[i] = 0;
00627                   _supersize[i] = this->size(i);
00628                 }
00629         }
00630  
00632         int offset (int i) const
00633         {
00634           return _offset[i];
00635         }
00636 
00638         const iTupel & offset () const
00639         {
00640           return _offset;
00641         }
00642 
00644         int supersize (int i) const
00645         {
00646           return _supersize[i];
00647         }
00648 
00650         const iTupel & supersize () const
00651         {
00652           return _supersize;
00653         }
00654 
00656         virtual SubYGrid<d,ct> intersection (const YGrid<d,ct>& r) const
00657         {
00658           // check if the two grids can be intersected, must have same mesh size and shift
00659           for (int i=0; i<d; i++)
00660                 if (fabs(this->meshsize(i)-r.meshsize(i))>Ytolerance) return SubYGrid<d,ct>();
00661           for (int i=0; i<d; i++)
00662                 if (fabs(this->shift(i)-r.shift(i))>Ytolerance) return SubYGrid<d,ct>();
00663 
00664           iTupel neworigin;
00665           iTupel newsize;
00666           iTupel offset;
00667 
00668           for (int i=0; i<d; ++i)
00669                 {
00670                   // intersect
00671                   neworigin[i] = std::max(this->min(i),r.min(i));
00672                   newsize[i] = std::min(this->max(i),r.max(i))-neworigin[i]+1;
00673                   if (newsize[i]<0) {
00674                         newsize[i] = 0;
00675                         neworigin[i] = this->min(i);
00676                   }
00677 
00678                   // offset to my supergrid
00679                   offset[i] = _offset[i]+neworigin[i]-this->origin(i);
00680                 }       
00681           return SubYGrid<d,ct>(neworigin,newsize,offset,_supersize,this->meshsize(),this->shift());
00682         }
00683 
00687         class SubIterator : public YGrid<d,ct>::Iterator {
00688         public:
00690           SubIterator (const SubYGrid<d,ct>& r) : YGrid<d,ct>::Iterator::Iterator (r)
00691           {
00693                 for (int i=0; i<d; ++i) _size[i] = r.size(i);
00694 
00695                 // compute superincrements 
00696                 int inc = 1;
00697                 for (int i=0; i<d; ++i)
00698                   {
00699                         _superincrement[i] = inc;
00700                         inc *= r.supersize(i);
00701                   }
00702 
00703                 // move superindex to first cell in subgrid
00704                 _superindex = 0;
00705                 for (int i=0; i<d; ++i)
00706                   _superindex += r.offset(i)*_superincrement[i];
00707           }
00708 
00710           SubIterator (const SubYGrid<d,ct>& r, const iTupel& coord) : YGrid<d,ct>::Iterator::Iterator (r,coord)
00711           {
00713                 for (int i=0; i<d; ++i) _size[i] = r.size(i);
00714 
00715                 // compute superincrements 
00716                 int inc = 1;
00717                 for (int i=0; i<d; ++i)
00718                   {
00719                         _superincrement[i] = inc;
00720                         inc *= r.supersize(i);
00721                   }
00722 
00723                 // move superindex to first cell in subgrid
00724                 _superindex = 0;
00725                 for (int i=0; i<d; ++i)
00726                   _superindex += (r.offset(i)+coord[i]-r.origin(i))*_superincrement[i];
00727           }
00728 
00730           SubIterator (const typename YGrid<d,ct>::Iterator& i) : YGrid<d,ct>::Iterator::Iterator(i)
00731         {}
00732 
00734           void reinit (const SubYGrid<d,ct>& r, const iTupel& coord) 
00735           {
00736                 YGrid<d,ct>::Iterator::reinit(r,coord);
00737 
00739                 for (int i=0; i<d; ++i) _size[i] = r.size(i);
00740 
00741                 // compute superincrements 
00742                 int inc = 1;
00743                 for (int i=0; i<d; ++i)
00744                   {
00745                         _superincrement[i] = inc;
00746                         inc *= r.supersize(i);
00747                   }
00748 
00749                 // move superindex to first cell in subgrid
00750                 _superindex = 0;
00751                 for (int i=0; i<d; ++i)
00752                   _superindex += (r.offset(i)+coord[i]-r.origin(i))*_superincrement[i];
00753           }
00754 
00756           int superindex () const
00757           {
00758                 return _superindex;
00759           }
00760 
00762           int superneighbor (int i, int dist) const
00763           {
00764                 return _superindex+dist*_superincrement[i];
00765           }
00766 
00768           int superdown (int i) const
00769           {
00770                 return _superindex-_superincrement[i];
00771           }
00772 
00774           int superup (int i) const
00775           {
00776                 return _superindex+_superincrement[i];
00777           }
00778 
00780           void move (int i, int dist)
00781           {
00782                 YGrid<d,ct>::Iterator::move(i,dist);    // move base iterator
00783                 _superindex += dist*_superincrement[i]; // move superindex
00784           }
00785 
00787           SubIterator& operator++ ()
00788           {
00789                 ++(this->_index);               // update consecutive index in grid
00790                 for (int i=0; i<d; i++)         // check for wrap around
00791                   {
00792                         _superindex += _superincrement[i]; // move on cell in direction i
00793                         if (++(this->_coord[i])<=this->_end[i])
00794                           return *this;
00795                         else 
00796                           { 
00797                                 this->_coord[i]=this->_origin[i];     // move back to origin in direction i
00798                                 _superindex -= _size[i]*_superincrement[i]; 
00799                           }
00800                   }
00801                 return *this;
00802           }
00803           
00805           void print (std::ostream& s) const
00806           {
00807                 YGrid<d,ct>::Iterator::print(s);
00808                 s << " super=" << superindex();
00809           }
00810 
00811         protected:
00812           int _superindex;        
00813           iTupel _superincrement; 
00814           iTupel _size;           
00815         };
00816 
00818         SubIterator subbegin () const { return SubIterator(*this); }
00819         
00821         SubIterator subend () const
00822         {
00823           iTupel last;
00824           for (int i=0; i<d; i++) last[i] = this->max(i);
00825           last[0] += 1;
00826           return SubIterator(*this,last);
00827         }
00828 
00833         class TransformingSubIterator : public SubIterator {
00834         public:
00836           TransformingSubIterator (const SubYGrid<d,ct>& r) : SubIterator(r)
00837           {
00838                 for (int i=0; i<d; ++i) _h[i] = r.meshsize(i);
00839                 for (int i=0; i<d; ++i) _begin[i] = r.origin(i)*r.meshsize(i)+r.shift(i);
00840                 for (int i=0; i<d; ++i) _position[i] = _begin[i];
00841           }
00842 
00844           TransformingSubIterator (const SubYGrid<d,ct>& r, const iTupel& coord) : SubIterator(r,coord)
00845           {
00846                 for (int i=0; i<d; ++i) _h[i] = r.meshsize(i);
00847                 for (int i=0; i<d; ++i) _begin[i] = r.origin(i)*r.meshsize(i)+r.shift(i);
00848                 for (int i=0; i<d; ++i) _position[i] = coord[i]*r.meshsize(i)+r.shift(i);
00849           }
00850 
00852       TransformingSubIterator (const SubIterator& i) :
00853         SubIterator(i)
00854         {}
00855 
00856       TransformingSubIterator (const TransformingSubIterator & t) :
00857         SubIterator(t), _h(t._h), _begin(t._begin), _position(t._position)
00858         {}
00859 
00861           void reinit (const SubYGrid<d,ct>& r, const iTupel& coord) 
00862           {
00863                 SubIterator::reinit(r,coord);
00864                 for (int i=0; i<d; ++i) _h[i] = r.meshsize(i);
00865                 for (int i=0; i<d; ++i) _begin[i] = r.origin(i)*r.meshsize(i)+r.shift(i);
00866                 for (int i=0; i<d; ++i) _position[i] = coord[i]*r.meshsize(i)+r.shift(i);
00867           }
00868 
00870           TransformingSubIterator& operator++ ()
00871           {
00872                 ++(this->_index);                       // update consecutive index in subgrid
00873                 for (int i=0; i<d; i++)         // check for wrap around
00874                   {
00875                         this->_superindex += this->_superincrement[i]; // move on cell in direction i
00876                         if (++(this->_coord[i])<=this->_end[i])
00877                           {
00878                                 _position[i] += _h[i];
00879                                 return *this;
00880                           }
00881                         else 
00882                           { 
00883                                 this->_coord[i]=this->_origin[i];     // move back to origin in direction i
00884                                 this->_superindex -= this->_size[i]*this->_superincrement[i]; 
00885                                 _position[i] = _begin[i];
00886                           }
00887                   }
00888                 return *this;
00889           }
00890 
00892           ct position (int i) const
00893           {
00894                 return _position[i];
00895           }
00896 
00898           const fTupel& position () const
00899           {
00900                 return _position;
00901           }
00902 
00904           ct meshsize (int i) const
00905           {
00906                 return _h[i];
00907           }
00908 
00910           const fTupel& meshsize () const
00911           {
00912                 return _h;
00913           }
00914 
00916           void move (int i, int dist)
00917           {
00918                 SubIterator::move(i,dist);
00919                 _position[i] += dist*_h[i];
00920           }
00921 
00923           void print (std::ostream& s) const
00924           {
00925                 SubIterator::print(s);
00926                 s << " [";
00927                 for (int i=0; i<d-1; i++) s << position(i) << ",";
00928                 s << position(d-1) << "]";
00929           }
00930 
00931         private:
00932           fTupel _h;        
00933           fTupel _begin;    
00934           fTupel _position; 
00935         };
00936 
00938         TransformingSubIterator tsubbegin () const
00939         {
00940           return TransformingSubIterator(*this);
00941         }
00942         
00944         TransformingSubIterator tsubbegin (iTupel& co) const
00945         {
00946           return TransformingSubIterator(*this,co);
00947         }
00948         
00950         TransformingSubIterator tsubend () const
00951         {
00952           SubIterator endit = subend();
00953           return TransformingSubIterator(endit);
00954         }
00955         
00956   private:
00957         iTupel _offset;    
00958         iTupel _supersize; 
00959   };
00960 
00961 
00963   template <int d, typename ct>
00964   inline std::ostream& operator<< (std::ostream& s, SubYGrid<d,ct> e)
00965   {
00966         YGrid<d,ct> x = e;
00967         s << x << " ofs=" << e.offset() << " ss=" << e.supersize();
00968         return s;
00969   }
00970 
00972   template <int d, typename ct>
00973   inline std::ostream& operator<< (std::ostream& s, typename SubYGrid<d,ct>::TransformingSubIterator& e)
00974   {
00975         e.print(s);
00976         return s;
00977   }
00978 
00979 
00993   template<int d>
00994   class Torus {
00995   public:
00997         typedef FieldVector<int, d> iTupel;
00998         typedef FieldVector<bool, d> bTupel;
00999 
01000 
01001   private:
01002         struct CommPartner {
01003           int rank;
01004           iTupel delta;
01005           int index;
01006         };
01007 
01008         struct CommTask {
01009           int rank;      // process to send to / receive from
01010           void *buffer;  // buffer to send / receive
01011           int size;      // size of buffer
01012 #if HAVE_MPI
01013           MPI_Request request; // used by MPI to handle request
01014 #else
01015           int request;
01016 #endif
01017           int flag;      // used by MPI 
01018         };
01019 
01020   public:
01022         Torus () 
01023         {  }
01024 
01026 #if HAVE_MPI
01027         Torus (MPI_Comm comm, int tag, iTupel size)
01028 #else
01029         Torus (int tag, iTupel size)
01030 #endif
01031         {
01032           // MPI stuff
01033 #if HAVE_MPI
01034           _comm = comm;
01035           MPI_Comm_size(comm,&_procs);
01036           MPI_Comm_rank(comm,&_rank);
01037 #else
01038           _procs=1; _rank=0;
01039 #endif
01040           _tag = tag;
01041 
01042           // determine dimensions
01043           iTupel dims;
01044           double opt=1E100;
01045           optimize_dims(d-1,size,_procs,dims,opt);
01046 //        if (_rank==0) std::cout << "Torus<" << d 
01047 //                                                        << ">: mapping " << _procs << " processes onto "
01048 //                                                        << _dims << " torus." << std::endl;
01049 
01050           // compute increments for lexicographic ordering
01051           int inc = 1;
01052           for (int i=0; i<d; i++)
01053                 {
01054                   _increment[i] = inc;
01055                   inc *= _dims[i];
01056                 }
01057 
01058           // make full schedule
01059           proclists();
01060         }
01061 
01063         int rank () const
01064         {
01065           return _rank;
01066         }
01067 
01069         iTupel coord () const
01070         {
01071           return rank_to_coord(_rank);
01072         }
01073 
01075         int procs () const
01076         {
01077           return _procs;
01078         }
01079 
01081         const iTupel & dims () const
01082         {
01083           return _dims;
01084         }
01085 
01087         int dims (int i) const
01088         {
01089           return _dims[i];
01090         }
01091 
01093 #if HAVE_MPI
01094         MPI_Comm comm () const
01095         {
01096           return _comm;
01097         }
01098 #endif
01099 
01101         int tag () const
01102         {
01103           return _tag;
01104         }
01105 
01107         bool inside (iTupel c) const
01108         {
01109           for (int i=d-1; i>=0; i--)
01110                 if (c[i]<0 || c[i]>=_dims[i]) return false;
01111           return true;
01112         }
01113 
01115         iTupel rank_to_coord (int rank) const
01116         {
01117           iTupel coord;
01118           rank = rank%_procs;
01119           for (int i=d-1; i>=0; i--)
01120                 {
01121                   coord[i] = rank/_increment[i];
01122                   rank = rank%_increment[i];
01123                 }
01124           return coord;
01125         }
01126 
01128         int coord_to_rank (iTupel coord) const
01129         {
01130           for (int i=0; i<d; i++) coord[i] = coord[i]%_dims[i];
01131           int rank = 0;
01132           for (int i=0; i<d; i++) rank += coord[i]*_increment[i];
01133           return rank;
01134         }
01135 
01137         int rank_relative (int rank, int dir, int cnt) const
01138         {
01139           iTupel coord = rank_to_coord(rank);
01140           coord[dir] = (coord[dir]+dims[dir]+cnt)%dims[dir];
01141           return coord_to_rank(coord);
01142         }
01143 
01145         int color (const iTupel & coord) const
01146         {
01147           int c = 0;
01148           int power = 1;
01149 
01150           // interior coloring
01151           for (int i=0; i<d; i++)
01152                 {
01153                   if (coord[i]%2==1) c += power;
01154                   power *= 2;
01155                 }
01156 
01157           // extra colors for boundary processes
01158           for (int i=0; i<d; i++)
01159                 {
01160                   if (_dims[i]>1 && coord[i]==_dims[i]-1) c += power;
01161                   power *= 2;
01162                 }
01163 
01164           return c;
01165         }
01166 
01168         int color (int rank) const
01169         {
01170           return color(rank_to_coord(rank));
01171         }
01172 
01174         int neighbors () const
01175         {
01176           int n=1;
01177           for (int i=0; i<d; ++i)
01178                 n *= 3;
01179           return n-1;
01180         }
01181 
01183         bool is_neighbor (iTupel delta, bTupel periodic) const
01184         {
01185           iTupel coord = rank_to_coord(_rank); // my own coordinate with 0 <= c_i < dims_i
01186 
01187 
01188           for (int i=0; i<d; i++)
01189                 {
01190                   if (delta[i]<0)
01191                         {
01192                           // if I am on the boundary and domain is not periodic => no neighbor
01193                           if (coord[i]==0 && periodic[i]==false) return false;
01194                         }
01195                   if (delta[i]>0)
01196                         {
01197                           // if I am on the boundary and domain is not periodic => no neighbor
01198                           if (coord[i]==_dims[i]-1 && periodic[i]==false) return false;
01199                         }
01200                 }
01201           return true;
01202         }
01203 
01205         double partition (int rank, iTupel origin_in, iTupel size_in, iTupel& origin_out, iTupel& size_out) const
01206         {
01207           iTupel coord = rank_to_coord(rank);
01208           double maxsize = 1;
01209           double sz = 1;
01210 
01211           // make a tensor product partition
01212           for (int i=0; i<d; i++)
01213                 {
01214                   // determine 
01215                   int m = size_in[i]/_dims[i];
01216                   int r = size_in[i]%_dims[i];
01217 
01218                   sz *= size_in[i];
01219 
01220                   if (coord[i]<_dims[i]-r)
01221                         {
01222                           origin_out[i] = origin_in[i] + coord[i]*m;
01223                           size_out[i] = m;
01224                           maxsize *= m;
01225                         }
01226                   else
01227                         {
01228                           origin_out[i] = origin_in[i] + (_dims[i]-r)*m + (coord[i]-(_dims[i]-r))*(m+1);
01229                           size_out[i] = m+1;
01230                           maxsize *= m+1;
01231                         }
01232                 }
01233           return maxsize/(sz/_procs);
01234         }
01235 
01242         class ProcListIterator {
01243         public:
01245           ProcListIterator (typename std::deque<CommPartner>::const_iterator iter)
01246           {
01247                 i = iter;
01248           }
01249 
01251           int rank () const
01252           {
01253                 return i->rank;
01254           }
01255 
01257           iTupel delta () const
01258           {
01259                 return i->delta;
01260           }
01261 
01263           int index () const
01264           {
01265                 return i->index;
01266           }
01267 
01269           int distance () const
01270           {
01271                 int dist = 0;
01272                 iTupel delta=i->delta;
01273                 for (int i=0; i<d; ++i)
01274                   dist += std::abs(delta[i]);
01275                 return dist;
01276           }
01277 
01279           bool operator== (const ProcListIterator& iter)
01280           {
01281                 return i == iter.i;
01282           }
01283 
01284 
01286           bool operator!= (const ProcListIterator& iter)
01287           {
01288                 return i != iter.i;
01289           }
01290         
01292           ProcListIterator& operator++ ()
01293           {
01294                 ++i;
01295                 return *this;
01296           }
01297           
01298         private:
01299           typename std::deque<CommPartner>::const_iterator i;
01300         };
01301 
01303         ProcListIterator sendbegin () const
01304         {
01305           return ProcListIterator(_sendlist.begin());
01306         }
01307 
01309         ProcListIterator sendend () const
01310         {
01311           return ProcListIterator(_sendlist.end());
01312         }
01313 
01315         ProcListIterator recvbegin () const
01316         {
01317           return ProcListIterator(_recvlist.begin());
01318         }
01319 
01321         ProcListIterator recvend () const
01322         {
01323           return ProcListIterator(_recvlist.end());
01324         }
01325 
01327         void send (int rank, void* buffer, int size) const
01328         {
01329           CommTask task;
01330           task.rank = rank;
01331           task.buffer = buffer;
01332           task.size = size;
01333           if (rank!=_rank)
01334                 _sendrequests.push_back(task);
01335           else
01336                 _localsendrequests.push_back(task);
01337         }
01338  
01340         void recv (int rank, void* buffer, int size) const
01341         {
01342           CommTask task;
01343           task.rank = rank;
01344           task.buffer = buffer;
01345           task.size = size;
01346           if (rank!=_rank)
01347                 _recvrequests.push_back(task);
01348           else
01349                 _localrecvrequests.push_back(task);
01350         }
01351  
01353         void exchange () const
01354         {
01355           // handle local requests first
01356           if (_localsendrequests.size()!=_localrecvrequests.size())
01357                 {
01358                   std::cout << "[" << rank() << "]: ERROR: local sends/receives do not match in exchange!" << std::endl;
01359                   return;
01360                 }
01361           for (unsigned int i=0; i<_localsendrequests.size(); i++)
01362                 {
01363                   if (_localsendrequests[i].size!=_localrecvrequests[i].size)
01364                         {
01365                           std::cout << "[" << rank() << "]: ERROR: size in local sends/receive does not match in exchange!" << std::endl;
01366                           return;
01367                         }
01368                   memcpy(_localrecvrequests[i].buffer,_localsendrequests[i].buffer,_localsendrequests[i].size);
01369                 }
01370           _localsendrequests.clear();
01371           _localrecvrequests.clear();
01372 
01373 #if HAVE_MPI
01374           // handle foreign requests
01375           int sends=0;
01376           int recvs=0;
01377 
01378           // issue sends to foreign processes
01379           for (unsigned int i=0; i<_sendrequests.size(); i++)
01380                 if (_sendrequests[i].rank!=rank())
01381                   {
01382 //                      std::cout << "[" << rank() << "]" << " send " << _sendrequests[i].size << " bytes " 
01383 //                                        << "to " << _sendrequests[i].rank << " p=" << _sendrequests[i].buffer << std::endl;
01384                         MPI_Isend(_sendrequests[i].buffer, _sendrequests[i].size, MPI_BYTE, 
01385                                           _sendrequests[i].rank, _tag, _comm, &(_sendrequests[i].request));
01386                         _sendrequests[i].flag = false;
01387                         sends++;
01388                   }
01389 
01390           // issue receives from foreign processes
01391           for (unsigned int i=0; i<_recvrequests.size(); i++)
01392                 if (_recvrequests[i].rank!=rank())
01393                   {
01394 //                      std::cout << "[" << rank() << "]"  << " recv " << _recvrequests[i].size << " bytes " 
01395 //                                        << "fm " << _recvrequests[i].rank << " p=" << _recvrequests[i].buffer << std::endl;
01396                         MPI_Irecv(_recvrequests[i].buffer, _recvrequests[i].size, MPI_BYTE, 
01397                                           _recvrequests[i].rank, _tag, _comm, &(_recvrequests[i].request));
01398                         _recvrequests[i].flag = false;
01399                         recvs++;
01400                   }
01401 
01402           // poll sends
01403           while (sends>0)
01404                 {
01405                   for (unsigned int i=0; i<_sendrequests.size(); i++)
01406                         if (!_sendrequests[i].flag)
01407                           {
01408                                 MPI_Status status;
01409                                 MPI_Test( &(_sendrequests[i].request), &(_sendrequests[i].flag), &status);
01410                                 if (_sendrequests[i].flag)
01411                                   {
01412                                         sends--;
01413 //                                      std::cout << "[" << rank() << "]"  << " send to " << _sendrequests[i].rank << " OK" << std::endl;
01414                                   }
01415                           }
01416                 }
01417 
01418           // poll receives
01419           while (recvs>0)
01420                 {
01421                   for (unsigned int i=0; i<_recvrequests.size(); i++)
01422                         if (!_recvrequests[i].flag)
01423                           {
01424                                 MPI_Status status;
01425                                 MPI_Test( &(_recvrequests[i].request), &(_recvrequests[i].flag), &status);
01426                                 if (_recvrequests[i].flag)
01427                                   {
01428                                         recvs--;
01429 //                                      std::cout << "[" << rank() << "]"  << " recv fm " << _recvrequests[i].rank << " OK" << std::endl;
01430                                   }
01431 
01432                           }
01433                 }
01434 
01435           // clear request buffers
01436           _sendrequests.clear();
01437           _recvrequests.clear();
01438 #endif
01439         }
01440 
01442         double global_sum (double x) const
01443         {
01444           double res;
01445         
01446           if (_procs==1) return x;
01447 #if HAVE_MPI
01448           MPI_Allreduce(&x,&res,1,MPI_DOUBLE,MPI_SUM,_comm);
01449           return res;
01450 #endif
01451         }
01452 
01454         double global_max (double x) const
01455         {
01456           double res;
01457         
01458           if (_procs==1) return x;
01459 #if HAVE_MPI
01460           MPI_Allreduce(&x,&res,1,MPI_DOUBLE,MPI_MAX,_comm);
01461           return res;
01462 #endif
01463         }
01464 
01466         double global_min (double x) const
01467         {
01468           double res;
01469         
01470           if (_procs==1) return x;
01471 #if HAVE_MPI
01472           MPI_Allreduce(&x,&res,1,MPI_DOUBLE,MPI_MIN,_comm);
01473           return res;
01474 #endif
01475         }
01476 
01477  
01479         void print (std::ostream& s) const
01480         {
01481           s << "[" << rank() <<  "]: Torus " << procs() << " processor(s) arranged as " << dims() << std::endl;
01482           for (ProcListIterator i=sendbegin(); i!=sendend(); ++i)
01483                 {
01484                   s << "[" << rank() <<  "]: send to   " 
01485                         << "rank=" << i.rank() 
01486                         << " index=" << i.index() 
01487                         << " delta=" << i.delta() << " dist=" << i.distance() << std::endl;
01488                 }
01489           for (ProcListIterator i=recvbegin(); i!=recvend(); ++i)
01490                 {
01491                   s << "[" << rank() <<  "]: recv from " 
01492                         << "rank=" << i.rank() 
01493                         << " index=" << i.index() 
01494                         << " delta=" << i.delta() << " dist=" << i.distance() << std::endl;
01495                 }
01496         }
01497 
01498   private:
01499 
01500         void optimize_dims (int i, iTupel& size, int P, iTupel& dims, double &opt )
01501         {
01502           if (i>0) // test all subdivisions recursively
01503                 {
01504                   for (int k=1; k<=P; k++)
01505                         if (P%k==0)
01506                           {
01507                                 // P divisible by k
01508                                 dims[i] = k;
01509                                 optimize_dims(i-1,size,P/k,dims,opt);
01510                           }
01511                 }
01512           else
01513                 {
01514                   // found a possible combination
01515                   dims[0] = P;
01516                   
01517                   // check for optimality
01518                   double m = -1.0;
01519 
01520                   for (int k=0; k<d; k++)
01521                         {
01522                       double mm=((double)size[k])/((double)dims[k]);
01523                       if (fmod((double)size[k],(double)dims[k])>0.0001) mm*=3;
01524                       if ( mm > m ) m = mm;
01525                         }
01526                   //              if (_rank==0) std::cout << "optimize_dims: " << size << " | " << dims << " norm=" << m << std::endl;  
01527                   if (m<opt) 
01528                         {
01529                           opt = m;
01530                           _dims = dims;
01531                         }
01532                 }
01533         }
01534 
01535         void proclists ()
01536         {
01537           // compile the full neighbor list
01538           CommPartner cp;
01539           iTupel delta;
01540 
01541           delta = -1;
01542           bool ready = false;
01543           iTupel me, nb;
01544           me = rank_to_coord(_rank);
01545           int index = 0;
01546           int last = neighbors()-1;
01547           while (!ready)
01548                 {
01549                   // find neighbors coordinates
01550                   for (int i=0; i<d; i++)
01551                         nb[i] = ( me[i]+_dims[i]+delta[i] ) % _dims[i];
01552 
01553                   // find neighbors rank
01554                   int nbrank = coord_to_rank(nb);
01555                   
01556                   // check if delta is not zero
01557                   for (int i=0; i<d; i++)
01558                         if (delta[i]!=0) 
01559                           {
01560                                 cp.rank = nbrank;
01561                                 cp.delta = delta;
01562                                 cp.index = index;
01563                                 _recvlist.push_back(cp);
01564                                 cp.index = last-index;
01565                                 _sendlist.push_front(cp);
01566                                 index++;
01567                                 break;
01568                           }
01569 
01570                   // next neighbor
01571                   ready = true;
01572                   for (int i=0; i<d; i++)
01573                         if (delta[i]<1)
01574                           {
01575                                 (delta[i])++;
01576                                 ready=false;
01577                                 break;
01578                           }
01579                         else
01580                           {
01581                                 delta[i] = -1;
01582                           }
01583                 }
01584 
01585         }
01586 
01587 #if HAVE_MPI
01588         MPI_Comm _comm;
01589 #endif
01590         int _rank;
01591         int _procs;
01592         iTupel _dims;
01593         iTupel _increment;
01594         int _tag;
01595         std::deque<CommPartner> _sendlist;
01596         std::deque<CommPartner> _recvlist;
01597  
01598         mutable std::vector<CommTask> _sendrequests;
01599         mutable std::vector<CommTask> _recvrequests;
01600         mutable std::vector<CommTask> _localsendrequests;
01601         mutable std::vector<CommTask> _localrecvrequests;
01602   };
01603 
01605   template <int d>
01606   inline std::ostream& operator<< (std::ostream& s, const Torus<d> & t)
01607   {
01608         t.print(s);
01609         return s;
01610   }
01611 
01612 
01615   template<int d, typename ct>
01616   class MultiYGrid {
01617   public:
01618         // some data types
01619         struct Intersection {
01620           SubYGrid<d,ct> grid; // the intersection as a subgrid of local grid
01621           int rank;            // rank of process where other grid is stored
01622           int distance;        // manhattan distance to other grid
01623         };
01624         
01625         struct YGridLevel {        // This stores all the information on one grid level 
01626           // cell (codim 0) data
01627           YGrid<d,ct> cell_global;         // the whole cell grid on that level
01628           SubYGrid<d,ct> cell_overlap;     // we have no ghost cells, so our part is overlap completely
01629           SubYGrid<d,ct> cell_interior;    // interior cells are a subgrid of all cells
01630           
01631           std::deque<Intersection> send_cell_overlap_overlap;  // each intersection is a subgrid of overlap
01632           std::deque<Intersection> recv_cell_overlap_overlap;  // each intersection is a subgrid of overlap
01633           
01634       std::deque<Intersection> send_cell_interior_overlap; // each intersection is a subgrid of overlap
01635       std::deque<Intersection> recv_cell_overlap_interior; // each intersection is a subgrid of overlap
01636           
01637           // vertex (codim dim) data
01638       YGrid<d,ct> vertex_global;           // the whole vertex grid on that level
01639       SubYGrid<d,ct> vertex_overlapfront;  // all our vertices are overlap and front
01640       SubYGrid<d,ct> vertex_overlap;       // subgrid containing only overlap
01641       SubYGrid<d,ct> vertex_interiorborder;// subgrid containing only interior and border
01642       SubYGrid<d,ct> vertex_interior;      // subgrid containing only interior
01643 
01644       std::deque<Intersection> send_vertex_overlapfront_overlapfront; // each intersection is a subgrid of overlapfront
01645       std::deque<Intersection> recv_vertex_overlapfront_overlapfront; // each intersection is a subgrid of overlapfront
01646 
01647       std::deque<Intersection> send_vertex_overlap_overlapfront; // each intersection is a subgrid of overlapfront
01648       std::deque<Intersection> recv_vertex_overlapfront_overlap; // each intersection is a subgrid of overlapfront
01649 
01650           std::deque<Intersection> send_vertex_interiorborder_interiorborder; // each intersection is a subgrid of overlapfront
01651           std::deque<Intersection> recv_vertex_interiorborder_interiorborder; // each intersection is a subgrid of overlapfront
01652         
01653           std::deque<Intersection> send_vertex_interiorborder_overlapfront; // each intersection is a subgrid of overlapfront
01654           std::deque<Intersection> recv_vertex_overlapfront_interiorborder; // each intersection is a subgrid of overlapfront
01655         
01656           // general
01657       MultiYGrid<d,ct>* mg;  // each grid level knows its multigrid
01658       int overlap;           // in mesh cells on this level
01659     };
01660 
01662         typedef FieldVector<int, d> iTupel;
01663         typedef FieldVector<ct, d> fTupel;
01664         typedef FieldVector<bool, d> bTupel;
01665 
01666         // communication tag used by multigrid
01667         enum { tag = 17 };
01668 
01670 #if HAVE_MPI
01671         MultiYGrid (MPI_Comm comm, fTupel L, iTupel s, bTupel periodic, int overlap) 
01672           : _torus(comm,tag,s) // torus gets s to compute procs/direction
01673         {
01674           // store parameters
01675           _LL = L;
01676           _s = s;
01677           _periodic = periodic;
01678           _overlap = overlap;
01679 
01680           // coarse cell interior  grid obtained through partitioning of global grid
01681           iTupel o_interior;
01682           iTupel s_interior;
01683           iTupel o = iTupel(0);
01684           double imbal = _torus.partition(_torus.rank(),o,s,o_interior,s_interior);
01685           imbal = _torus.global_max(imbal);
01686 
01687           // add level
01688           _maxlevel = 0;
01689           _levels[_maxlevel] = makelevel(L,s,periodic,o_interior,s_interior,overlap);
01690 
01691           // output
01692 //        if (_torus.rank()==0) std::cout << "MultiYGrid<" << d // changed dinfo to cout
01693 //                                                                        << ">: coarse grid with size " << s 
01694 //                                                                        << " imbalance=" << (imbal-1)*100 << "%" << std::endl;
01695           //      print(std::cout);
01696         }
01697 #else
01698         MultiYGrid (fTupel L, iTupel s, bTupel periodic, int overlap) 
01699           : _torus(tag,s) // torus gets s to compute procs/direction
01700         {
01701           // store parameters
01702           _LL = L;
01703           _s = s;
01704           _periodic = periodic;
01705           _overlap = overlap;
01706 
01707           // coarse cell interior  grid obtained through partitioning of global grid
01708           iTupel o = iTupel(0);
01709           iTupel o_interior(o);
01710           iTupel s_interior(s);
01711 
01712           // add level
01713           _maxlevel = 0;
01714           _levels[_maxlevel] = makelevel(L,s,periodic,o_interior,s_interior,overlap);
01715           // output
01716 //        if (_torus.rank()==0) std::cout << "MultiYGrid<" << d // changed dinfo to cout
01717 //                                                                        << ">: coarse grid with size " << s 
01718 //                                                                        << " imbalance=" << (imbal-1)*100 << "%" << std::endl;
01719           //      print(std::cout);
01720         }
01721 #endif
01722 
01724         void refine (bool keep_overlap)
01725         {
01726           // access to coarser grid level
01727           YGridLevel& cg = _levels[maxlevel()];
01728 
01729           // compute size of new global grid
01730           iTupel s;
01731           for (int i=0; i<d; i++) s[i] = 2*cg.cell_global.size(i);
01732 
01733           // compute overlap
01734           int overlap;
01735           if (keep_overlap) overlap = 2*cg.overlap; else overlap = cg.overlap;
01736 
01737           // output
01738 //        if (_torus.rank()==0) std::cout << "MultiYGrid<" // changed dinfo to cout
01739 //                                                                        << d << ">: refined to size " 
01740 //                                                                        << s << std::endl;
01741 
01742           // the cell interior grid obtained from coarse cell interior grid
01743           iTupel o_interior;
01744           iTupel s_interior;
01745           for (int i=0; i<d; i++) o_interior[i] = 2*cg.cell_interior.origin(i);
01746           for (int i=0; i<d; i++) s_interior[i] = 2*cg.cell_interior.size(i);
01747           
01748           // add level
01749           _maxlevel++;
01750           _levels[_maxlevel] = makelevel(_LL,s,_periodic,o_interior,s_interior,overlap);
01751         }
01752 
01754         const Torus<d>& torus () const
01755         {
01756           return _torus;
01757         }
01758 
01760         int maxlevel () const
01761         {
01762           return _maxlevel;
01763         }
01764 
01766         bool periodic (int i) const
01767         {
01768           return _periodic[i];
01769         }
01770 
01772         class YGridLevelIterator {
01773         private:
01774           int l;
01775           const YGridLevel* i;
01776         public:
01778           YGridLevelIterator ()
01779           {
01780           }
01781 
01783           YGridLevelIterator (const YGridLevel* start, int level)
01784           {
01785                 i=start; l=level;
01786           }
01787 
01789           YGridLevelIterator (const YGridLevelIterator & it)
01790         : l(it.l), i(it.i)
01791       {}
01792 
01794           int level () const
01795           {
01796                 return l;
01797           }
01798 
01800           int overlap () const
01801           {
01802                 return i->overlap;
01803           }
01804 
01806           const MultiYGrid<d,ct>* mg () const
01807           {
01808                 return i->mg;
01809           }
01810 
01812           bool operator== (const YGridLevelIterator& iter) const
01813           {
01814                 return i == iter.i;
01815           }
01816 
01818           bool operator!= (const YGridLevelIterator& iter) const
01819           {
01820                 return i != iter.i;
01821           }
01822         
01824           YGridLevelIterator& operator++ ()
01825           {
01826                 ++i; // assumes built-in array
01827                 ++l;
01828                 return *this;
01829           }
01830 
01832           YGridLevelIterator& operator-- ()
01833           {
01834                 --i;
01835                 --l;
01836                 return *this;
01837           }
01838 
01840           YGridLevelIterator finer () const
01841           {
01842                 return YGridLevelIterator(i+1,l+1);
01843           }
01844 
01846           YGridLevelIterator coarser () const
01847           {
01848                 return YGridLevelIterator(i-1,l-1);
01849           }
01850 
01852           const YGrid<d,ct>& cell_global () const
01853           {
01854                 return i->cell_global;
01855           }
01856 
01858           const SubYGrid<d,ct>& cell_overlap () const
01859           {
01860                 return i->cell_overlap;
01861           }
01862 
01864           const SubYGrid<d,ct>& cell_interior () const
01865           {
01866                 return i->cell_interior;
01867           }
01868 
01869 
01871           const std::deque<Intersection>& send_cell_overlap_overlap () const
01872           {
01873                 return i->send_cell_overlap_overlap;
01874           }
01875           const std::deque<Intersection>& recv_cell_overlap_overlap () const
01876           {
01877                 return i->recv_cell_overlap_overlap;
01878           }
01879           const std::deque<Intersection>& send_cell_interior_overlap () const
01880           {
01881                 return i->send_cell_interior_overlap;
01882           }
01883           const std::deque<Intersection>& recv_cell_overlap_interior () const 
01884           {
01885                 return i->recv_cell_overlap_interior;
01886           }
01887 
01888 
01890           const YGrid<d,ct>& vertex_global () const
01891           {
01892                 return i->vertex_global;
01893           }
01895           const SubYGrid<d,ct>& vertex_overlapfront () const
01896           {
01897                 return i->vertex_overlapfront;
01898           }
01900           const SubYGrid<d,ct>& vertex_overlap () const
01901           {
01902                 return i->vertex_overlap;
01903           }
01905           const SubYGrid<d,ct>& vertex_interiorborder () const
01906           {
01907                 return i->vertex_interiorborder;
01908           }
01910           const SubYGrid<d,ct>& vertex_interior () const
01911           {
01912                 return i->vertex_interior;
01913           }
01914 
01916           const std::deque<Intersection>& send_vertex_overlapfront_overlapfront () const
01917           {
01918                 return i->send_vertex_overlapfront_overlapfront;
01919           }
01920           const std::deque<Intersection>& recv_vertex_overlapfront_overlapfront () const
01921           {
01922                 return i->recv_vertex_overlapfront_overlapfront;
01923           }
01924           const std::deque<Intersection>& send_vertex_overlap_overlapfront () const
01925           {
01926                 return i->send_vertex_overlap_overlapfront;
01927           }
01928           const std::deque<Intersection>& recv_vertex_overlapfront_overlap () const
01929           {
01930                 return i->recv_vertex_overlapfront_overlap;
01931           }
01932           const std::deque<Intersection>& send_vertex_interiorborder_interiorborder () const
01933           {
01934                 return i->send_vertex_interiorborder_interiorborder;
01935           }
01936           const std::deque<Intersection>& recv_vertex_interiorborder_interiorborder () const
01937           {
01938                 return i->recv_vertex_interiorborder_interiorborder;
01939           }
01940           const std::deque<Intersection>& send_vertex_interiorborder_overlapfront () const
01941           {
01942                 return i->send_vertex_interiorborder_overlapfront;
01943           }
01944           const std::deque<Intersection>& recv_vertex_overlapfront_interiorborder () const
01945           {
01946                 return i->recv_vertex_overlapfront_interiorborder;
01947           }
01948         };
01949 
01951         YGridLevelIterator begin () const
01952         {
01953           return YGridLevelIterator(_levels,0);
01954         }
01955 
01957         YGridLevelIterator begin (int i) const
01958         {
01959           if (i<0 || i>maxlevel())
01960                 DUNE_THROW(GridError, "level not existing");
01961           return YGridLevelIterator(_levels+i,i);
01962         }
01963 
01965         YGridLevelIterator end () const
01966         {
01967           return YGridLevelIterator(_levels+(_maxlevel+1),_maxlevel+1);
01968         }
01969 
01971         YGridLevelIterator rbegin () const
01972         {
01973           return YGridLevelIterator(_levels+_maxlevel,_maxlevel);
01974         }
01975 
01977         YGridLevelIterator rend () const
01978         {
01979           return YGridLevelIterator(_levels-1,-1);
01980         }
01981 
01983         inline void print (std::ostream& s) const
01984         {
01985           int rank = torus().rank();
01986 
01987           s << "[" << rank << "]:" << " MultiYGrid maxlevel=" << maxlevel() << std::endl;
01988 
01989           for (YGridLevelIterator g=begin(); g!=end(); ++g)
01990                 {
01991                   s << "[" << rank << "]:   " << std::endl;
01992                   s << "[" << rank << "]:   " << "==========================================" << std::endl;
01993                   s << "[" << rank << "]:   " << "level=" << g.level() << std::endl;
01994                   s << "[" << rank << "]:   " << "cell_global=" << g.cell_global() << std::endl;
01995                   s << "[" << rank << "]:   " << "cell_overlap=" << g.cell_overlap() << std::endl;
01996                   s << "[" << rank << "]:   " << "cell_interior=" << g.cell_interior() << std::endl;          
01997                   for (typename std::deque<Intersection>::const_iterator i=g.send_cell_overlap_overlap().begin();
01998                            i!=g.send_cell_overlap_overlap().end(); ++i)
01999                         {
02000                           s << "[" << rank << "]:    " << " s_c_o_o "
02001                                 << i->rank << " " << i->grid << std::endl;
02002                         }
02003                   for (typename std::deque<Intersection>::const_iterator i=g.recv_cell_overlap_overlap().begin();
02004                            i!=g.recv_cell_overlap_overlap().end(); ++i)
02005                         {
02006                           s << "[" << rank << "]:    " << " r_c_o_o "
02007                                 << i->rank << " " << i->grid << std::endl;
02008                         }
02009                   for (typename std::deque<Intersection>::const_iterator i=g.send_cell_interior_overlap().begin();
02010                            i!=g.send_cell_interior_overlap().end(); ++i)
02011                         {
02012                           s << "[" << rank << "]:    " << " s_c_i_o "
02013                                 << i->rank << " " << i->grid << std::endl;
02014                         }
02015                   for (typename std::deque<Intersection>::const_iterator i=g.recv_cell_overlap_interior().begin();
02016                            i!=g.recv_cell_overlap_interior().end(); ++i)
02017                         {
02018                           s << "[" << rank << "]:    " << " r_c_o_i "
02019                                 << i->rank << " " << i->grid << std::endl;
02020                         }
02021 
02022                   s << "[" << rank << "]:   " << "-----------------------------------------------"  << std::endl;
02023                   s << "[" << rank << "]:   " << "vertex_global="         << g.vertex_global() << std::endl;
02024                   s << "[" << rank << "]:   " << "vertex_overlapfront="   << g.vertex_overlapfront() << std::endl;
02025                   s << "[" << rank << "]:   " << "vertex_overlap="        << g.vertex_overlap() << std::endl;
02026                   s << "[" << rank << "]:   " << "vertex_interiorborder=" << g.vertex_interiorborder() << std::endl;
02027                   s << "[" << rank << "]:   " << "vertex_interior="       << g.vertex_interior() << std::endl;
02028                   for (typename std::deque<Intersection>::const_iterator i=g.send_vertex_overlapfront_overlapfront().begin();
02029                            i!=g.send_vertex_overlapfront_overlapfront().end(); ++i)
02030                         {
02031                           s << "[" << rank << "]:    " << " s_v_of_of "
02032                                 << i->rank << " " << i->grid << std::endl;
02033                         }
02034                   for (typename std::deque<Intersection>::const_iterator i=g.recv_vertex_overlapfront_overlapfront().begin();
02035                            i!=g.recv_vertex_overlapfront_overlapfront().end(); ++i)
02036                         {
02037                           s << "[" << rank << "]:    " << " r_v_of_of "
02038                                 << i->rank << " " << i->grid << std::endl;
02039                         }
02040                   for (typename std::deque<Intersection>::const_iterator i=g.send_vertex_overlap_overlapfront().begin();
02041                            i!=g.send_vertex_overlap_overlapfront().end(); ++i)
02042                         {
02043                           s << "[" << rank << "]:    " << " s_v_o_of "
02044                                 << i->rank << " " << i->grid << std::endl;
02045                         }
02046                   for (typename std::deque<Intersection>::const_iterator i=g.recv_vertex_overlapfront_overlap().begin();
02047                            i!=g.recv_vertex_overlapfront_overlap().end(); ++i)
02048                         {
02049                           s << "[" << rank << "]:    " << " r_v_of_o "
02050                                 << i->rank << " " << i->grid << std::endl;
02051                         }
02052                   for (typename std::deque<Intersection>::const_iterator i=g.send_vertex_interiorborder_interiorborder().begin();
02053                            i!=g.send_vertex_interiorborder_interiorborder().end(); ++i)
02054                         {
02055                           s << "[" << rank << "]:    " << " s_v_ib_ib "
02056                                 << i->rank << " " << i->grid << std::endl;
02057                         }
02058                   for (typename std::deque<Intersection>::const_iterator i=g.recv_vertex_interiorborder_interiorborder().begin();
02059                            i!=g.recv_vertex_interiorborder_interiorborder().end(); ++i)
02060                         {
02061                           s << "[" << rank << "]:    " << " r_v_ib_ib "
02062                                 << i->rank << " " << i->grid << std::endl;
02063                         }
02064                   for (typename std::deque<Intersection>::const_iterator i=g.send_vertex_interiorborder_overlapfront().begin();
02065                            i!=g.send_vertex_interiorborder_overlapfront().end(); ++i)
02066                         {
02067                           s << "[" << rank << "]:    " << " s_v_ib_of "
02068                                 << i->rank << " " << i->grid << std::endl;
02069                         }
02070                   for (typename std::deque<Intersection>::const_iterator i=g.recv_vertex_overlapfront_interiorborder().begin();
02071                            i!=g.recv_vertex_overlapfront_interiorborder().end(); ++i)
02072                         {
02073                           s << "[" << rank << "]:    " << " s_v_of_ib "
02074                                 << i->rank << " " << i->grid << std::endl;
02075                         }
02076                 }
02077 
02078           s << std::endl;
02079         }
02080 
02081   private:
02082     // make a new YGridLevel structure. For that we need
02083     // L           size of the whole domain in each direction
02084     // s           number of cells in each direction
02085     // periodic    boolean indication periodicity in each direction
02086     // o_interior  origin of interior (non-overlapping) cell decomposition
02087     // s_interior  size of interior cell decomposition
02088     // overlap     to be used on this grid level
02089         YGridLevel makelevel (fTupel L, iTupel s, bTupel periodic, iTupel o_interior, iTupel s_interior, int overlap)
02090     {
02091           // first, lets allocate a new structure
02092           YGridLevel g;
02093           g.overlap = overlap;
02094           g.mg = this;
02095 
02096           // the global cell grid
02097           iTupel o = iTupel(0); // logical origin is always 0, that is not a restriction
02098           fTupel h;
02099           fTupel r;
02100           for (int i=0; i<d; i++) h[i] = L[i]/s[i]; // the mesh size in each direction
02101           for (int i=0; i<d; i++) r[i] = 0.5*h[i];  // the shift for cell centers
02102           g.cell_global = YGrid<d,ct>(o,s,h,r);     // this is the global cell grid
02103 
02104           // extend the cell interior grid by overlap considering periodicity
02105           iTupel o_overlap;
02106           iTupel s_overlap;
02107           for (int i=0; i<d; i++)
02108                 {
02109                   if (periodic[i])
02110                         {
02111                           // easy case, extend by 2 overlaps in total
02112                           o_overlap[i] = o_interior[i]-overlap;  // Note: origin might be negative now 
02113                           s_overlap[i] = s_interior[i]+2*overlap;// Note: might be larger than global size
02114                         }
02115                   else
02116                         {
02117                           // nonperiodic case, intersect with global size
02118                           int min = std::max(0,o_interior[i]-overlap);
02119                           int max = std::min(s[i]-1,o_interior[i]+s_interior[i]-1+overlap);
02120                           o_overlap[i] = min;
02121                           s_overlap[i] = max-min+1;
02122                         }
02123                 }
02124           g.cell_overlap = SubYGrid<d,ct>(YGrid<d,ct>(o_overlap,s_overlap,h,r));
02125 
02126           // now make the interior grid a subgrid of the overlapping grid
02127           iTupel offset;
02128           for (int i=0; i<d; i++) offset[i] = o_interior[i]-o_overlap[i];
02129           g.cell_interior = SubYGrid<d,ct>(o_interior,s_interior,offset,s_overlap,h,r);
02130 
02131           // compute cell intersections 
02132           intersections(g.cell_overlap,g.cell_overlap,g.cell_global.size(),g.send_cell_overlap_overlap,g.recv_cell_overlap_overlap);
02133           intersections(g.cell_interior,g.cell_overlap,g.cell_global.size(),g.send_cell_interior_overlap,g.recv_cell_overlap_interior);
02134 
02135           // now we can do the vertex grids. They are derived completely from the cell grids
02136           iTupel o_vertex_global, s_vertex_global;
02137           for (int i=0; i<d; i++) r[i] = 0.0;  // the shift for vertices is zero, and the mesh size is same as for cells
02138 
02139           // first let's make the global grid
02140           for (int i=0; i<d; i++) o_vertex_global[i] = g.cell_global.origin(i);
02141           for (int i=0; i<d; i++) s_vertex_global[i] = g.cell_global.size(i)+1; // one more vertices than cells ...
02142           g.vertex_global = YGrid<d,ct>(o_vertex_global,s_vertex_global,h,r);
02143 
02144           // now the local grid stored in this processor. All other grids are subgrids of this
02145           iTupel o_vertex_overlapfront;
02146           iTupel s_vertex_overlapfront;
02147           for (int i=0; i<d; i++) o_vertex_overlapfront[i] = g.cell_overlap.origin(i);
02148           for (int i=0; i<d; i++) s_vertex_overlapfront[i] = g.cell_overlap.size(i)+1; // one more vertices than cells ...
02149           g.vertex_overlapfront = SubYGrid<d,ct>(YGrid<d,ct>(o_vertex_overlapfront,s_vertex_overlapfront,h,r));
02150 
02151           // now overlap only (i.e. without front), is subgrid of overlapfront
02152           iTupel o_vertex_overlap;
02153           iTupel s_vertex_overlap;
02154           for (int i=0; i<d; i++) 
02155                 {
02156                   o_vertex_overlap[i] = g.cell_overlap.origin(i);
02157                   s_vertex_overlap[i] = g.cell_overlap.size(i)+1;
02158 
02159                   if (!periodic[i] && g.cell_overlap.origin(i)>g.cell_global.origin(i))
02160                   {
02161                         // not at the lower boundary
02162                         o_vertex_overlap[i] += 1;
02163                         s_vertex_overlap[i] -= 1;
02164                   }
02165 
02166                   if (!periodic[i] && g.cell_overlap.origin(i)+g.cell_overlap.size(i)<g.cell_global.origin(i)+g.cell_global.size(i))
02167                   {
02168                         // not at the upper boundary
02169                         s_vertex_overlap[i] -= 1;
02170                   }
02171 
02172 
02173                   offset[i] = o_vertex_overlap[i]-o_vertex_overlapfront[i];
02174                 }
02175           g.vertex_overlap = SubYGrid<d,ct>(o_vertex_overlap,s_vertex_overlap,offset,s_vertex_overlapfront,h,r);
02176 
02177           // now interior with border
02178           iTupel o_vertex_interiorborder;
02179           iTupel s_vertex_interiorborder;
02180           for (int i=0; i<d; i++) o_vertex_interiorborder[i] = g.cell_interior.origin(i);
02181           for (int i=0; i<d; i++) s_vertex_interiorborder[i] = g.cell_interior.size(i)+1;
02182           for (int i=0; i<d; i++) offset[i] = o_vertex_interiorborder[i]-o_vertex_overlapfront[i];
02183           g.vertex_interiorborder = SubYGrid<d,ct>(o_vertex_interiorborder,s_vertex_interiorborder,offset,s_vertex_overlapfront,h,r);
02184 
02185           // now only interior
02186           iTupel o_vertex_interior;
02187           iTupel s_vertex_interior;
02188           for (int i=0; i<d; i++) 
02189                 {
02190                   o_vertex_interior[i] = g.cell_interior.origin(i);
02191                   s_vertex_interior[i] = g.cell_interior.size(i)+1;
02192 
02193                   if (!periodic[i] && g.cell_interior.origin(i)>g.cell_global.origin(i))
02194                   {
02195                         // not at the lower boundary
02196                         o_vertex_interior[i] += 1;
02197                         s_vertex_interior[i] -= 1;
02198                   }
02199 
02200                   if (!periodic[i] && g.cell_interior.origin(i)+g.cell_interior.size(i)<g.cell_global.origin(i)+g.cell_global.size(i))
02201                   {
02202                         // not at the upper boundary
02203                         s_vertex_interior[i] -= 1;
02204                   }
02205 
02206                   offset[i] = o_vertex_interior[i]-o_vertex_overlapfront[i];
02207                 }
02208           g.vertex_interior = SubYGrid<d,ct>(o_vertex_interior,s_vertex_interior,offset,s_vertex_overlapfront,h,r);
02209 
02210           // compute vertex intersections
02211           intersections(g.vertex_overlapfront,g.vertex_overlapfront,g.cell_global.size(),
02212                                         g.send_vertex_overlapfront_overlapfront,g.recv_vertex_overlapfront_overlapfront);
02213           intersections(g.vertex_overlap,g.vertex_overlapfront,g.cell_global.size(),
02214                                         g.send_vertex_overlap_overlapfront,g.recv_vertex_overlapfront_overlap);
02215           intersections(g.vertex_interiorborder,g.vertex_interiorborder,g.cell_global.size(),
02216                                         g.send_vertex_interiorborder_interiorborder,g.recv_vertex_interiorborder_interiorborder);
02217           intersections(g.vertex_interiorborder,g.vertex_overlapfront,g.cell_global.size(),
02218                                         g.send_vertex_interiorborder_overlapfront,g.recv_vertex_overlapfront_interiorborder);
02219 
02220           // return the whole thing
02221           return g;
02222         }
02223 
02224 
02225         // construct list of intersections with neighboring processors:
02226     //   recvgrid: the grid stored in this processor
02227     //   sendgrid:  the subgrid to be sent to neighboring processors
02228     //   size: needed to shift local grid in periodic case
02229     //   returns two lists: Intersections to be sent and Intersections to be received
02230         // Note: sendgrid/recvgrid may be SubYGrids. Since intersection method is virtual it should work properly
02231         void intersections (const SubYGrid<d,ct>& sendgrid, const SubYGrid<d,ct>& recvgrid, const iTupel& size, 
02232                                                 std::deque<Intersection>& sendlist, std::deque<Intersection>& recvlist)
02233         {
02234           // the exchange buffers
02235           std::vector<YGrid<d,ct> > send_recvgrid(_torus.neighbors()); 
02236           std::vector<YGrid<d,ct> > recv_recvgrid(_torus.neighbors());
02237           std::vector<YGrid<d,ct> > send_sendgrid(_torus.neighbors()); 
02238           std::vector<YGrid<d,ct> > recv_sendgrid(_torus.neighbors());
02239 
02240           // fill send buffers; iterate over all neighboring processes
02241           // non-periodic case is handled automatically because intersection will be zero
02242           for (typename Torus<d>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
02243                 {
02244                   // determine if we communicate with this neighbor (and what)
02245                   bool skip = false;
02246                   iTupel coord = _torus.coord(); // my coordinates
02247                   iTupel delta = i.delta();      // delta to neighbor
02248                   iTupel nb = coord;             // the neighbor
02249                   for (int k=0; k<d; k++) nb[k] += delta[k];
02250                   iTupel v = iTupel(0);                  // grid movement
02251 
02252                   for (int k=0; k<d; k++)
02253                         {
02254                           if (nb[k]<0) 
02255                                 {
02256                                   if (_periodic[k])
02257                                         v[k] += size[k];
02258                                   else
02259                                         skip = true;
02260                                 }
02261                           if (nb[k]>=_torus.dims(k))
02262                                 {
02263                                   if (_periodic[k])
02264                                         v[k] -= size[k];
02265                                   else
02266                                         skip = true;
02267                                 }
02268                           // neither might be true, then v=0
02269                         }
02270 
02271                   // store moved grids in send buffers
02272                   if (!skip)
02273                         {
02274                           send_sendgrid[i.index()] = sendgrid.move(v);
02275                           send_recvgrid[i.index()] = recvgrid.move(v);
02276                         }
02277                   else
02278                         {
02279                           send_sendgrid[i.index()] = YGrid<d,ct>(iTupel(0),iTupel(0),fTupel(0.0),fTupel(0.0));
02280                           send_recvgrid[i.index()] = YGrid<d,ct>(iTupel(0),iTupel(0),fTupel(0.0),fTupel(0.0));
02281                         }
02282                 }
02283 
02284           // issue send requests for sendgrid being sent to all neighbors
02285           for (typename Torus<d>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
02286                 _torus.send(i.rank(), &send_sendgrid[i.index()], sizeof(YGrid<d,ct>));
02287 
02288           // issue recv requests for sendgrids of neighbors
02289           for (typename Torus<d>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
02290                 _torus.recv(i.rank(), &recv_sendgrid[i.index()], sizeof(YGrid<d,ct>));
02291           
02292           // exchange the sendgrids
02293           _torus.exchange();
02294 
02295           // issue send requests for recvgrid being sent to all neighbors
02296           for (typename Torus<d>::ProcListIterator i=_torus.sendbegin(); i!=_torus.sendend(); ++i)
02297                 _torus.send(i.rank(), &send_recvgrid[i.index()], sizeof(YGrid<d,ct>));
02298 
02299           // issue recv requests for recvgrid of neighbors
02300           for (typename Torus<d>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
02301                 _torus.recv(i.rank(), &recv_recvgrid[i.index()], sizeof(YGrid<d,ct>));
02302           
02303           // exchange the recvgrid
02304           _torus.exchange();
02305 
02306           // process receive buffers and compute intersections
02307           for (typename Torus<d>::ProcListIterator i=_torus.recvbegin(); i!=_torus.recvend(); ++i)
02308                 {
02309                   // what must be sent to this neighbor
02310                   Intersection send_intersection;
02311                   send_intersection.grid = sendgrid.intersection(recv_recvgrid[i.index()]);
02312 //                std::cout << "[" << _torus.rank() << "]:   " << "sendgrid=" << sendgrid << std::endl;
02313 //                std::cout << "[" << _torus.rank() << "]:   " << "recved recvgrid=" << recv_recvgrid[i.index()] << std::endl;
02314 //                std::cout << "[" << _torus.rank() << "]:   " << "intersection=" << send_intersection.grid << std::endl;
02315                   send_intersection.rank = i.rank();
02316                   send_intersection.distance = i.distance();
02317                   if (!send_intersection.grid.empty()) sendlist.push_front(send_intersection);
02318 
02319                   Intersection recv_intersection;
02320                   recv_intersection.grid = recvgrid.intersection(recv_sendgrid[i.index()]);
02321 //                std::cout << "[" << _torus.rank() << "]:   " << "recvgrid=" << recvgrid << std::endl;
02322 //                std::cout << "[" << _torus.rank() << "]:   " << "recved sendgrid=" << recv_sendgrid[i.index()] << std::endl;
02323 //                std::cout << "[" << _torus.rank() << "]:   " << "intersection=" << recv_intersection.grid << std::endl;
02324                   recv_intersection.rank = i.rank();
02325                   recv_intersection.distance = i.distance();
02326                   if(!recv_intersection.grid.empty()) recvlist.push_back(recv_intersection);
02327                 }
02328         }
02329 
02330         // private data of multigrid
02331         fTupel _LL;
02332         iTupel _s;
02333         bTupel _periodic;
02334         int _maxlevel;
02335         YGridLevel _levels[32];
02336         int _overlap;
02337         Torus<d> _torus;
02338   };
02339 
02341   template <int d, class ct>
02342   inline std::ostream& operator<< (std::ostream& s, MultiYGrid<d,ct>& mg)
02343   {
02344         mg.print(s);
02345         s << std::endl;
02346         return s;
02347   }
02348 
02349 } // namespace Dune
02350 
02351 #endif

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].