00001 #ifndef DUNE_ALBERTAGRID_IMP_HH
00002 #define DUNE_ALBERTAGRID_IMP_HH
00003
00009 #if HAVE_ALBERTA
00010
00011 #include <iostream>
00012 #include <fstream>
00013 #include <dune/common/deprecated.hh>
00014
00015 #include <vector>
00016 #include <assert.h>
00017 #include <algorithm>
00018
00019
00020 #include <dune/common/misc.hh>
00021 #include <dune/common/interfaces.hh>
00022 #include <dune/common/fvector.hh>
00023 #include <dune/common/fmatrix.hh>
00024 #include <dune/common/stdstreams.hh>
00025 #include <dune/common/collectivecommunication.hh>
00026
00027 #include <dune/grid/common/grid.hh>
00028 #include <dune/grid/common/adaptcallback.hh>
00029 #include <dune/grid/common/defaultindexsets.hh>
00030 #include <dune/grid/common/sizecache.hh>
00031 #include <dune/grid/common/intersectioniteratorwrapper.hh>
00032 #include <dune/grid/common/defaultgridview.hh>
00033
00034
00035
00036 #include "albertaheader.hh"
00037
00038
00039 #include <dune/grid/utility/grapedataioformattypes.hh>
00040
00041
00042 #define CALC_COORD 0
00043
00044
00045 #include "albertaextra.hh"
00046
00047 #include <dune/grid/albertagrid/misc.hh>
00048 #include <dune/grid/albertagrid/capabilities.hh>
00049
00050
00051 #include "agmemory.hh"
00052
00053 #include <dune/grid/albertagrid/coordcache.hh>
00054 #include <dune/grid/albertagrid/level.hh>
00055
00056 #include "indexsets.hh"
00057 #include "geometry.hh"
00058 #include "entity.hh"
00059 #include "entitypointer.hh"
00060 #include "hierarchiciterator.hh"
00061 #include "treeiterator.hh"
00062 #include "leveliterator.hh"
00063 #include "leafiterator.hh"
00064 #include "intersection.hh"
00065
00066 namespace Dune
00067 {
00068
00069
00070
00071
00072 template< int dim, int dimworld >
00073 class AlbertaGrid;
00074
00075
00076
00077
00078
00079
00080 template <int dim, int dimworld>
00081 struct AlbertaGridFamily
00082 {
00083 typedef AlbertaGrid<dim,dimworld> GridImp;
00084
00085 typedef DefaultLevelIndexSet< AlbertaGrid<dim,dimworld> > LevelIndexSetImp;
00086 typedef DefaultLeafIndexSet< AlbertaGrid<dim,dimworld> > LeafIndexSetImp;
00087
00088 typedef AlbertaGridIdSet< dim, dimworld > IdSetImp;
00089 typedef unsigned int IdType;
00090
00091 struct Traits
00092 {
00093 typedef GridImp Grid;
00094
00095 typedef Dune :: Intersection< const GridImp, LeafIntersectionIteratorWrapper > LeafIntersection;
00096 typedef Dune :: Intersection< const GridImp, LeafIntersectionIteratorWrapper > LevelIntersection;
00097 typedef Dune::IntersectionIterator<const GridImp, LeafIntersectionIteratorWrapper, LeafIntersectionIteratorWrapper > LeafIntersectionIterator;
00098 typedef Dune::IntersectionIterator<const GridImp, LeafIntersectionIteratorWrapper, LeafIntersectionIteratorWrapper > LevelIntersectionIterator;
00099
00100 typedef Dune::HierarchicIterator<const GridImp, AlbertaGridHierarchicIterator> HierarchicIterator;
00101
00102 typedef IdType GlobalIdType;
00103 typedef IdType LocalIdType;
00104
00105 template <int cd>
00106 struct Codim
00107 {
00108
00109 typedef Dune::Geometry<dim-cd, dimworld, const GridImp, AlbertaGridGeometry> Geometry;
00110 typedef Dune::Geometry<dim-cd, dim, const GridImp, AlbertaGridGeometry> LocalGeometry;
00111
00112 typedef Dune::Entity< cd, dim, const GridImp, AlbertaGridEntity > Entity;
00113
00114 typedef AlbertaGridEntityPointer< cd, const GridImp > EntityPointerImpl;
00115 typedef Dune::EntityPointer< const GridImp, EntityPointerImpl > EntityPointer;
00116
00117 template <PartitionIteratorType pitype>
00118 struct Partition
00119 {
00120 typedef Dune::LevelIterator<cd,pitype,const GridImp,AlbertaGridLevelIterator> LevelIterator;
00121 typedef Dune::LeafIterator<cd,pitype,const GridImp,AlbertaGridLeafIterator> LeafIterator;
00122 };
00123
00124 typedef typename Partition< All_Partition >::LevelIterator LevelIterator;
00125 typedef typename Partition< All_Partition >::LeafIterator LeafIterator;
00126 };
00127
00128 template <PartitionIteratorType pitype>
00129 struct Partition
00130 {
00131 typedef Dune::GridView<DefaultLevelGridViewTraits<const GridImp,pitype> >
00132 LevelGridView;
00133 typedef Dune::GridView<DefaultLeafGridViewTraits<const GridImp,pitype> >
00134 LeafGridView;
00135 };
00136
00137 typedef IndexSet<GridImp,LevelIndexSetImp,DefaultLevelIteratorTypes<GridImp> > LevelIndexSet;
00138 typedef IndexSet<GridImp,LeafIndexSetImp,DefaultLeafIteratorTypes<GridImp> > LeafIndexSet;
00139 typedef AlbertaGridHierarchicIndexSet< dim, dimworld > HierarchicIndexSet;
00140 typedef IdSet<GridImp,IdSetImp,IdType> GlobalIdSet;
00141 typedef IdSet<GridImp,IdSetImp,IdType> LocalIdSet;
00142
00143 typedef Dune::CollectiveCommunication< int > CollectiveCommunication;
00144 };
00145 };
00146
00147
00148
00149
00150
00151
00192 template< int dim, int dimworld = Alberta::dimWorld >
00193 class AlbertaGrid
00194 : public GridDefaultImplementation
00195 < dim, dimworld, Alberta::Real, AlbertaGridFamily< dim, dimworld > >,
00196 public HasHierarchicIndexSet
00197 {
00198 typedef AlbertaGrid< dim, dimworld > This;
00199 typedef GridDefaultImplementation
00200 < dim, dimworld, Alberta::Real, AlbertaGridFamily< dim, dimworld > >
00201 Base;
00202
00203
00204 template< class, class > friend class Conversion;
00205
00206 template< int, int, class > friend class AlbertaGridEntity;
00207
00208 friend class AlbertaGridHierarchicIterator<AlbertaGrid<dim,dimworld> >;
00209
00210 friend class AlbertaGridIntersectionIterator<AlbertaGrid<dim,dimworld> >;
00211 friend class AlbertaGridIntersectionIterator<const AlbertaGrid<dim,dimworld> >;
00212
00213 friend class AlbertaMarkerVector< dim, dimworld >;
00214 friend class AlbertaGridHierarchicIndexSet<dim,dimworld>;
00215
00216 public:
00217 typedef Alberta::Real ctype;
00218
00219 static const int dimension = dim;
00220 static const int dimensionworld = dimworld;
00221
00223 typedef AlbertaGridFamily< dim, dimworld > GridFamily;
00224
00225
00226 typedef typename AlbertaGridFamily< dim, dimworld >::Traits Traits;
00227
00229 typedef typename Traits::HierarchicIndexSet HierarchicIndexSet;
00230
00232 typedef typename Traits::CollectiveCommunication CollectiveCommunication;
00233
00234 private:
00236 typedef typename Traits::template Codim<0>::LeafIterator LeafIterator;
00237
00239 typedef typename GridFamily:: LevelIndexSetImp LevelIndexSetImp;
00240 typedef typename GridFamily:: LeafIndexSetImp LeafIndexSetImp;
00241
00243 typedef typename Traits :: LeafIndexSet LeafIndexSet;
00244
00246 typedef AlbertaGridIdSet<dim,dimworld> IdSetImp;
00247 typedef typename Traits :: GlobalIdSet GlobalIdSet;
00248 typedef typename Traits :: LocalIdSet LocalIdSet;
00249
00250 struct AdaptationState;
00251
00252 template< class DataHandler >
00253 struct AdaptationCallback;
00254
00255 public:
00257 typedef typename ALBERTA AlbertHelp::AlbertLeafData<dimworld,dim+1> LeafDataType;
00258
00259 private:
00260
00261 static const int MAXL = 64;
00262
00263 typedef Alberta::ElementInfo< dimension > ElementInfo;
00264 typedef Alberta::MeshPointer< dimension > MeshPointer;
00265 typedef Alberta::HierarchyDofNumbering< dimension > DofNumbering;
00266 typedef AlbertaGridLevelProvider< dimension > LevelProvider;
00267
00268
00269 AlbertaGrid ( const This & );
00270 This &operator= ( const This & );
00271
00272 public:
00274 AlbertaGrid ();
00275
00281 AlbertaGrid ( const Alberta::MacroData< dimension > ¯oData,
00282 const std::string &gridName = "AlbertaGrid" );
00283
00289 AlbertaGrid ( const std::string ¯oGridFileName,
00290 const std::string &gridName = "AlbertaGrid" );
00291
00293 ~AlbertaGrid ();
00294
00297 int maxLevel () const;
00298
00300 template<int cd, PartitionIteratorType pitype>
00301 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
00302 lbegin (int level) const;
00303
00305 template<int cd, PartitionIteratorType pitype>
00306 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
00307 lend (int level) const;
00308
00310 template< int codim >
00311 typename Traits::template Codim< codim >::LevelIterator
00312 lbegin ( int level ) const;
00313
00315 template< int codim >
00316 typename Traits::template Codim< codim >::LevelIterator
00317 lend ( int level ) const;
00318
00320 template< int codim, PartitionIteratorType pitype >
00321 typename Traits
00322 ::template Codim< codim >::template Partition< pitype >::LeafIterator
00323 leafbegin () const;
00324
00326 template< int codim, PartitionIteratorType pitype >
00327 typename Traits
00328 ::template Codim< codim >::template Partition< pitype >::LeafIterator
00329 leafend () const;
00330
00332 template< int codim >
00333 typename Traits::template Codim< codim >::LeafIterator
00334 leafbegin () const;
00335
00337 template< int codim >
00338 typename Traits::template Codim< codim >::LeafIterator
00339 leafend () const;
00340
00345 int size (int level, int codim) const;
00346
00348 int size (int level, GeometryType type) const;
00349
00351 int size (int codim) const;
00352
00354 int size (GeometryType type) const;
00355
00356 public:
00357
00358
00359
00360 using Base::getMark;
00361 using Base::mark;
00362
00364 int getMark ( const typename Traits::template Codim< 0 >::Entity &e ) const;
00365
00367 bool mark ( int refCount, const typename Traits::template Codim< 0 >::Entity &e );
00368
00370 bool globalRefine ( int refCount );
00371
00372 template< class DataHandle >
00373 bool globalRefine ( int refCount, AdaptDataHandleInterface< This, DataHandle > &handle );
00374
00376 bool adapt ();
00377
00378 template< class DataHandle >
00379 bool adapt ( AdaptDataHandleInterface< This, DataHandle > &handle );
00380
00382 template< class DofManager, class RestrictProlongOperator >
00383 bool DUNE_DEPRECATED
00384 adapt ( DofManager &, RestrictProlongOperator &, bool verbose = false );
00385
00387 bool preAdapt ();
00388
00390 void postAdapt();
00391
00394 const CollectiveCommunication &comm () const
00395 {
00396 return comm_;
00397 }
00398
00399 static std::string typeName ()
00400 {
00401 std::ostringstream s;
00402 s << "AlbertaGrid< " << dim << ", " << dimworld << " >";
00403 return s.str();
00404 }
00405
00407 std::string name () const
00408 {
00409 return mesh_.name();
00410 };
00411
00412
00413
00414
00416 template< GrapeIOFileFormatType ftype >
00417 bool writeGrid( const std::string &filename, ctype time ) const;
00418
00420 template< GrapeIOFileFormatType ftype >
00421 bool readGrid( const std::string &filename, ctype &time );
00422
00423
00424 const HierarchicIndexSet & hierarchicIndexSet () const { return hIndexSet_; }
00425
00427 const typename Traits :: LevelIndexSet & levelIndexSet (int level) const;
00428
00430 const typename Traits :: LeafIndexSet & leafIndexSet () const;
00431
00433 const GlobalIdSet &globalIdSet () const
00434 {
00435 return idSet_;
00436 }
00437
00439 const LocalIdSet &localIdSet () const
00440 {
00441 return idSet_;
00442 }
00443
00444
00445 ALBERTA MESH* getMesh () const
00446 {
00447 return mesh_;
00448 };
00449
00450 const MeshPointer &meshPointer () const
00451 {
00452 return mesh_;
00453 }
00454
00455 const DofNumbering &dofNumbering () const
00456 {
00457 return dofNumbering_;
00458 }
00459
00460 const LevelProvider &levelProvider () const
00461 {
00462 return levelProvider_;
00463 }
00464
00465 int dune2alberta ( int codim, int i ) const
00466 {
00467 return numberingMap_.dune2alberta( codim, i );
00468 }
00469
00470 int alberta2dune ( int codim, int i ) const
00471 {
00472 return numberingMap_.alberta2dune( codim, i );
00473 }
00474
00475 private:
00476 using Base::getRealImplementation;
00477
00478 typedef std::vector<int> ArrayType;
00479
00480 void setup ();
00481
00482
00483
00484 void calcExtras();
00485
00486
00487 bool writeGridXdr ( const std::string &filename, ctype time ) const;
00488
00490 bool readGridXdr ( const std::string &filename, ctype &time );
00491
00492 #if 0
00494 bool readGridAscii ( const std::string &filename, ctype &time );
00495 #endif
00496
00497
00498 void removeMesh();
00499
00500
00501 MeshPointer mesh_;
00502
00503
00504 CollectiveCommunication comm_;
00505
00506
00507 int maxlevel_;
00508
00509
00510
00511
00512 typedef MakeableInterfaceObject< typename Traits::template Codim< 0 >::Entity >
00513 EntityObject;
00514
00515 public:
00516 typedef AGMemoryProvider< EntityObject > EntityProvider;
00517
00518 typedef AlbertaGridIntersectionIterator< const This > IntersectionIteratorImp;
00519 typedef IntersectionIteratorImp LeafIntersectionIteratorImp;
00520 typedef AGMemoryProvider< LeafIntersectionIteratorImp > LeafIntersectionIteratorProviderType;
00521 friend class LeafIntersectionIteratorWrapper< const This >;
00522
00523 typedef LeafIntersectionIteratorWrapper< const This >
00524 AlbertaGridIntersectionIteratorType;
00525
00526 LeafIntersectionIteratorProviderType & leafIntersetionIteratorProvider() const { return leafInterItProvider_; }
00527
00528 private:
00529 mutable EntityProvider entityProvider_;
00530 mutable LeafIntersectionIteratorProviderType leafInterItProvider_;
00531
00532 public:
00533 template< class IntersectionInterfaceType >
00534 const typename Base
00535 :: template ReturnImplementationType< IntersectionInterfaceType >
00536 :: ImplementationType & DUNE_DEPRECATED
00537 getRealIntersectionIterator ( const IntersectionInterfaceType &iterator ) const
00538 {
00539 return this->getRealImplementation( iterator );
00540 }
00541
00542 template< class IntersectionType >
00543 const typename Base
00544 :: template ReturnImplementationType< IntersectionType >
00545 :: ImplementationType &
00546 getRealIntersection ( const IntersectionType &intersection ) const
00547 {
00548 return this->getRealImplementation( intersection );
00549 }
00550
00551
00552 template< int codim >
00553 MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *
00554 getNewEntity () const;
00555
00556
00557 template <int codim>
00558 void freeEntity ( MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *entity ) const;
00559
00560 public:
00561
00562 const Alberta::GlobalVector &
00563 getCoord ( const ElementInfo &elementInfo, int vertex ) const;
00564
00565 private:
00566
00567 Alberta::NumberingMap< dimension > numberingMap_;
00568
00569 DofNumbering dofNumbering_;
00570
00571 LevelProvider levelProvider_;
00572
00573
00574 HierarchicIndexSet hIndexSet_;
00575
00576
00577 IdSetImp idSet_;
00578
00579
00580
00581 mutable std::vector< LevelIndexSetImp * > levelIndexVec_;
00582
00583
00584
00585 mutable LeafIndexSetImp* leafIndexSet_;
00586
00587 typedef SingleTypeSizeCache< This > SizeCacheType;
00588 SizeCacheType * sizeCache_;
00589
00590 typedef AlbertaMarkerVector< dim, dimworld > MarkerVector;
00591
00592
00593 mutable MarkerVector leafMarkerVector_;
00594
00595
00596 mutable std::vector< MarkerVector > levelMarkerVector_;
00597
00598 #if !CALC_COORD
00599 Alberta::CoordCache< dimension > coordCache_;
00600 #endif
00601
00602
00603 AdaptationState adaptationState_;
00604 };
00605
00606
00607
00608
00609
00610
00611 template< int dim, int dimworld >
00612 struct AlbertaGrid< dim, dimworld >::AdaptationState
00613 {
00614 enum Phase { ComputationPhase, PreAdaptationPhase, PostAdaptationPhase };
00615
00616 private:
00617 Phase phase_;
00618 int coarsenMarked_;
00619 int refineMarked_;
00620
00621 public:
00622 AdaptationState ()
00623 : phase_( ComputationPhase ),
00624 coarsenMarked_( 0 ),
00625 refineMarked_( 0 )
00626 {}
00627
00628 void mark ( int count )
00629 {
00630 if( count < 0 )
00631 ++coarsenMarked_;
00632 if( count > 0 )
00633 refineMarked_ += (2 << count);
00634 }
00635
00636 void unmark ( int count )
00637 {
00638 if( count < 0 )
00639 --coarsenMarked_;
00640 if( count > 0 )
00641 refineMarked_ -= (2 << count);
00642 }
00643
00644 bool coarsen () const
00645 {
00646 return (coarsenMarked_ > 0);
00647 }
00648
00649 int refineMarked () const
00650 {
00651 return refineMarked_;
00652 }
00653
00654 void preAdapt ()
00655 {
00656 if( phase_ != ComputationPhase )
00657 error( "preAdapt may only be called in computation phase." );
00658 phase_ = PreAdaptationPhase;
00659 }
00660
00661 void adapt ()
00662 {
00663 if( phase_ != PreAdaptationPhase )
00664 error( "adapt may only be called in preadapdation phase." );
00665 phase_ = PostAdaptationPhase;
00666 }
00667
00668 void postAdapt ()
00669 {
00670 if( phase_ != PostAdaptationPhase )
00671 error( "postAdapt may only be called in postadaptation phase." );
00672 phase_ = ComputationPhase;
00673
00674 coarsenMarked_ = 0;
00675 refineMarked_ = 0;
00676 }
00677
00678 private:
00679 void error ( const std::string &message )
00680 {
00681 DUNE_THROW( InvalidStateException, message );
00682 }
00683 };
00684
00685 }
00686
00687 #include "agmemory.hh"
00688 #include "albertagrid.cc"
00689
00690
00691 #undef DIM
00692 #undef DIM_OF_WORLD
00693 #undef CALC_COORD
00694
00695 #ifdef _ABS_NOT_DEFINED_
00696 #undef ABS
00697 #endif
00698
00699 #ifdef _MIN_NOT_DEFINED_
00700 #undef MIN
00701 #endif
00702
00703 #ifdef _MAX_NOT_DEFINED_
00704 #undef MAX
00705 #endif
00706
00707 #if DUNE_ALBERTA_VERSION >= 0x201
00708 #ifdef obstack_chunk_alloc
00709 #undef obstack_chunk_alloc
00710 #endif
00711 #ifdef obstack_chunk_free
00712 #undef obstack_chunk_free
00713 #endif
00714 #include <dune/grid/albertagrid/undefine-2.1.hh>
00715 #elif DUNE_ALBERTA_VERSION == 0x200
00716 #include <dune/grid/albertagrid/undefine-2.0.hh>
00717 #else
00718 #include <dune/grid/albertagrid/undefine-1.2.hh>
00719 #endif
00720
00721 #define _ALBERTA_H_
00722
00723 #endif // HAVE_ALBERTA
00724
00725 #endif