agcommunicator.hh

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

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