agcommunicator.hh

00001 #ifndef DUNE_AGCOMMUNICATOR_HH
00002 #define DUNE_AGCOMMUNICATOR_HH
00003 
00004 #include <cstring>
00005 
00006 // use this define to control if Albert should use the found MPI
00007 
00008 // do not use at the moment
00009 #undef ALBERTA_USES_MPI  
00010 
00011 #if defined(HAVE_MPI) && defined(ALBERTA_USES_MPI)
00012 #include <mpi.h>
00013 #endif
00014 
00015 #include <dune/common/dlist.hh>
00016 
00017 #define _ANSI_HEADER
00018 
00019 #if defined(HAVE_MPI) && defined(ALBERTA_USES_MPI)
00020 
00021 #include <dune/grid/bsgrid/systemincludes.hh>
00022 namespace BernhardSchuppGrid {
00023 #include <dune/grid/bsgrid/serialize.h>
00024 //#include <dune/grid/bsgrid/mpAccess_MPI.h>
00025 #include <dune/grid/bsgrid/mpAccess_MPI.cc>
00026 #include "loadbalancer.cc"
00027 }
00028 typedef BernhardSchuppGrid :: ObjectStream ObjectStreamType;
00029 typedef BernhardSchuppGrid :: ObjectStream AlbertaObjectStream;
00030 #endif
00031 
00032 #include <map>
00033 
00034 namespace Dune {
00035 
00036 
00037   enum ObjectId { BEGINELEMENT = -665 , ENDOFSTREAM = -666 , REFINEEL = 1 , STOPHERE = 0 };
00038 
00039   static const int COMMUNICATOR_COMM_TAG = 457;
00040 
00047 template <class DofManagerType>  
00048 class CommunicatorInterface 
00049 {
00050 public: 
00051   virtual bool firstMark() = 0;
00052   virtual bool secondMark() = 0;
00053   virtual bool thirdMark() = 0;
00054   
00055   virtual bool markFirstLevel() = 0;
00056   virtual bool markNextLevel () = 0;
00057   
00058   virtual bool xtractData (DofManagerType & dm) = 0;
00059   
00060   virtual bool repartition (DofManagerType & dm) = 0;
00061   virtual bool communicate (DofManagerType & dm) = 0;
00062   virtual bool consistencyGhosts () = 0;
00063 };
00064 
00065 #if defined(HAVE_MPI) && defined(ALBERTA_USES_MPI)
00066 typedef BernhardSchuppGrid :: ObjectStream ObjectStreamType;
00067 typedef BernhardSchuppGrid :: LoadBalancer LoadBalancer;
00068 typedef BernhardSchuppGrid :: MpAccessMPI MpAccessMPI;
00069 template <class T>
00070 struct MPIType
00071 {
00072   enum { MpiType = MPI_BYTE };
00073 };
00074   
00075 template <> struct MPIType<double> { enum { MpiType = MPI_DOUBLE }; };
00076 template <> struct MPIType<int>    { enum { MpiType = MPI_INT    }; };
00077 
00078 static int cycle_ = 0;
00079  
00080 template <class GridType, class DofManagerType>
00081 class AlbertGridCommunicator : public CommunicatorInterface<DofManagerType>
00082 {
00083 public: 
00084 
00086   AlbertGridCommunicator(MPI_Comm mpiComm, GridType &grid, int mySize) 
00087     : grid_ (grid) , myRank_ (grid.myRank() ) , _ldbOver(1.2) , _ldbUnder(0.0) 
00088     , mpAccess_ ( mpiComm )
00089     , elmap_ (mySize)
00090     , elmap2_(mySize)
00091     , osv_ (0), interiorEls_ (0) , ghostEls_(0)
00092     , pSize_(-1)
00093   {
00094     int test = MPI_Comm_dup (mpiComm, &_mpiComm) ;
00095     assert (test == MPI_SUCCESS) ;
00096 
00097     int i ;
00098     test = MPI_Comm_size (_mpiComm, & i) ;
00099     assert (test == MPI_SUCCESS) ;
00100     const_cast< int & > (pSize_) = i ;
00101     createLinkage();
00102     //sprintf(name,"data_%d/grid",myRank_);
00103     sprintf(name,"data/p_%d/grid",myRank_);
00104     cyc2 = -1;
00105   }
00106 
00107   char name [1024];
00108  
00109   void createLinkage () 
00110   {
00111     mpAccess_.removeLinkage();
00112     std::set < int, std::less < int > > s ;
00113     secondScan (s);
00114     mpAccess_.insertRequestSymetric ( s );
00115     //std::cout << "Proc " <<myRank_ << " ="; 
00116     mpAccess_.printLinkage(std::cout);
00117   }
00118 
00119   template <class EntityType> 
00120   int unmarkAllChildren ( EntityType & en )
00121   {
00122     typedef typename EntityType :: Traits :: HierarchicIterator HierIt;
00123     int mxl = grid_.maxlevel();
00124     if(en.isLeaf()) return; 
00125 
00126     en.mark( -1 );
00127     HierIt endit = en.hend( mxl );
00128     for(HierIt it = en.hbegin( mxl ); it != endit ; ++it )
00129     {
00130       it->mark( -1 );
00131     }
00132   }
00133   
00134   typedef std :: map < int , int > OlderElsMap;
00135   std::vector < OlderElsMap > * interiorEls_; 
00136   std::vector < OlderElsMap > * ghostEls_; 
00137  
00138 
00139   std::vector < std:: map < int , ObjectStreamType > > elmap_;
00140   std::vector < std:: map < int , std::map < int , int > > > elmap2_;
00141   
00142   void xtractRefinementTree ()
00143   {
00144     const int nl = mpAccess_.nlinks();
00145 
00146     assert(osv_);
00147     assert(osv_->size() == nl);
00148     
00149     for(int link=0; link<nl; link++)
00150     {
00151     
00152       elmap_[link].clear();
00153       elmap2_[link].clear();
00154     
00155       std:: map < int , ObjectStreamType > & elmap = elmap_[link];
00156       typedef std::map < int , ObjectStreamType > LevelMap;
00157 
00158       ObjectStreamType & os = (*osv_)[link];
00159     
00160       int buff;
00161       int newmxl = 0;
00162       os.readObject( buff );
00163       //std::cout << buff << " Read buff \n";
00164       if(buff == ENDOFSTREAM ) return ;
00165       else
00166       {
00167         assert(buff == BEGINELEMENT );
00168         while( buff == BEGINELEMENT )
00169         {
00170           os.readObject( buff ); // read elnum 
00171           int elnum = buff;
00172           //std::cout << "Unpack el = " << elnum << "\n";
00173           os.readObject(buff); // read refine info  
00174           if(buff == BEGINELEMENT ) continue;
00175           if(buff == ENDOFSTREAM  ) break;
00176           if(buff == 1) // means that macro element has children 
00177           {
00178             //std::cout << "Read level info = " << buff << "\n";
00179             if(elmap.find(elnum) == elmap.end())
00180             {
00181               ObjectStreamType elstr;
00182               elmap[elnum] = elstr;
00183             }
00184             ObjectStreamType & elstr = elmap[elnum];
00185           
00186             os.readObject(buff); // read level 
00187             while((buff != ENDOFSTREAM) && (buff != BEGINELEMENT ))
00188             {
00189               if(buff < 0) newmxl = std::max( newmxl, std::abs( buff ));
00190               elstr.writeObject( buff );
00191               os.readObject( buff );
00192             }
00193           }
00194         }
00195       }
00196       oldmxl_ = std::max( newmxl ,oldmxl_ );
00197     }
00198   }
00199 
00200   bool markFirstLevel () 
00201   {
00202     bool marked = false;
00203     const int nl = mpAccess_.nlinks();
00204     for(int link = 0; link<nl; link++) 
00205     {
00206       typedef std :: map < int , int > HierMap ;
00207       std:: map < int , std::map < int , int > > & elmap2 = elmap2_[link];
00208       std:: map < int , ObjectStreamType > & elmap = elmap_[link];
00209       typedef std::map < int , ObjectStreamType > LevelMap;
00210       {
00211         // now refine grid 
00212         typedef typename GridType :: template Codim<0>::LevelIterator LevelIteratorType;
00213         LevelIteratorType endit  = grid_.template lend<0>  (0);
00214         for(LevelIteratorType it = grid_.template lbegin<0>(0);
00215             it != endit ; ++it )
00216         {
00217           int id = it->globalIndex();
00218           if(elmap.find(id) != elmap.end())
00219           {
00220             std::map < int , int > hiertree;
00221             elmap2[id] = hiertree;
00222             marked = true;
00223             //std::cout << "On p=" << myRank_ << "mark El= " << id << "\n";
00224             if(it->isLeaf()) 
00225             {
00226               grid_.mark(1,(*it));
00227             }
00228           }
00229         }
00230       }
00231     }
00232     nextlevel_ = 1;
00233     return marked;
00234   }
00235 
00236   bool markNextLevel () 
00237   {
00238     if(nextlevel_ > oldmxl_) return false;
00239     bool marked = false;
00240     const int nl = mpAccess_.nlinks();
00241     typedef std :: map < int , int > HierMap ;
00242     for(int link=0; link < nl; link++) 
00243     {
00244       typedef std :: map < int , int > HierMap ;
00245       std:: map < int , std::map < int , int > > & elmap2 = elmap2_[link];
00246       std:: map < int , ObjectStreamType > & elmap = elmap_[link];
00247     
00248       // std::cout << "Begin on Level l = " << mxl << "\n";
00249       // now refine grid 
00250       typedef typename GridType :: template Codim<0>::LevelIterator LevelIteratorType;
00251       LevelIteratorType endit  = grid_.template lend<0>  (0);
00252       for(LevelIteratorType it = grid_.template lbegin<0>(0);
00253           it != endit ; ++it )
00254       {
00255         int id = it->globalIndex();
00256         //std::cout << "Begin LevelIter it = " << id << "\n";
00257         if(elmap.find(id) != elmap.end())
00258         {
00259           int mxl = nextlevel_; 
00260           int buff;
00261           // build a new hier tree 
00262           ObjectStreamType & levstr = elmap[id];
00263           try { 
00264             levstr.readObject( buff );
00265           }
00266           catch (ObjectStreamType :: EOFException) 
00267           {
00268             continue;
00269           }
00270           assert( buff < 0);
00271           assert( std::abs( buff ) == mxl );
00272           
00273           HierMap  & hiertree = elmap2[id];
00274           typedef typename GridType ::template Codim<0> :: Entity :: HierarchicIterator HierIt;
00275 
00276           // Hier muss eine ineinandergeschateltes HierarchiIt kommen.
00277 
00278           typedef typename GridType :: template Codim<0> :: Entity EntityType;
00279           typedef typename GridType :: template Codim<0> :: EntityPointer EntityPointer;
00280           hiertree[id] = 1;
00281 
00282           HierIt hendit = it->hend(mxl);
00283           for(HierIt hit = it->hbegin(mxl); hit != hendit; ++hit)
00284           {
00285             if(hit->level() != mxl) continue;
00286             //std::cout << "Next " << hit->globalIndex() << " on p = " << myRank_ << "\n";
00287             //std::cout << mxl <<  " " << hit->level() << "\n";
00288             // if father isnt in tree then we dont do anything here
00289             EntityPointer vati = hit->father();
00290             //std::cout << "fathe is " << vati.globalIndex() << "\n";
00291             if( hiertree.find( vati->globalIndex() ) ==  hiertree.end()) continue;
00292 
00293             int mark = 0;
00294             try {
00295               levstr.readObject ( mark );
00296             }
00297             catch (ObjectStreamType :: EOFException) 
00298             {
00299               std::cout << "Read " << mark << " from ObjectStream \n";
00300               assert(false);
00301             }
00302 
00303             if(mark == 1)
00304             {
00305               //std::cout << "Mark el = "<<hit->globalIndex()<<" on p="<<myRank_  << "\n";
00306               hiertree[hit->globalIndex()] = mark;
00307               marked = true;
00308               if(hit->isLeaf())
00309               {
00310                 grid_.mark(1,(*hit));
00311               }
00312             }
00313           }
00314         }
00315       }
00316     }
00317 
00318     nextlevel_ ++ ;
00319     return marked;
00320   }
00321 
00322   bool xtractData (DofManagerType & dm) 
00323   {
00324     dm.resize();
00325     
00326     const int nl = mpAccess_.nlinks();
00327     assert(osv_); 
00328     assert(osv_->size() == nl);
00329     for(int link=0; link < nl; link++) 
00330     {
00331       //std::cout << "Read link " << link << "\n";
00332       ObjectStreamType & os = (*osv_)[link]; 
00333       int id = 0;
00334       os.readObject( id );
00335       //std::cout << "first On p=" <<myRank_ << " ident = " << id << "\n";
00336       while ( id != ENDOFSTREAM )
00337       {
00338         typedef std :: map < int , int > HierMap ;
00339         std:: map < int , std::map < int , int > > & elmap2 = elmap2_[link];
00340         //std:: map < int , ObjectStreamType > & elmap = elmap_[link];
00341 
00342         // now refine grid 
00343         typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator LevelIteratorType;
00344         LevelIteratorType endit  = grid_.template lend<0,Interior_Partition>  (0);
00345         for(LevelIteratorType it = grid_.template lbegin<0,Interior_Partition>(0);
00346             it != endit ; ++it )
00347         {
00348           if( id == it->globalIndex())
00349           {
00350             //std::cout << "Read macro el = " << it->globalIndex() << " on p = " << myRank_ << "\n";
00351             // read macro data 
00352             dm.gather( os, *it );
00353             int count =1;
00354 
00355             int mxl = grid_.maxlevel();
00356             HierMap & hiertree = elmap2[id];
00357             typedef typename GridType :: template Codim<0> :: Entity EntityType;
00358             typedef typename GridType :: template Codim<0> :: EntityPointer EntityPointer;
00359               
00360             typedef typename GridType :: template Codim<0> :: Entity :: HierarchicIterator HierIt;
00361             HierIt hendit  = it->hend(mxl);
00362             for(HierIt hit = it->hbegin(mxl); hit != hendit; ++hit)
00363             {
00364               //if(hit->level() != l) continue;
00365               EntityPointer vati = hit->father();
00366               if(hiertree.find(vati->globalIndex()) == hiertree.end()) continue;
00367               //std::cout << "Read el = " << hit->globalIndex() << " on p = " << myRank_ << "\n";
00368               dm.gather(os, *hit );
00369               count ++;
00370             }
00371             //std::cout << "On p="<<myRank_ << " read els = " << count << "\n";
00372           }
00373         }
00374         os.readObject( id );
00375         //std::cout << "On p=" <<myRank_ << " ident = " << id << "\n";
00376       }
00377     }
00378 
00379     std::cout << "Done Start xtract DAta \n";
00380 
00381     coarsenNotNeeded();
00382 
00383     for(int l=0; l < elmap_.size(); l++) 
00384     {
00385       elmap_[l].clear();
00386       elmap2_[l].clear();
00387     }
00388     delete osv_; osv_ = 0;
00389   }
00390 
00391   void coarsenNotNeeded () 
00392   {
00393     std::cout << "Check coarsen on p= " << myRank_ << "\n";
00394     bool notDone = true;
00395     int mxl = grid_.maxlevel();
00396     for(int i=mxl; i > 0 ; i--)
00397     {
00398       typedef typename GridType :: LeafIterator LeafIteratorType;
00399       LeafIteratorType endit  = grid_.leafend  (mxl);
00400       for(LeafIteratorType it = grid_.leafbegin(mxl);
00401           it != endit ; ++it )
00402       {
00403         if((grid_.owner(*it) != grid_.myRank()) && (it->partitionType() != GhostEntity))
00404         {
00405           grid_.mark(-1, (*it)); 
00406         }
00407       }
00408     
00409       notDone = grid_.adapt();
00410     }
00411   }
00412 
00413   template <class EntityType>
00414   int checkRefineStatus ( EntityType & en , PartitionType pitype )
00415   {
00416     typedef typename EntityType :: Traits :: HierarchicIterator HierIt; 
00417     
00418     // if we don have children, do nothing
00419     if(en.isLeaf()) return 0;
00420     
00421     int mxl = grid_.maxlevel();
00422     int count = en.level(); 
00423     HierIt endit = en.hend(mxl); 
00424     for(HierIt it = en.hbegin(mxl); it != endit; ++it)
00425     {
00426       if(it->partitionType() == pitype)
00427         count = std::max(count, (*it).level()); 
00428     }
00429     // return deep below entity
00430     return count - en.level();
00431   }
00432 
00433   template <class EntityType>
00434   void writeChildren ( ObjectStreamType & os, EntityType & en)
00435   {
00436     typedef typename EntityType :: HierarchicIterator HierIt; 
00437     assert( !en.isLeaf() );
00438     
00439     int mxl = en.level() + 1; 
00440     HierIt endit = en.hend(mxl); 
00441     for(HierIt it = en.hbegin(mxl); it != endit; ++it)
00442     {
00443       //if(it->partitionType() != BorderEntity ) continue;
00444       //std::cout << "write State of " << it->globalIndex() << " on p = "<< myRank_ << "\n";
00445       os.writeObject( (it->isLeaf()) ? 0 : 1 );
00446     }
00447     // return deep below entity
00448   }
00449   
00450   template <class EntityType>
00451   void readChildren ( ObjectStreamType & os, EntityType & en )
00452   {
00453     typedef typename EntityType ::  HierarchicIterator HierIt; 
00454     assert( !en.isLeaf() );
00455     
00456     int mxl = en.level() + 1; 
00457     HierIt endit = en.hend(mxl); 
00458     for(HierIt it = en.hbegin(mxl); it != endit; ++it)
00459     {
00460       //if(it->partitionType() != GhostEntity ) continue;
00461       int m = 0;
00462       os.readObject( m );
00463       assert( m != ENDOFSTREAM );
00464       if( m > 0 ) grid_.mark(1, (*it));
00465     }
00466     // return deep below entity
00467   }
00468 
00469   template <class EntityType>
00470   void markChildren ( EntityType & en, int m )
00471   {
00472     typedef typename EntityType ::  HierarchicIterator HierIt; 
00473     if(en.isLeaf()) 
00474     { 
00475       grid_.mark(1,en); 
00476       return ; 
00477     }
00478     
00479     int mxl = grid_.maxlevel();
00480     HierIt endit = en.hend(mxl); 
00481     for(HierIt it = en.hbegin(mxl); it != endit; ++it)
00482     {
00483       if((it->isLeaf()) && (it->partitionType() == GhostEntity ))
00484       {
00485         //std::cout << "on p=" << myRank_ << " mark child " << it->globalIndex() << "\n";
00486         grid_.mark(1,(*it)); 
00487       }
00488     }
00489     // return deep below entity
00490   }
00491   
00492   void unmarkNotOwned()
00493   {
00494     int mxl = grid_.maxlevel();
00495     typedef typename GridType:: LeafIterator LeafIteratorType;
00496     LeafIteratorType endit  = grid_.leafend  (mxl);
00497     for(LeafIteratorType it = grid_.leafbegin(mxl); 
00498         it != endit; ++it)
00499     {
00500       if(grid_.owner(*it) != grid_.myRank() && (it->partitionType() != GhostEntity))
00501       {
00502         grid_.mark(-1,(*it));
00503       }
00504     }
00505   }
00506     
00507   bool firstMark ()
00508   {
00509     cycle_++;
00510     int ts = (cycle_ * 100) + 1;
00511     //GrapeDataIO < GridType > dataIO;
00512     //dataIO.writeGrid( grid_, xdr , name, (double) ts, ts);
00513     
00514     const int nl = mpAccess_.nlinks();
00515     std::vector < ObjectStreamType > osv (nl) ;
00516     std::vector < int > d =  mpAccess_.dest();
00517 
00518     if(interiorEls_) delete interiorEls_;
00519     if(ghostEls_ ) delete ghostEls_;
00520     
00521     interiorEls_ = new std::vector < OlderElsMap > ( nl );
00522     ghostEls_    = new std::vector < OlderElsMap > ( nl );
00523 
00524     //unmarkNotOwned();
00525         
00526     int checkmxl =0;
00527     {
00528       int mxl = grid_.maxlevel();
00529       for(int link=0; link<nl; link++)
00530       {
00531         OlderElsMap & interiorEls = (*interiorEls_)[link];
00532         interiorEls.clear();
00533         int count = 0; 
00534         osv[link].writeObject( mxl );
00535         {
00536           //typedef typename GridType::template 
00537           //  LeafIteratorDef<InteriorBorder_Partition>::LeafIteratorType IBLevelIteratorType;
00538           typedef typename GridType:: LeafIterator IBLevelIteratorType;
00539           IBLevelIteratorType endit  = grid_.template leafend<InteriorBorder_Partition>  ( mxl, d[link] );
00540           for(IBLevelIteratorType it = grid_.template leafbegin<InteriorBorder_Partition>( mxl, d[link] ); 
00541               it != endit; ++it)
00542           {
00543             int id = it->globalIndex();
00544             checkmxl = std::max(checkmxl , it->level());
00545           
00546             int m = std::max ( grid_.getMark(*it) , 0 );
00547             grid_.mark(m,*it);
00548                        
00549             // if element is not marked, check afterwards 
00550             if(interiorEls.find(id) == interiorEls.end()) interiorEls[id] = count;
00551             count ++ ;
00552 
00553             osv[link].writeObject( m );
00554           }
00555         }
00556       }
00557     }
00558 
00559     for(int link=0; link<nl; link++)
00560     {
00561       osv[link].writeObject( checkmxl );
00562       
00563       OlderElsMap & interiorEls = (*interiorEls_)[link];
00564       int s = interiorEls.size();
00565       osv[link].writeObject( s );
00566 
00567       osv[link].writeObject( ENDOFSTREAM );
00568     }
00569    
00570     osv = mpAccess_.exchange(osv);
00571 
00572     // unpack 
00573     int oldmxl = 0;
00574     {
00575       int mxl = grid_.maxlevel();
00576       //assert(mxl == checkmxl);
00577       for(int link=0; link<nl; link++)
00578       {
00579         osv[link].readObject( oldmxl ); 
00580         oldmxl = std::max ( mxl , oldmxl );
00581         
00582         OlderElsMap & ghostEls = (*ghostEls_)[link];
00583         ghostEls.clear();
00584         {
00585           int count = 0;
00586           typedef typename GridType::LeafIterator GLeafIteratorType;
00587           //typedef typename GridType::template LeafIteratorDef<Ghost_Partition>::LeafIteratorType GLeafIteratorType;
00588           GLeafIteratorType endit  = grid_.template leafend<Ghost_Partition> ( mxl, d[link] );
00589           for(GLeafIteratorType it = grid_.template leafbegin<Ghost_Partition> (mxl, d[link] );
00590               it != endit; ++it)
00591           {
00592             int id = it->globalIndex();
00593             if(ghostEls.find(id) == ghostEls.end()) ghostEls[id] = count;
00594             
00595             int m = 0;
00596             osv[link].readObject(m);
00597             assert(m != ENDOFSTREAM );
00598 
00599             if( m == -1 ) 
00600             {
00601               if(ghostEls.find(-id) == ghostEls.end()) ghostEls[-id] = count;
00602             }
00603             m = std::max( m , 0 ); // do not mark for coarsen, first check if coarsen ok 
00604             grid_.mark(m,(*it));
00605             count ++ ;
00606           }
00607         }
00608       }
00609     }
00610     
00611     for(int link=0; link<nl; link++)
00612     {
00613       OlderElsMap & interiorEls = (*interiorEls_)[link];
00614       int buff = 0;
00615       osv[link].readObject( buff );
00616       osv[link].readObject( buff );
00617       //assert(buff == interiorEls.size() );
00618     }
00619   
00620     oldmxl_ = checkmxl;
00621 
00622     return true;   
00623   }
00624 
00625   bool secondMark ()
00626   {
00627     int ts = (cycle_*100) + 2;
00628     //GrapeDataIO < GridType > dataIO;
00629     //dataIO.writeGrid( grid_, xdr , name, (double) ts, ts );
00630 
00631     int oldmxl = oldmxl_;
00632 
00633     const int nl = mpAccess_.nlinks();
00634     std::vector < ObjectStreamType > osv (nl) ;
00635     std::vector < int > d =  mpAccess_.dest();
00636 
00637     //unmarkNotOwned();
00638 
00639     //std::cout << "Starte secondAdapt on proc = " << grid_.myRank() << " with intEls = " << intEls.size() << "\n";
00640 
00641     {
00642       int mxl = oldmxl;
00643       //std::cout << "Iterate over level = " << mxl << "\n";
00644       for(int link=0; link<nl; link++)
00645       {
00646         OlderElsMap & interiorEls = (*interiorEls_)[link];
00647         OlderElsMap & ghostEls    = (*ghostEls_)   [link];
00648         for(int l=0; l<=mxl; l++) 
00649         {
00650           // over all levels 
00651           //typedef typename GridType::template Codim<0>::InteriorBorderLevelIterator IBLevelIteratorType;
00652           typedef typename GridType::template Codim<0>::template Partition<InteriorBorder_Partition>::LevelIterator IBLevelIteratorType;
00653           IBLevelIteratorType endit  = grid_.template lend  <0,InteriorBorder_Partition> ( l, d[link] );
00654           for(IBLevelIteratorType it = grid_.template lbegin<0,InteriorBorder_Partition> ( l, d[link] ); 
00655               it != endit; ++it)
00656           {
00657             int id = it->globalIndex();
00658             // if element is not in former leaf list, then go on 
00659             if( interiorEls.find(id) == interiorEls.end() ) continue;
00660             
00661             int mak = 0;
00662             if( ! (it->isLeaf() ) ) 
00663             {
00664               // this means this element was refined, so send 1 to other
00665               // side
00666               mak = 1;
00667             }
00668               
00669             osv[link].writeObject( interiorEls[id] );
00670             osv[link].writeObject( mak );
00671           }
00672         }
00673 
00674         for(int l=0; l<=mxl; l++) 
00675         {
00676           // over all levels 
00677           typedef typename GridType::template Codim<0>::template Partition<Ghost_Partition>::LevelIterator GLevelIteratorType;
00678           GLevelIteratorType endit  = grid_.template lend  <0,Ghost_Partition> ( l, d[link] );
00679           for(GLevelIteratorType it = grid_.template lbegin<0,Ghost_Partition> ( l, d[link] ); 
00680               it != endit; ++it)
00681           {
00682             int id = it->globalIndex();
00683             
00684             // if element is not in former leaf list, then go on 
00685             if( ghostEls.find(id) == ghostEls.end() ) continue;
00686             
00687             int mak = 0;
00688             if( !(it->isLeaf()) ) 
00689             {
00690               mak = 1;
00691               //std::cout << it->globalIndex() << " write of entity = " << mak << " on p =" <<myRank_ << "\n";
00692             }
00693               
00694             osv[link].writeObject( ghostEls[id] );
00695             osv[link].writeObject( mak );
00696           }
00697         }
00698       }
00699     }
00700     
00701     //std::cout << "SecondAD " << intEls.size () << " " << count << " InteriorBorder Els on p=" << myRank_ << "\n";
00702 
00703     for(int link=0; link<nl; link++)
00704     {
00705       osv[link].writeObject( ENDOFSTREAM );
00706     }
00707     
00708     osv = mpAccess_.exchange(osv);
00709 
00710     std::vector < std::vector < int > > markerIB (nl);
00711     std::vector < std::vector < int > > markerGH (nl);
00712     for(int l=0; l<nl; l++) 
00713     {
00714       OlderElsMap & interiorEls = (*interiorEls_)[l];  
00715       OlderElsMap & ghostEls    = (*ghostEls_)[l];
00716       markerIB[l].resize(interiorEls.size());
00717       markerGH[l].resize(ghostEls.size());
00718       for(int i=0; i<markerIB[l].size(); i++) markerIB[l][i] = -2;
00719       for(int i=0; i<markerGH[l].size(); i++) markerGH[l][i] = -2;
00720     }
00721     
00722     for(int link=0; link<nl; link++)
00723     {
00724       int buff; 
00725       osv[link].readObject( buff );
00726       while (buff != ENDOFSTREAM )
00727       {
00728         int id = buff;
00729         osv[link].readObject( buff );
00730         assert( buff != ENDOFSTREAM );
00731         if(markerGH[link][id] == -2) markerGH[link][id] = buff;
00732         else markerIB[link][id] = buff;
00733         osv[link].readObject( buff );
00734       }
00735     }
00736 
00737     {
00738       int mxl = oldmxl;
00739       for(int link=0; link<nl; link++)
00740       {
00741         OlderElsMap & interiorEls = (*interiorEls_)[link];
00742         OlderElsMap & ghostEls    = (*ghostEls_)[link];
00743         for(int l=0; l<=mxl;l++)
00744         {
00745         typedef typename GridType::template Codim<0>::template Partition<Ghost_Partition>::LevelIterator GLevelIteratorType;
00746         //typedef typename GridType::template LeafIteratorDef<Ghost_Partition>::LeafIteratorType GLevelIteratorType;
00747         GLevelIteratorType endit  = grid_.template lend  <0,Ghost_Partition>(l, d[link] );
00748         for(GLevelIteratorType it = grid_.template lbegin<0,Ghost_Partition>(l, d[link] );
00749             it != endit; ++it)
00750         {
00751           int id = it->globalIndex(); 
00752           if(ghostEls.find(id) == ghostEls.end()) continue;
00753 
00754           int m = markerGH[link][ghostEls[id]];
00755           if((m<=0) && (it->isLeaf())) 
00756           {
00757             ghostEls.erase( id );
00758             continue;
00759           }
00760           else 
00761           {
00762             //if(it->isLeaf())
00763             {
00764               //std::cout << "Mark el =" << it->globalIndex() << " on p = "<< myRank_ << "\n";
00765               grid_.mark(m,(*it));
00766             }
00767           }
00768           
00769           //std::cout << count2 << " on el \n";
00770         }
00771         }
00772 
00773         for(int l=0; l<=mxl;l++)
00774         {
00775           
00776         typedef typename GridType::template Codim<0>::template Partition<InteriorBorder_Partition>::LevelIterator IBLevelIteratorType;
00777         //typedef typename GridType::template LeafIteratorDef<Ghost_Partition>::LeafIteratorType GLevelIteratorType;
00778         IBLevelIteratorType endit  = grid_.template lend  <0,InteriorBorder_Partition>(l, d[link] );
00779         for(IBLevelIteratorType it = grid_.template lbegin<0,InteriorBorder_Partition>(l, d[link] );
00780             it != endit; ++it)
00781         {
00782           int id = it->globalIndex(); 
00783           if(interiorEls.find(id) == interiorEls.end()) continue;
00784 
00785           int m = markerIB[link][interiorEls[id]];
00786           if((m<=0) && (it->isLeaf())) 
00787           {
00788             interiorEls.erase( id );
00789             continue;
00790           }
00791           else 
00792           {
00793             grid_.mark(m,(*it));
00794           }
00795         }
00796         }
00797       }
00798     }
00799    
00800     return true;
00801   }
00802 
00803   bool thirdMark ()
00804   {
00805     int ts = (cycle_*100) + 3;
00806     //GrapeDataIO < GridType > dataIO;
00807     //dataIO.writeGrid( grid_, xdr , name, (double) ts, ts );
00808     
00809     int oldmxl = oldmxl_;
00810     const int nl = mpAccess_.nlinks();
00811     std::vector < ObjectStreamType > osv (nl) ;
00812     std::vector < int > d =  mpAccess_.dest();
00813 
00814     //std::cout << "Starte thirdAdapt on proc = " << grid_.myRank() << " with intEls = " << intEls.size() << "\n";
00815 
00816     int count = 0;
00817     {
00818       int mxl = oldmxl;
00819       //std::cout << "Iterate over level = " << mxl << "\n";
00820       for(int link=0; link<nl; link++)
00821       {
00822         OlderElsMap & interiorEls = (*interiorEls_)[link];
00823         OlderElsMap & ghostEls    = (*ghostEls_)[link];
00824         for(int l=0; l<=mxl; l++) 
00825         {
00826           // over all levels 
00827           typedef typename GridType::template Codim<0>::template Partition<InteriorBorder_Partition>::LevelIterator IBLevelIteratorType;
00828           IBLevelIteratorType endit  = grid_.template lend  <0,InteriorBorder_Partition> ( l, d[link] );
00829           for(IBLevelIteratorType it = grid_.template lbegin<0,InteriorBorder_Partition> ( l, d[link] ); 
00830               it != endit; ++it)
00831           {
00832             int id = it->globalIndex();
00833             
00834             // if element is not in former leaf list, then go on 
00835             if( interiorEls.find(id) == interiorEls.end()) continue;
00836            
00837             writeChildren ( osv[link] , *it );
00838             count++;
00839           }
00840         }
00841       }
00842     }
00843     
00844     //std::cout << "SecondAD " << intEls.size () << " " << count << " InteriorBorder Els on p=" << myRank_ << "\n";
00845 
00846     for(int link=0; link<nl; link++)
00847       osv[link].writeObject( ENDOFSTREAM );
00848     
00849     osv = mpAccess_.exchange(osv);
00850 
00851     // unpack 
00852     int count2 = 0;
00853     {
00854       int mxl = oldmxl;
00855       for(int link=0; link<nl; link++)
00856       {
00857         OlderElsMap & interiorEls = (*interiorEls_)[link];
00858         OlderElsMap & ghostEls    = (*ghostEls_)[link];
00859         for(int l=0; l<=mxl;l++)
00860         {
00861           typedef typename GridType::template Codim<0>::template Partition<Ghost_Partition>::LevelIterator GLevelIteratorType;
00862           GLevelIteratorType endit  = grid_.template lend  <0,Ghost_Partition>(l, d[link] );
00863           for(GLevelIteratorType it = grid_.template lbegin<0,Ghost_Partition>(l, d[link] );
00864               it != endit; ++it)
00865           {
00866             int id = it->globalIndex(); 
00867             if(ghostEls.find(id) == ghostEls.end()) continue;
00868 
00869             //ghostEls.erase( id );
00870             readChildren( osv[link], *it );
00871           }
00872         }
00873       }
00874     }
00875    
00876     return true;
00877   }
00878 
00879   int cyc2; 
00881   bool communicate( DofManagerType & dm)
00882   {
00883     if(cyc2 == cycle_) 
00884     {
00885        cycle_++; 
00886     }
00887     cyc2 = cycle_;
00888     int ts = (cycle_*100) + 4;
00889     //GrapeDataIO < GridType > dataIO;
00890     //dataIO.writeGrid( grid_, xdr , name, (double) ts, ts );
00891 
00892     const int nl = mpAccess_.nlinks();
00893     std::vector < ObjectStreamType > osv (nl) ;
00894     std::vector < int > d =  mpAccess_.dest();
00895 
00896     // pack 
00897 #ifndef NDEBUG
00898     {
00899       for(int link=0; link<nl; link++)
00900       {
00901         int count = 0;
00902         typedef typename GridType::LeafIterator IBLevelIteratorType;
00903         //typedef typename GridType::template LeafIteratorDef<InteriorBorder_Partition>::LeafIteratorType IBLevelIteratorType;
00904         IBLevelIteratorType endit  = grid_.template leafend<InteriorBorder_Partition>   ( grid_.maxlevel(), d[link] );
00905         for(IBLevelIteratorType it = grid_.template leafbegin<InteriorBorder_Partition> ( grid_.maxlevel(), d[link] ); 
00906             it != endit; ++it)
00907         {
00908           count ++ ;
00909         }
00910         osv[link].writeObject( count );
00911       }
00912     }
00913 #endif
00914 
00915     {
00916       for(int link=0; link<nl; link++)
00917       {
00918         typedef typename GridType::LeafIterator  IBLevelIteratorType;
00919         //typedef typename GridType::template LeafIteratorDef<InteriorBorder_Partition>::LeafIteratorType IBLevelIteratorType;
00920         IBLevelIteratorType endit  = grid_.template leafend<InteriorBorder_Partition>   ( grid_.maxlevel(), d[link] );
00921         for(IBLevelIteratorType it = grid_.template leafbegin<InteriorBorder_Partition> ( grid_.maxlevel(), d[link] ); 
00922             it != endit; ++it)
00923         {
00924           dm.scatter(osv[link],*it);
00925         }
00926       }
00927     }
00928 
00929     osv = mpAccess_.exchange(osv);
00930 
00931     dm.resize();
00932 
00933 #ifndef NDEBUG
00934     // unpack 
00935     {
00936       for(int link=0; link<nl; link++)
00937       {
00938         int s = 0;
00939         int count = 0;
00940         osv[link].readObject( s );
00941         typedef typename GridType::LeafIterator GLevelIteratorType;
00942         //typedef typename GridType::template LeafIteratorDef<Ghost_Partition>::LeafIteratorType GLevelIteratorType;
00943         GLevelIteratorType endit  = grid_.template leafend<Ghost_Partition>   ( grid_.maxlevel(), d[link] );
00944         for(GLevelIteratorType it = grid_.template leafbegin<Ghost_Partition> ( grid_.maxlevel(), d[link] );
00945             it != endit; ++it)
00946         {
00947           count ++ ;
00948         }
00949         assert( s == count );
00950       }
00951     }
00952 #endif
00953     
00954     {
00955       for(int link=0; link<nl; link++)
00956       {
00957         typedef typename GridType::LeafIterator GLevelIteratorType;
00958         //typedef typename GridType::template LeafIteratorDef<Ghost_Partition>::LeafIteratorType GLevelIteratorType;
00959         GLevelIteratorType endit  = grid_.template leafend<Ghost_Partition>   ( grid_.maxlevel(), d[link] );
00960         for(GLevelIteratorType it = grid_.template leafbegin<Ghost_Partition> ( grid_.maxlevel(), d[link] );
00961             it != endit; ++it)
00962         {
00963           dm.gather( osv[link], *it);
00964         }
00965       }
00966     }
00967     return true; 
00968   }
00969  
00970   template <class T> 
00971   T globalMax (T val) const 
00972   {
00973     return mpAccess_.gmax(val); 
00974   }
00975     
00976   template <class T> 
00977   T globalMin (T val) const 
00978   {
00979     return mpAccess_.gmin(val); 
00980   }
00981     
00982   template <class T> 
00983   T globalSum (T val) const 
00984   {
00985     return mpAccess_.gsum(val);
00986   }
00987    
00988   template <class T>
00989   void globalSum (T *send, int s , T *recv) const 
00990   {
00991     MPI_Allreduce(send,recv, s, MPIType<T>:: MpiType , MPI_SUM, _mpiComm); 
00992     return ;
00993   }
00994 
00995 private:
00996   int psize () const 
00997   {
00998     return pSize_;
00999   }
01000 
01001 public:
01002   template <class EntityType> 
01003   void ldbUpdateVertex ( EntityType & en , LoadBalancer :: DataBase & db )
01004   {
01005     int weight = 1; // a least weight 1 
01006     {
01007       int mxl = grid_.maxlevel();
01008       typedef typename EntityType :: HierarchicIterator HierIt; 
01009       HierIt endit = en.hend(mxl); 
01010       for(HierIt it = en.hbegin(mxl); it != endit; ++it)
01011       {
01012         weight++;
01013       }
01014       double center[3] = {0.0,0.0,0.0};
01015       db.vertexUpdate ( LoadBalancer :: 
01016           GraphVertex (en.globalIndex(),weight, center) ) ;
01017     }
01018     
01019     {
01020       typedef typename EntityType :: IntersectionIterator InterIt;
01021       InterIt endit = en.iend(); 
01022       for(InterIt nit = en.ibegin(); nit != endit; ++nit)
01023       {
01024         if(nit.neighbor())
01025         if(en.globalIndex() < (*nit).globalIndex())
01026         {
01027           db.edgeUpdate ( LoadBalancer :: GraphEdge ( en.globalIndex(),
01028                 (*nit).globalIndex(), weight) );
01029         }
01030       }
01031     }
01032   }
01033 
01034   bool calcRepartition (DofManagerType & dm) 
01035   { 
01036     std::vector <int> procPart (grid_.size(0,0));
01037     for(int i=0; i<procPart.size() ;i++) procPart[i] = -1;
01038     LoadBalancer :: DataBase db ;
01039     // build up loadbalance data base with macro vertices and edges 
01040     {
01041       typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator InteriorLevelIterator;
01042       InteriorLevelIterator endit  = grid_.template lend   <0,Interior_Partition> (0);
01043       for(InteriorLevelIterator it = grid_.template lbegin <0,Interior_Partition> (0); 
01044           it != endit; ++it )
01045       {
01046         ldbUpdateVertex ( *it , db );
01047       }
01048     }
01049 
01050     //db.printLoad ();
01051     bool neu = false ;
01052     {
01053       const int np = psize ();
01054       // Kriterium, wann eine Lastneuverteilung vorzunehmen ist:
01055       // 
01056       // load  - eigene ElementLast
01057       // mean  - mittlere ElementLast
01058       // nload - Lastverh"altnis
01059 
01060       double load = db.accVertexLoad () ;
01061       std :: vector < double > v (mpAccess_.gcollect (load)) ;
01062       double mean = std::accumulate (v.begin (), v.end (), 0.0) / double (np) ;
01063 
01064       for (std :: vector < double > :: iterator i = v.begin () ; i != v.end () ; i ++)
01065         neu |= (*i > mean ? (*i > (_ldbOver * mean) ? true : false) : 
01066             (*i < (_ldbUnder * mean) ? true : false)) ;
01067 
01068       //std::cout << load << " Load \n";
01069     }
01070     int val = (neu) ? 1 : 0; 
01071     int v2  = mpAccess_.gmax(val); 
01072     neu = (v2 == 1) ? true : false;
01073     
01074     if (neu) 
01075     {
01076       db.repartition( mpAccess_ ,  LoadBalancer :: DataBase :: METIS_PartGraphRecursive );
01077       {
01078         int countMyEls = 0;
01079         int firstEl = 0;
01080         typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator InteriorLevelIterator;
01081         InteriorLevelIterator endit= grid_.template lend  <0,Interior_Partition> (0);
01082         InteriorLevelIterator it   = grid_.template lbegin<0,Interior_Partition> (0); 
01083         if(it != endit) firstEl = it->globalIndex(); 
01084         for(; it != endit; ++it )
01085         {
01086           int id = (*it).globalIndex();
01087           procPart[id] = db.getDestination ( id ) ;
01088           if(procPart[id] == grid_.myRank()) countMyEls++;
01089         }
01090         //assert(countMyEls > 0);
01091         if(countMyEls == 0) procPart[firstEl] = grid_.myRank();
01092 
01093       }
01094       
01095       //std::cout << "Do repartition of MacroGrid on p= "<<myRank_ <<"\n";
01096       repartitionMacroGrid( procPart , dm );
01097       return true;
01098     }
01099 
01100     return false;
01101   }
01102 
01103   void  secondScan ( std::set < int, std::less < int > >  & s ) 
01104   {
01105     {
01106       // fake , all to all linkage at the moment 
01107       {
01108         std :: vector < int > l (pSize_-1);
01109         int count =0;
01110         for(int i=0; i<pSize_; i++)
01111         {
01112           if(i!=myRank_) 
01113           {
01114             l[count] = i;
01115             count ++;
01116           }
01117         }
01118         
01119         for (std::vector < int > :: const_iterator i = l.begin () ; i != l.end () ; s.insert (* i ++)) ;
01120         int pos = *(s.find(1-myRank_)); 
01121         //std::cout << " link = " << pos << "\n";
01122       }
01123     }
01124   }  
01125 
01126   std::vector < ObjectStreamType > * osv_ ;
01127   
01128   void repartitionMacroGrid ( std::vector<int> & procPart , DofManagerType & dm )
01129   {
01130     for( int l = grid_.maxlevel(); l>0; l-- )
01131     {
01132       unmarkNotOwned();
01133       grid_.adapt();
01134     }
01135 
01136     
01137     std::vector < int > d = mpAccess_.dest(); 
01138     const int nlinks = mpAccess_.nlinks ();
01139 
01140     if(osv_) delete osv_;
01141     osv_ = new std::vector < ObjectStreamType > (nlinks);
01142     std::vector < ObjectStreamType > & osv = *osv_;
01143 
01144     {
01145       // write new element owners, write the number only if you are the owner
01146       {
01147         typedef typename GridType::template Codim<0> :: LevelIterator LevelIteratorType;
01148         LevelIteratorType it    = grid_.template lbegin<0> ( 0 );
01149         LevelIteratorType endit = grid_.template lend<0>   ( 0 );
01150         for( ; it != endit; ++it)
01151         {
01152           for(int p=0;p<nlinks; p++)
01153           {
01154             if(grid_.owner(*it) == grid_.myRank())
01155             {
01156               osv[p].writeObject( procPart[ (*it).globalIndex()] );
01157             }
01158             else
01159             {
01160               osv[p].writeObject( -1 );
01161             }
01162           }
01163         }
01164       }
01165       
01166       {
01167         // walk over my interior macro elements
01168         typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator InteriorLevelIteratorType;
01169         InteriorLevelIteratorType it    = grid_.template lbegin<0,Interior_Partition> ( 0 );
01170         InteriorLevelIteratorType endit = grid_.template lend<0,Interior_Partition> ( 0 );
01171         for( ; it != endit; ++it)
01172         {
01173           if(grid_.owner(*it) == grid_.myRank())
01174           {
01175             int id = it->globalIndex();
01176             if(procPart[id] != grid_.myRank() && (procPart[id] != -1))
01177             {
01178               //std::cout << grid_.myRank() << " mr|pp " << procPart[id] <<"\n";
01179               //std::cout << "Pack Element ! = " << id << "\n";
01180               int newProc = mpAccess_.link ( procPart[id] );
01181               grid_.packAll ( osv[newProc] , (*it) );
01182             }
01183           }
01184         }
01185       }
01186       
01187       // set new owners 
01188       for(int p=0; p<nlinks; p++)
01189       {
01190         osv[p].writeObject(ENDOFSTREAM); // end stream
01191       }
01192 
01193       // pack data 
01194       {
01195         // walk over my interior macro elements
01196         typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator InteriorLevelIteratorType;
01197         InteriorLevelIteratorType it    = grid_.template lbegin<0,Interior_Partition> ( 0 );
01198         InteriorLevelIteratorType endit = grid_.template lend<0,Interior_Partition> ( 0 );
01199         for( ; it != endit; ++it)
01200         {
01201           if(grid_.owner(*it) == grid_.myRank())
01202           {
01203             int id = it->globalIndex();
01204             if(procPart[id] != grid_.myRank() && (procPart[id] != -1))
01205             {
01206               int link = mpAccess_.link ( procPart[id] );
01207               grid_.partition( procPart[id] , (*it) );
01208 
01209               osv[link].writeObject( id );
01210               packAllData ( osv[link], dm , *it );
01211               //dm.inlineData( osv[link], *it ); 
01212             }
01213           }
01214         }
01215       }
01216 
01217       for(int p=0; p<nlinks; p++)
01218       {
01219         osv[p].writeObject(ENDOFSTREAM); // end stream
01220       }
01221     }
01222 
01223     osv = mpAccess_.exchange( osv );
01224 
01225     { 
01226       typedef typename GridType::template Codim<0> :: LevelIterator LevelIteratorType;
01227       LevelIteratorType it    = grid_.template lbegin<0> ( 0 );
01228       LevelIteratorType endit = grid_.template lend<0>   ( 0 );
01229       for( ; it != endit; ++it)
01230       {
01231         for(int p=0; p<nlinks; p++)
01232         {
01233           int proc = d[p];
01234           int np = 0 ;
01235           osv[p].readObject( np );
01236           //std::cout << myRank_ << "mr | " <<np << "\n";
01237           //std::cout << it->owner() << "own | " << proc << "\n";
01238 
01239           // the new number comes only from the old owner 
01240           if((grid_.owner(*it) == proc) && (np >= 0)) 
01241           {
01242             grid_.partition ( np , (*it)); 
01243           }
01244         }
01245       }
01246     }
01247 
01248     { 
01249       typedef typename GridType::template Codim<0>::template Partition<Interior_Partition>::LevelIterator InteriorLevelIteratorType;
01250       InteriorLevelIteratorType it    = grid_.template lbegin<0,Interior_Partition> ( 0 );
01251       InteriorLevelIteratorType endit = grid_.template lend<0,Interior_Partition>   ( 0 );
01252       for( ; it != endit; ++it)
01253       {
01254         // the new number comes only from the old owner 
01255         if((grid_.owner(*it) == grid_.myRank())) 
01256         {
01257           int id = (*it).globalIndex();
01258           //std::cout << "Setze part von el = " << id << " auf Proc = " << myRank_ << "\n";
01259           //std::cout << procPart[id] << " Procpart = \n";
01260           if(procPart[id] >= 0)
01261             grid_.partition ( procPart[ id ] , (*it) ); 
01262         }
01263       }
01264     }
01265     
01266     xtractRefinementTree();
01267     //ts++;
01268     //dataIO.writeGrid( grid_, xdr , name, (double) ts, ts );
01269   }
01270   
01271   template <class EntityType>
01272   void packAllData ( ObjectStreamType & os, DofManagerType & dm , EntityType & en )
01273   {
01274     typedef typename EntityType :: HierarchicIterator HierIt; 
01275     dm.scatter( os , en );
01276     //std::cout << "Write Macro el = " << en.globalIndex()<<" on p=" <<myRank_ <<"\n";
01277 
01278     int count = 1;
01279     if(! (en.isLeaf() ))
01280     {
01281       //int mxl = en.level() + 1;
01282       int mxl = grid_.maxlevel();
01283       HierIt endit  = en.hend(mxl); 
01284       for(HierIt it = en.hbegin(mxl); it != endit; ++it)
01285       {
01286         //std::cout << "Pack el = " << it->globalIndex() << " on p=" << myRank_ <<"\n";
01287         dm.scatter( os, *it );
01288         count ++;
01289       }
01290     }
01291     //std::cout << "Packed "<<count<< " El on p=" << myRank_ <<"\n";
01292   }
01293 
01294   bool consistencyGhosts ()
01295   {
01296     std::vector < int > d = mpAccess_.dest(); 
01297     const int nlinks = mpAccess_.nlinks ();
01298     std::vector < ObjectStreamType > osv (nlinks) ;
01299 
01300     {
01301       for(int link=0; link<nlinks; link++)
01302       {
01303         // walk over my interior macro elements
01304         typedef typename GridType::template Codim<0>::template Partition<InteriorBorder_Partition>::LevelIterator IBLevelIteratorType;
01305         IBLevelIteratorType it    = grid_.template lbegin<0,InteriorBorder_Partition>( 0, d[link] );
01306         IBLevelIteratorType endit = grid_.template lend  <0,InteriorBorder_Partition>( 0, d[link] );
01307         for( ; it != endit; ++it)
01308         {
01309           grid_.packBorder ( osv[link] , (*it) );
01310         }
01311       }
01312     } 
01313     
01314     // set new owners 
01315     for(int p=0; p<nlinks; p++)
01316     {
01317       osv[p].writeObject(ENDOFSTREAM); // end stream
01318     }
01319 
01320     osv = mpAccess_.exchange( osv );
01321 
01322     {
01323       for(int link=0; link<nlinks; link++)
01324       {
01325         grid_.unpackAll ( osv[link] ); 
01326       }
01327     }
01328     
01329     grid_.createGhosts ();
01330     return true;
01331   }
01332 
01333   
01334   bool consistencyAdapt ()
01335   {
01336     assert(false);
01337     std::vector < int > d = mpAccess_.dest(); 
01338     const int nlinks = mpAccess_.nlinks ();
01339     std::vector < ObjectStreamType > osv (nlinks) ;
01340 
01341     {
01342       for(int link=0; link<nlinks; link++)
01343       {
01344         int mxl = grid_.maxlevel(); 
01345         // walk over my interior macro elements
01346         typedef typename GridType::LeafIterator IBLeafIteratorType;
01347         IBLeafIteratorType it    = grid_.template leafbegin<InteriorBorder_Partition>( mxl, d[link] );
01348         IBLeafIteratorType endit = grid_.template leafend  <InteriorBorder_Partition>( mxl, d[link] );
01349         for( ; it != endit; ++it)
01350         {
01351           osv[link].writeObject( it->level() );
01352         }
01353       }
01354     } 
01355     
01356     // set new owners 
01357     for(int p=0; p<nlinks; p++)
01358     {
01359       osv[p].writeObject(ENDOFSTREAM); // end stream
01360     }
01361 
01362     {
01363       for(int link=0; link<nlinks; link++)
01364       {
01365         int mxl = grid_.maxlevel(); 
01366         // walk over my interior macro elements
01367         typedef typename GridType::LeafIterator GLeafIteratorType;
01368         GLeafIteratorType it    = grid_.template leafbegin<Ghost_Partition>( mxl, d[link] );
01369         GLeafIteratorType endit = grid_.template leafend  <Ghost_Partition>( mxl, d[link] );
01370         for( ; it != endit; ++it)
01371         {
01372           int lev = 0;
01373           osv[link].readObject( lev );
01374           assert( lev != ENDOFSTREAM );
01375           if(lev > it->level()) it->mark(1);
01376         }
01377       }
01378     } 
01379    
01380     grid_.adapt();
01381     return true;
01382   }
01383 
01384   bool repartition (DofManagerType & dm) 
01385   {
01386     return calcRepartition (dm);
01387   } 
01388 
01389   // reference to corresponding grid 
01390   GridType & grid_;
01391 
01392   MPI_Comm _mpiComm;
01393 
01394   // rank of my thread 
01395   const int myRank_;
01396   const int pSize_;
01397 
01398   int oldmxl_; 
01399   int nextlevel_; 
01400 
01401   MpAccessMPI mpAccess_;
01402 
01403   double _ldbOver; 
01404   double _ldbUnder; 
01405 
01406 };
01407 #else 
01408 
01409 class ObjectStream 
01410 {
01411   public:
01412     class EOFException {} ;
01413     template <class T>
01414     void readObject (T &) {}
01415     void readObject (int) {}
01416     void readObject (double) {}
01417     template <class T> 
01418     void writeObject (T &) {}
01419     void writeObject (int) {}
01420     void writeObject (double) {}
01421     
01422     template <class T>
01423     void read (T &) const {}
01424     template <class T> 
01425     void write (const T &) {}
01426 };
01427 
01428 typedef ObjectStream ObjectStreamType;
01429 typedef ObjectStream AlbertaObjectStream;
01430 
01431 template <class GridType, class DofManagerType>
01432 class AlbertGridCommunicator : public CommunicatorInterface<DofManagerType>
01433 {
01434 public: 
01435 
01436 #if defined(HAVE_MPI) && defined(ALBERTA_USES_MPI)
01437   AlbertGridCommunicator(MPI_Comm mpiComm, GridType &grid, int anzp) {}
01438 #else 
01439   AlbertGridCommunicator(GridType &grid) {}
01440 #endif 
01441 
01442   template <class T> 
01443   T globalMax (T val) const 
01444   {
01445     return val;
01446   }
01447     
01448   template <class T> 
01449   T globalMin (T val) const 
01450   {
01451     return val;
01452   }
01453     
01454   template <class T> 
01455   T globalSum (T val) const 
01456   {
01457     return val;
01458   }
01459    
01460   template <class T>
01461   void globalSum (T *send, int s , T *recv) const 
01462   {
01463     std::memcpy(recv,send,s*sizeof(T));
01464     return ;
01465   }
01466   
01467   bool loadBalance ( DofManagerType & dm) { return false; }
01468   void loadBalance ()  { return false; }
01469   
01470   bool firstMark() { return false; }
01471   bool secondMark() { return false; }
01472   bool thirdMark() { return false;}
01473   
01474   bool markFirstLevel() { return false; }
01475   bool markNextLevel () { return false; }
01476   
01477   bool xtractData (DofManagerType & dm) { return false; }
01478   
01479   bool repartition (DofManagerType & dm) { return false; }
01480   bool communicate (DofManagerType & dm) { return false; }
01481   bool consistencyGhosts () { return false; }
01482   
01483 };
01484 #endif
01485 
01486 template <class GridType, class CritType>
01487 void makeParallelGrid (GridType &grid, CritType &crit)
01488 {
01489   for(int l=0; l <= grid.maxlevel(); l++)
01490   {
01491     typedef typename GridType::template Codim<0>::LevelIterator LevelIteratorType;
01492     
01493     LevelIteratorType it    = grid.template lbegin<0> (l);
01494     LevelIteratorType endit = grid.template lend  <0> (l);
01495  
01496     for( ; it != endit; ++it )
01497     {
01498       crit.classify( *it ); 
01499     }
01500   }
01501 }
01502 
01503 } // end namespace Dune 
01504 #endif

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