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 template< int, class > friend class AlbertaGridEntityPointer;
00208
00209 friend class AlbertaGridHierarchicIterator<AlbertaGrid<dim,dimworld> >;
00210
00211 friend class AlbertaGridIntersectionIterator<AlbertaGrid<dim,dimworld> >;
00212 friend class AlbertaGridIntersectionIterator<const AlbertaGrid<dim,dimworld> >;
00213
00214 friend class AlbertaMarkerVector< dim, dimworld >;
00215 friend class AlbertaGridHierarchicIndexSet<dim,dimworld>;
00216
00217 public:
00218 typedef Alberta::Real ctype;
00219
00220 static const int dimension = dim;
00221 static const int dimensionworld = dimworld;
00222
00224 typedef AlbertaGridFamily< dim, dimworld > GridFamily;
00225
00226
00227 typedef typename AlbertaGridFamily< dim, dimworld >::Traits Traits;
00228
00230 typedef typename Traits::HierarchicIndexSet HierarchicIndexSet;
00231
00233 typedef typename Traits::CollectiveCommunication CollectiveCommunication;
00234
00235 private:
00237 typedef typename Traits::template Codim<0>::LeafIterator LeafIterator;
00238
00240 typedef typename GridFamily:: LevelIndexSetImp LevelIndexSetImp;
00241 typedef typename GridFamily:: LeafIndexSetImp LeafIndexSetImp;
00242
00244 typedef typename Traits :: LeafIndexSet LeafIndexSet;
00245
00247 typedef AlbertaGridIdSet<dim,dimworld> IdSetImp;
00248 typedef typename Traits :: GlobalIdSet GlobalIdSet;
00249 typedef typename Traits :: LocalIdSet LocalIdSet;
00250
00251 struct AdaptationState;
00252
00253 template< class DataHandler >
00254 struct AdaptationCallback;
00255
00256 public:
00258 typedef typename ALBERTA AlbertHelp::AlbertLeafData<dimworld,dim+1> LeafDataType;
00259
00260 private:
00261
00262 static const int MAXL = 64;
00263
00264 typedef Alberta::ElementInfo< dimension > ElementInfo;
00265 typedef Alberta::MeshPointer< dimension > MeshPointer;
00266 typedef Alberta::HierarchyDofNumbering< dimension > DofNumbering;
00267 typedef AlbertaGridLevelProvider< dimension > LevelProvider;
00268
00269
00270 AlbertaGrid ( const This & );
00271 This &operator= ( const This & );
00272
00273 public:
00275 AlbertaGrid ();
00276
00282 AlbertaGrid ( const Alberta::MacroData< dimension > ¯oData,
00283 const std::string &gridName = "AlbertaGrid" );
00284
00290 AlbertaGrid ( const std::string ¯oGridFileName,
00291 const std::string &gridName = "AlbertaGrid" );
00292
00294 ~AlbertaGrid ();
00295
00298 int maxLevel () const;
00299
00301 template<int cd, PartitionIteratorType pitype>
00302 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
00303 lbegin (int level) const;
00304
00306 template<int cd, PartitionIteratorType pitype>
00307 typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
00308 lend (int level) const;
00309
00311 template< int codim >
00312 typename Traits::template Codim< codim >::LevelIterator
00313 lbegin ( int level ) const;
00314
00316 template< int codim >
00317 typename Traits::template Codim< codim >::LevelIterator
00318 lend ( int level ) const;
00319
00321 template< int codim, PartitionIteratorType pitype >
00322 typename Traits
00323 ::template Codim< codim >::template Partition< pitype >::LeafIterator
00324 leafbegin () const;
00325
00327 template< int codim, PartitionIteratorType pitype >
00328 typename Traits
00329 ::template Codim< codim >::template Partition< pitype >::LeafIterator
00330 leafend () const;
00331
00333 template< int codim >
00334 typename Traits::template Codim< codim >::LeafIterator
00335 leafbegin () const;
00336
00338 template< int codim >
00339 typename Traits::template Codim< codim >::LeafIterator
00340 leafend () const;
00341
00346 int size (int level, int codim) const;
00347
00349 int size (int level, GeometryType type) const;
00350
00352 int size (int codim) const;
00353
00355 int size (GeometryType type) const;
00356
00357 public:
00358
00359
00360
00361 using Base::getMark;
00362 using Base::mark;
00363
00365 int getMark ( const typename Traits::template Codim< 0 >::Entity &e ) const;
00366
00368 bool mark ( int refCount, const typename Traits::template Codim< 0 >::Entity &e );
00369
00371 void globalRefine ( int refCount );
00372
00373 template< class DataHandle >
00374 void globalRefine ( int refCount, AdaptDataHandleInterface< This, DataHandle > &handle );
00375
00377 bool adapt ();
00378
00383 template< class DataHandle >
00384 bool adapt ( AdaptDataHandleInterface< This, DataHandle > &handle );
00385
00386 template< class DofManager, class RestrictProlongOperator >
00387 bool DUNE_DEPRECATED
00388 adapt ( DofManager &, RestrictProlongOperator &, bool verbose = false );
00389
00391 bool preAdapt ();
00392
00394 void postAdapt();
00395
00398 const CollectiveCommunication &comm () const
00399 {
00400 return comm_;
00401 }
00402
00403 static std::string typeName ()
00404 {
00405 std::ostringstream s;
00406 s << "AlbertaGrid< " << dim << ", " << dimworld << " >";
00407 return s.str();
00408 }
00409
00411 std::string name () const
00412 {
00413 return mesh_.name();
00414 };
00415
00416
00417
00418
00420 template< GrapeIOFileFormatType ftype >
00421 bool writeGrid( const std::string &filename, ctype time ) const;
00422
00424 template< GrapeIOFileFormatType ftype >
00425 bool readGrid( const std::string &filename, ctype &time );
00426
00427
00428 const HierarchicIndexSet & hierarchicIndexSet () const { return hIndexSet_; }
00429
00431 const typename Traits :: LevelIndexSet & levelIndexSet (int level) const;
00432
00434 const typename Traits :: LeafIndexSet & leafIndexSet () const;
00435
00437 const GlobalIdSet &globalIdSet () const
00438 {
00439 return idSet_;
00440 }
00441
00443 const LocalIdSet &localIdSet () const
00444 {
00445 return idSet_;
00446 }
00447
00448
00449 ALBERTA MESH* getMesh () const
00450 {
00451 return mesh_;
00452 };
00453
00454 const MeshPointer &meshPointer () const
00455 {
00456 return mesh_;
00457 }
00458
00459 const DofNumbering &dofNumbering () const
00460 {
00461 return dofNumbering_;
00462 }
00463
00464 const LevelProvider &levelProvider () const
00465 {
00466 return levelProvider_;
00467 }
00468
00469 int dune2alberta ( int codim, int i ) const
00470 {
00471 return numberingMap_.dune2alberta( codim, i );
00472 }
00473
00474 int alberta2dune ( int codim, int i ) const
00475 {
00476 return numberingMap_.alberta2dune( codim, i );
00477 }
00478
00479 private:
00480 using Base::getRealImplementation;
00481
00482 typedef std::vector<int> ArrayType;
00483
00484 void setup ();
00485
00486
00487
00488 void calcExtras();
00489
00490
00491 bool writeGridXdr ( const std::string &filename, ctype time ) const;
00492
00494 bool readGridXdr ( const std::string &filename, ctype &time );
00495
00496 #if 0
00498 bool readGridAscii ( const std::string &filename, ctype &time );
00499 #endif
00500
00501
00502 void removeMesh();
00503
00504
00505 MeshPointer mesh_;
00506
00507
00508 CollectiveCommunication comm_;
00509
00510
00511 int maxlevel_;
00512
00513
00514
00515
00516 typedef MakeableInterfaceObject< typename Traits::template Codim< 0 >::Entity >
00517 EntityObject;
00518
00519 public:
00520 typedef AGMemoryProvider< EntityObject > EntityProvider;
00521
00522 typedef AlbertaGridIntersectionIterator< const This > IntersectionIteratorImp;
00523 typedef IntersectionIteratorImp LeafIntersectionIteratorImp;
00524 typedef AGMemoryProvider< LeafIntersectionIteratorImp > LeafIntersectionIteratorProviderType;
00525 friend class LeafIntersectionIteratorWrapper< const This >;
00526
00527 typedef LeafIntersectionIteratorWrapper< const This >
00528 AlbertaGridIntersectionIteratorType;
00529
00530 LeafIntersectionIteratorProviderType & leafIntersetionIteratorProvider() const { return leafInterItProvider_; }
00531
00532 private:
00533 mutable EntityProvider entityProvider_;
00534 mutable LeafIntersectionIteratorProviderType leafInterItProvider_;
00535
00536 public:
00537 template< class IntersectionInterfaceType >
00538 const typename Base
00539 :: template ReturnImplementationType< IntersectionInterfaceType >
00540 :: ImplementationType & DUNE_DEPRECATED
00541 getRealIntersectionIterator ( const IntersectionInterfaceType &iterator ) const
00542 {
00543 return this->getRealImplementation( iterator );
00544 }
00545
00546 template< class IntersectionType >
00547 const typename Base
00548 :: template ReturnImplementationType< IntersectionType >
00549 :: ImplementationType &
00550 getRealIntersection ( const IntersectionType &intersection ) const
00551 {
00552 return this->getRealImplementation( intersection );
00553 }
00554
00555
00556 template< int codim >
00557 MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *
00558 getNewEntity () const;
00559
00560
00561 template <int codim>
00562 void freeEntity ( MakeableInterfaceObject< typename Traits::template Codim< codim >::Entity > *entity ) const;
00563
00564 public:
00565
00566 const Alberta::GlobalVector &
00567 getCoord ( const ElementInfo &elementInfo, int vertex ) const;
00568
00569 private:
00570
00571 Alberta::NumberingMap< dimension > numberingMap_;
00572
00573 DofNumbering dofNumbering_;
00574
00575 LevelProvider levelProvider_;
00576
00577
00578 HierarchicIndexSet hIndexSet_;
00579
00580
00581 IdSetImp idSet_;
00582
00583
00584
00585 mutable std::vector< LevelIndexSetImp * > levelIndexVec_;
00586
00587
00588
00589 mutable LeafIndexSetImp* leafIndexSet_;
00590
00591 typedef SingleTypeSizeCache< This > SizeCacheType;
00592 SizeCacheType * sizeCache_;
00593
00594 typedef AlbertaMarkerVector< dim, dimworld > MarkerVector;
00595
00596
00597 mutable MarkerVector leafMarkerVector_;
00598
00599
00600 mutable std::vector< MarkerVector > levelMarkerVector_;
00601
00602 #if !CALC_COORD
00603 Alberta::CoordCache< dimension > coordCache_;
00604 #endif
00605
00606
00607 AdaptationState adaptationState_;
00608 };
00609
00610
00611
00612
00613
00614
00615 template< int dim, int dimworld >
00616 struct AlbertaGrid< dim, dimworld >::AdaptationState
00617 {
00618 enum Phase { ComputationPhase, PreAdaptationPhase, PostAdaptationPhase };
00619
00620 private:
00621 Phase phase_;
00622 int coarsenMarked_;
00623 int refineMarked_;
00624
00625 public:
00626 AdaptationState ()
00627 : phase_( ComputationPhase ),
00628 coarsenMarked_( 0 ),
00629 refineMarked_( 0 )
00630 {}
00631
00632 void mark ( int count )
00633 {
00634 if( count < 0 )
00635 ++coarsenMarked_;
00636 if( count > 0 )
00637 refineMarked_ += (2 << count);
00638 }
00639
00640 void unmark ( int count )
00641 {
00642 if( count < 0 )
00643 --coarsenMarked_;
00644 if( count > 0 )
00645 refineMarked_ -= (2 << count);
00646 }
00647
00648 bool coarsen () const
00649 {
00650 return (coarsenMarked_ > 0);
00651 }
00652
00653 int refineMarked () const
00654 {
00655 return refineMarked_;
00656 }
00657
00658 void preAdapt ()
00659 {
00660 if( phase_ != ComputationPhase )
00661 error( "preAdapt may only be called in computation phase." );
00662 phase_ = PreAdaptationPhase;
00663 }
00664
00665 void adapt ()
00666 {
00667 if( phase_ != PreAdaptationPhase )
00668 error( "adapt may only be called in preadapdation phase." );
00669 phase_ = PostAdaptationPhase;
00670 }
00671
00672 void postAdapt ()
00673 {
00674 if( phase_ != PostAdaptationPhase )
00675 error( "postAdapt may only be called in postadaptation phase." );
00676 phase_ = ComputationPhase;
00677
00678 coarsenMarked_ = 0;
00679 refineMarked_ = 0;
00680 }
00681
00682 private:
00683 void error ( const std::string &message )
00684 {
00685 DUNE_THROW( InvalidStateException, message );
00686 }
00687 };
00688
00689 }
00690
00691 #include "agmemory.hh"
00692 #include "albertagrid.cc"
00693
00694
00695 #undef DIM
00696 #undef DIM_OF_WORLD
00697 #undef CALC_COORD
00698
00699 #ifdef _ABS_NOT_DEFINED_
00700 #undef ABS
00701 #endif
00702
00703 #ifdef _MIN_NOT_DEFINED_
00704 #undef MIN
00705 #endif
00706
00707 #ifdef _MAX_NOT_DEFINED_
00708 #undef MAX
00709 #endif
00710
00711 #if DUNE_ALBERTA_VERSION >= 0x201
00712 #ifdef obstack_chunk_alloc
00713 #undef obstack_chunk_alloc
00714 #endif
00715 #ifdef obstack_chunk_free
00716 #undef obstack_chunk_free
00717 #endif
00718 #include <dune/grid/albertagrid/undefine-2.1.hh>
00719 #elif DUNE_ALBERTA_VERSION == 0x200
00720 #include <dune/grid/albertagrid/undefine-2.0.hh>
00721 #else
00722 #include <dune/grid/albertagrid/undefine-1.2.hh>
00723 #endif
00724
00725 #define _ALBERTA_H_
00726
00727 #endif // HAVE_ALBERTA
00728
00729 #endif