Skip to content
Snippets Groups Projects
Commit ad54d973 authored by Oliver Sander's avatar Oliver Sander
Browse files

fixed the parallel lbegin/lend methods

[[Imported from SVN: r1662]]
parent e61521f9
No related branches found
No related tags found
No related merge requests found
......@@ -86,14 +86,12 @@ namespace Dune
Please send any questions, suggestions, or bug reports to
sander@math.fu-berlin.de
@{
*/
/** \brief The type used by UG to store coordinates */
typedef double UGCtype;
/** @} end documentation group */
// forward declarations
template<int dim, int dimworld> class UGGrid;
......@@ -104,10 +102,6 @@ namespace Dune
template<class GridImp> class UGGridHierarchicIterator;
template<class GridImp> class UGGridIntersectionIterator;
// singleton holding reference elements
//template<int dim> struct UGGridReferenceElement;
} // namespace Dune
#include "uggrid/uggridgeometry.hh"
......@@ -128,7 +122,15 @@ namespace Dune {
/** \brief The UG %Grid class
* \ingroup UGGrid
*
* \todo Please doc me!
This is the implementation of the grid interface
using the UG grid management system. UG provides conforming grids
in two and three space dimensions. The grids can be mixed, i.e.
2d grids can contain triangles and quadrilaterals and 3d grid can
contain tetrahedra and hexahedra and also pyramids and prisms.
The grid refinement rules are very flexible. Local adaptive red/green
refinement is the default, but a special method in the UGGrid class
allows you to directly access a number of anisotropic refinements rules.
Last but not least, the UG grid manager is completely parallelized.
*/
template <int dim, int dimworld>
class UGGrid : public GridDefault <dim, dimworld,UGCtype, UGGrid<dim,dimworld> >
......@@ -200,11 +202,11 @@ namespace Dune {
//! Iterator to first entity of given codim on level
template<int codim, PartitionIteratorType PiType>
UGGridLevelIterator<codim,PiType,UGGrid<dim,dimworld> > lbegin (int level) const;
typename Traits::template codim<codim>::template partition<PiType>::LevelIterator lbegin (int level) const;
//! one past the end on this level
template<int codim, PartitionIteratorType PiType>
UGGridLevelIterator<codim,PiType,UGGrid<dim,dimworld> > lend (int level) const;
typename Traits::template codim<codim>::template partition<PiType>::LevelIterator lend (int level) const;
/** \brief Number of grid entities per level and codim
*/
......@@ -287,9 +289,11 @@ namespace Dune {
void makeNewUGMultigrid();
/** \brief Does one uniform refinement step
/** \brief Does uniform refinement
*
* \param n Number of uniform refinement steps
*/
void globalRefine();
void globalRefine(int n);
// UG multigrid, which contains the data
typename UGTypes<dimworld>::MultiGridType* multigrid_;
......
......@@ -307,7 +307,7 @@ UGGrid<dim, dimworld>::lbegin (int level) const
template<int dim, int dimworld>
template<int codim, PartitionIteratorType PiType>
inline UGGridLevelIterator<codim,PiType,UGGrid<dim,dimworld> >
inline typename UGGrid<dim,dimworld>::Traits::template codim<codim>::template partition<PiType>::LevelIterator
UGGrid<dim, dimworld>::lbegin (int level) const
{
assert(multigrid_);
......@@ -316,7 +316,7 @@ UGGrid<dim, dimworld>::lbegin (int level) const
if (!theGrid)
DUNE_THROW(GridError, "LevelIterator in nonexisting level " << level << " requested!");
return UGGridLevelIteratorFactory<codim,PiType,UGGrid<dim,dimworld> >::getIterator(theGrid, level);
return UGGridLevelIteratorFactory<codim,PiType, const UGGrid<dim,dimworld> >::getIterator(theGrid, level);
}
template < int dim, int dimworld >
......@@ -324,17 +324,15 @@ template<int codim>
typename UGGrid<dim,dimworld>::Traits::template codim<codim>::LevelIterator
UGGrid < dim, dimworld >::lend (int level) const
{
UGGridLevelIterator<codim,All_Partition, const UGGrid<dim,dimworld> > it(level);
return it;
return UGGridLevelIterator<codim,All_Partition, const UGGrid<dim,dimworld> >(level);
}
template < int dim, int dimworld >
template<int codim, PartitionIteratorType PiType>
inline UGGridLevelIterator<codim,PiType,UGGrid<dim,dimworld> >
inline typename UGGrid<dim,dimworld>::Traits::template codim<codim>::template partition<PiType>::LevelIterator
UGGrid < dim, dimworld >::lend (int level) const
{
UGGridLevelIterator<codim,PiType,UGGrid<dim,dimworld> > it(level);
return it;
return UGGridLevelIterator<codim,PiType, const UGGrid<dim,dimworld> >(level);
}
template < int dim, int dimworld >
......@@ -543,18 +541,22 @@ void UGGrid <dim, dimworld>::adaptWithoutClosure()
}
template < int dim, int dimworld >
void UGGrid < dim, dimworld >::globalRefine()
void UGGrid < dim, dimworld >::globalRefine(int n)
{
// mark all entities for grid refinement
typename Traits::template codim<0>::LevelIterator iIt = lbegin<0>(maxlevel());
typename Traits::template codim<0>::LevelIterator iEndIt = lend<0>(maxlevel());
for (int i=0; i<n; i++) {
for (; iIt!=iEndIt; ++iIt)
mark(1, iIt);
// mark all entities for grid refinement
typename Traits::template codim<0>::LevelIterator iIt = lbegin<0>(maxlevel());
typename Traits::template codim<0>::LevelIterator iEndIt = lend<0>(maxlevel());
this->preAdapt();
adapt();
this->postAdapt();
for (; iIt!=iEndIt; ++iIt)
mark(1, iIt);
this->preAdapt();
adapt();
this->postAdapt();
}
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment