Skip to content
Snippets Groups Projects
Commit 3c8bca7e authored by graeser's avatar graeser
Browse files

Add new example demonstrating grid views

parent 4b649ec2
No related branches found
No related tags found
No related merge requests found
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifdef HAVE_CONFIG_H
# include "config.h"
#endif
// included standard library headers
#include <iostream>
#include <array>
// included dune-common headers
#include <dune/common/parallel/mpihelper.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/fvector.hh>
#include <dune/common/stringutility.hh>
// included dune-geometry headers
#include <dune/geometry/quadraturerules.hh>
// included dune-grid headers
#include <dune/grid/io/file/vtk/vtkwriter.hh>
#include <dune/grid/utility/structuredgridfactory.hh>
#include <dune/grid/yaspgrid.hh>
#include <dune/grid/uggrid.hh>
// included dune-fu-tutorial headers
template<class Grid>
auto createCubeGrid()
{
static const int dim = Grid::dimension;
auto lowerLeft = Dune::FieldVector<double,dim>(0);
auto upperRight = Dune::FieldVector<double,dim>(1);
auto elementsPerDirection = std::array<unsigned int,dim>();
for(auto& e : elementsPerDirection)
e = 2;
return Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, elementsPerDirection);
}
template<class Grid>
auto createSimplexGrid()
{
static const int dim = Grid::dimension;
auto lowerLeft = Dune::FieldVector<double,dim>(0);
auto upperRight = Dune::FieldVector<double,dim>(1);
auto elementsPerDirection = std::array<unsigned int,dim>();
for(auto& e : elementsPerDirection)
e = 2;
return Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, elementsPerDirection);
}
template<class GridView>
void writeGridView(const GridView& gridView, std::string postFix)
{
Dune::VTKWriter<GridView> vtkWriter(gridView);
vtkWriter.write(std::string("04-gridviews-")+postFix);
}
template<class Grid>
void writeAllGridViews(const Grid& grid, std::string gridName)
{
for(int level = 0; level <= grid.maxLevel(); ++level)
writeGridView(grid.levelGridView(level), Dune::formatString(gridName+"-level-%02d", level));
writeGridView(grid.leafGridView(), gridName+"-leaf");
}
int main(int argc, char** argv)
{
try{
// Maybe initialize MPI
Dune::MPIHelper& helper = Dune::MPIHelper::instance(argc, argv);
// Print process rank
if(Dune::MPIHelper::isFake)
std::cout<< "This is a sequential program." << std::endl;
else
std::cout<<"I am rank "<<helper.rank()<<" of "<<helper.size()
<<" processes!"<<std::endl;
{
using Grid = Dune::YaspGrid<2>;
auto gridPointer = createCubeGrid<Grid>();
Grid& grid = *gridPointer;
grid.globalRefine(2);
writeAllGridViews(grid, "2d_yaspgrid");
}
{
using Grid = Dune::YaspGrid<3>;
auto gridPointer = createCubeGrid<Grid>();
Grid& grid = *gridPointer;
grid.globalRefine(2);
writeAllGridViews(grid, "3d_yaspgrid");
}
{
using Grid = Dune::UGGrid<2>;
auto gridPointer = createSimplexGrid<Grid>();
Grid& grid = *gridPointer;
grid.globalRefine(2);
for(int i=0; i<2; ++i)
{
auto gridView = grid.levelGridView(grid.maxLevel());
auto it = Dune::elements(gridView).begin();
grid.mark(1, *it);
grid.adapt();
}
writeAllGridViews(grid, "2d_uggrid");
}
return 0;
}
catch (Dune::Exception &e){
std::cerr << "Dune reported error: " << e << std::endl;
}
catch (...){
std::cerr << "Unknown exception thrown!" << std::endl;
}
}
......@@ -9,3 +9,6 @@ target_link_dune_default_libraries("02-basic-grid-interface")
add_executable("03-function-integration" 03-function-integration.cc)
target_link_dune_default_libraries("03-function-integration")
add_executable("04-gridviews" 04-gridviews.cc)
target_link_dune_default_libraries("04-gridviews")
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment