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< AlbertaGrid<dim,dimworld> > 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
00382 template< class DataHandle >
00383 bool adapt ( AdaptDataHandleInterface< This, DataHandle > &handle );
00384
00385 template< class DofManager, class RestrictProlongOperator >
00386 bool DUNE_DEPRECATED
00387 adapt ( DofManager &, RestrictProlongOperator &, bool verbose = false );
00388
00390 bool preAdapt ();
00391
00393 void postAdapt();
00394
00397 const CollectiveCommunication &comm () const
00398 {
00399 return comm_;
00400 }
00401
00402 static std::string typeName ()
00403 {
00404 std::ostringstream s;
00405 s << "AlbertaGrid< " << dim << ", " << dimworld << " >";
00406 return s.str();
00407 }
00408
00410 std::string name () const
00411 {
00412 return mesh_.name();
00413 };
00414
00415
00416
00417
00419 template< GrapeIOFileFormatType ftype >
00420 bool writeGrid( const std::string &filename, ctype time ) const;
00421
00423 template< GrapeIOFileFormatType ftype >
00424 bool readGrid( const std::string &filename, ctype &time );
00425
00426
00427 const HierarchicIndexSet & hierarchicIndexSet () const { return hIndexSet_; }
00428
00430 const typename Traits :: LevelIndexSet & levelIndexSet (int level) const;
00431
00433 const typename Traits :: LeafIndexSet & leafIndexSet () const;
00434
00436 const GlobalIdSet &globalIdSet () const
00437 {
00438 return idSet_;
00439 }
00440
00442 const LocalIdSet &localIdSet () const
00443 {
00444 return idSet_;
00445 }
00446
00447
00448 ALBERTA MESH* getMesh () const
00449 {
00450 return mesh_;
00451 };
00452
00453 const MeshPointer &meshPointer () const
00454 {
00455 return mesh_;
00456 }
00457
00458 const DofNumbering &dofNumbering () const
00459 {
00460 return dofNumbering_;
00461 }
00462
00463 const LevelProvider &levelProvider () const
00464 {
00465 return levelProvider_;
00466 }
00467
00468 int dune2alberta ( int codim, int i ) const
00469 {
00470 return numberingMap_.dune2alberta( codim, i );
00471 }
00472
00473 int alberta2dune ( int codim, int i ) const
00474 {
00475 return numberingMap_.alberta2dune( codim, i );
00476 }
00477
00478 private:
00479 using Base::getRealImplementation;
00480
00481 typedef std::vector<int> ArrayType;
00482
00483 void setup ();
00484
00485
00486
00487 void calcExtras();
00488
00489
00490 bool writeGridXdr ( const std::string &filename, ctype time ) const;
00491
00493 bool readGridXdr ( const std::string &filename, ctype &time );
00494
00495 #if 0
00497 bool readGridAscii ( const std::string &filename, ctype &time );
00498 #endif
00499
00500
00501 void removeMesh();
00502
00503
00504 MeshPointer mesh_;
00505
00506
00507 CollectiveCommunication comm_;
00508
00509
00510 int maxlevel_;
00511
00512
00513
00514
00515 typedef MakeableInterfaceObject< typename Traits::template Codim< 0 >::Entity >
00516 EntityObject;
00517
00518 public:
00519 typedef AGMemoryProvider< EntityObject > EntityProvider;
00520
00521 typedef AlbertaGridIntersectionIterator< const This > IntersectionIteratorImp;
00522 typedef IntersectionIteratorImp LeafIntersectionIteratorImp;
00523 typedef AGMemoryProvider< LeafIntersectionIteratorImp > LeafIntersectionIteratorProviderType;
00524 friend class LeafIntersectionIteratorWrapper< const This >;
00525
00526 typedef LeafIntersectionIteratorWrapper< const This >
00527 AlbertaGridIntersectionIteratorType;
00528
00529 LeafIntersectionIteratorProviderType & leafIntersetionIteratorProvider() const { return leafInterItProvider_; }
00530
00531 private:
00532 mutable EntityProvider entityProvider_;
00533 mutable LeafIntersectionIteratorProviderType leafInterItProvider_;
00534
00535 public:
00536 template< class IntersectionInterfaceType >
00537 const typename Base
00538 :: template ReturnImplementationType< IntersectionInterfaceType >
00539 :: ImplementationType & DUNE_DEPRECATED
00540 getRealIntersectionIterator ( const IntersectionInterfaceType &iterator ) const
00541 {
00542 return this->getRealImplementation( iterator );
00543 }
00544
00545 template< class IntersectionType >
00546 const typename Base
00547 :: template ReturnImplementationType< IntersectionType >
00548 :: ImplementationType &
00549 getRealIntersection ( const IntersectionType &intersection ) const
00550 {
00551 return this->getRealImplementation( intersection );
00552 }
00553
00554
00555 template< int codim >
00556 MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *
00557 getNewEntity () const;
00558
00559
00560 template <int codim>
00561 void freeEntity ( MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *entity ) const;
00562
00563 public:
00564
00565 const Alberta::GlobalVector &
00566 getCoord ( const ElementInfo &elementInfo, int vertex ) const;
00567
00568 private:
00569
00570 Alberta::NumberingMap< dimension > numberingMap_;
00571
00572 DofNumbering dofNumbering_;
00573
00574 LevelProvider levelProvider_;
00575
00576
00577 HierarchicIndexSet hIndexSet_;
00578
00579
00580 IdSetImp idSet_;
00581
00582
00583
00584 mutable std::vector< LevelIndexSetImp * > levelIndexVec_;
00585
00586
00587
00588 mutable LeafIndexSetImp* leafIndexSet_;
00589
00590 typedef SingleTypeSizeCache< This > SizeCacheType;
00591 SizeCacheType * sizeCache_;
00592
00593 typedef AlbertaMarkerVector< dim, dimworld > MarkerVector;
00594
00595
00596 mutable MarkerVector leafMarkerVector_;
00597
00598
00599 mutable std::vector< MarkerVector > levelMarkerVector_;
00600
00601 #if !CALC_COORD
00602 Alberta::CoordCache< dimension > coordCache_;
00603 #endif
00604
00605
00606 AdaptationState adaptationState_;
00607 };
00608
00609
00610
00611
00612
00613
00614 template< int dim, int dimworld >
00615 struct AlbertaGrid< dim, dimworld >::AdaptationState
00616 {
00617 enum Phase { ComputationPhase, PreAdaptationPhase, PostAdaptationPhase };
00618
00619 private:
00620 Phase phase_;
00621 int coarsenMarked_;
00622 int refineMarked_;
00623
00624 public:
00625 AdaptationState ()
00626 : phase_( ComputationPhase ),
00627 coarsenMarked_( 0 ),
00628 refineMarked_( 0 )
00629 {}
00630
00631 void mark ( int count )
00632 {
00633 if( count < 0 )
00634 ++coarsenMarked_;
00635 if( count > 0 )
00636 refineMarked_ += (2 << count);
00637 }
00638
00639 void unmark ( int count )
00640 {
00641 if( count < 0 )
00642 --coarsenMarked_;
00643 if( count > 0 )
00644 refineMarked_ -= (2 << count);
00645 }
00646
00647 bool coarsen () const
00648 {
00649 return (coarsenMarked_ > 0);
00650 }
00651
00652 int refineMarked () const
00653 {
00654 return refineMarked_;
00655 }
00656
00657 void preAdapt ()
00658 {
00659 if( phase_ != ComputationPhase )
00660 error( "preAdapt may only be called in computation phase." );
00661 phase_ = PreAdaptationPhase;
00662 }
00663
00664 void adapt ()
00665 {
00666 if( phase_ != PreAdaptationPhase )
00667 error( "adapt may only be called in preadapdation phase." );
00668 phase_ = PostAdaptationPhase;
00669 }
00670
00671 void postAdapt ()
00672 {
00673 if( phase_ != PostAdaptationPhase )
00674 error( "postAdapt may only be called in postadaptation phase." );
00675 phase_ = ComputationPhase;
00676
00677 coarsenMarked_ = 0;
00678 refineMarked_ = 0;
00679 }
00680
00681 private:
00682 void error ( const std::string &message )
00683 {
00684 DUNE_THROW( InvalidStateException, message );
00685 }
00686 };
00687
00688 }
00689
00690 #include "agmemory.hh"
00691 #include "albertagrid.cc"
00692
00693
00694 #undef DIM
00695 #undef DIM_OF_WORLD
00696 #undef CALC_COORD
00697
00698 #ifdef _ABS_NOT_DEFINED_
00699 #undef ABS
00700 #endif
00701
00702 #ifdef _MIN_NOT_DEFINED_
00703 #undef MIN
00704 #endif
00705
00706 #ifdef _MAX_NOT_DEFINED_
00707 #undef MAX
00708 #endif
00709
00710 #if DUNE_ALBERTA_VERSION >= 0x201
00711 #ifdef obstack_chunk_alloc
00712 #undef obstack_chunk_alloc
00713 #endif
00714 #ifdef obstack_chunk_free
00715 #undef obstack_chunk_free
00716 #endif
00717 #include <dune/grid/albertagrid/undefine-2.1.hh>
00718 #elif DUNE_ALBERTA_VERSION == 0x200
00719 #include <dune/grid/albertagrid/undefine-2.0.hh>
00720 #else
00721 #include <dune/grid/albertagrid/undefine-1.2.hh>
00722 #endif
00723
00724 #define _ALBERTA_H_
00725
00726 #endif // HAVE_ALBERTA
00727
00728 #endif