Skip to content
Snippets Groups Projects
Commit 94a50424 authored by Kilian Weishaupt's avatar Kilian Weishaupt Committed by Timo Koch
Browse files

[io][grid] Remove vtk output from subgrid manager

parent 9e5a44f8
No related branches found
No related tags found
1 merge request!1594[io][grid] Remove vtk output from subgrid manager
......@@ -29,12 +29,15 @@
#include <memory>
#include <dune/subgrid/subgrid.hh>
#include <dune/grid/io/file/vtk.hh>
#include <dune/grid/io/file/dgfparser/dgfwriter.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/boundaryflag.hh>
// TODO: remove this after 3.1 is released
#include <dune/grid/io/file/vtk.hh>
// TODO: remove this after 3.1 is released
#include <dune/grid/io/file/dgfparser/dgfwriter.hh>
namespace Dumux {
/*!
......@@ -71,24 +74,32 @@ public:
// choose which elements to add to the subgrid.
auto hostGridView = subgridPtr->getHostGrid().leafGridView();
for (const auto& e : elements(hostGridView))
if(selector(e))
if (selector(e))
elementsForSubgrid.insert(globalIDset.template id<0>(e));
subgridPtr->insertSetPartial(elementsForSubgrid);
subgridPtr->createEnd();
// TODO: remove this after 3.1 is released
// If desired, write out the final subgrid as a dgf file.
if(getParamFromGroup<bool>(modelParamGroup, "Grid.WriteSubGridToDGF", false))
if (getParamFromGroup<bool>(modelParamGroup, "Grid.WriteSubGridToDGF", false))
{
std::cerr << "Deprecation warning: SubGridManager: Grid.WriteSubGridToDGF is deprecated."
<< "Use Dune::VTKWriter to write out your grid manually." << std::endl;
const auto postfix = getParamFromGroup<std::string>(modelParamGroup, "Problem.Name", "");
const std::string name = postfix == "" ? "subgrid" : "subgrid_" + postfix;
Dune::DGFWriter<typename Grid::LeafGridView> writer(subgridPtr->leafGridView());
writer.write(name + ".dgf");
}
// TODO: remove this after 3.1 is released
// If desired, write out the hostgrid as vtk file.
if(getParamFromGroup<bool>(modelParamGroup, "Grid.WriteSubGridToVtk", false))
if (getParamFromGroup<bool>(modelParamGroup, "Grid.WriteSubGridToVtk", false))
{
std::cerr << "Deprecation warning: SubGridManager: Grid.WriteSubGridToVtk is deprecated."
<< "Use Dune::VTKWriter to write out your grid manually." << std::endl;
const auto postfix = getParamFromGroup<std::string>(modelParamGroup, "Problem.Name", "");
const std::string name = postfix == "" ? "subgrid" : "subgrid_" + postfix;
Dune::VTKWriter<typename Grid::LeafGridView> vtkWriter(subgridPtr->leafGridView());
......
......@@ -110,15 +110,29 @@ int main(int argc, char** argv) try
CircleSelector<GlobalPosition> elementSelectorThree(center);
// Create three different subgrids from the same hostgrid.
auto subgridPtrOne = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorOne, "SubGridOne");
auto subgridPtrTwo = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorTwo, "SubGridTwo");
auto subgridPtrThree = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorThree, "SubGridThree");
auto subgridPtrOne = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorOne);
auto subgridPtrTwo = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorTwo);
auto subgridPtrThree = SubgridManager<HostGrid>::makeGrid(hostGrid, elementSelectorThree);
std::cout << "Constructing a host grid and three subgrids took " << timer.elapsed() << " seconds.\n";
// Write out the host grid and the subgrids.
Dune::VTKWriter<HostGrid::LeafGridView> vtkWriter(hostGrid.leafGridView());
vtkWriter.write("hostgrid");
{
Dune::VTKWriter<HostGrid::LeafGridView> vtkWriter(hostGrid.leafGridView());
vtkWriter.write("hostgrid");
}
{
Dune::VTKWriter<SubgridManager<HostGrid>::Grid::LeafGridView> vtkWriter(subgridPtrOne->leafGridView());
vtkWriter.write("subgrid_one");
}
{
Dune::VTKWriter<SubgridManager<HostGrid>::Grid::LeafGridView> vtkWriter(subgridPtrTwo->leafGridView());
vtkWriter.write("subgrid_two");
}
{
Dune::VTKWriter<SubgridManager<HostGrid>::Grid::LeafGridView> vtkWriter(subgridPtrThree->leafGridView());
vtkWriter.write("subgrid_three");
}
return 0;
}
......
......@@ -8,5 +8,3 @@ Positions1 = 0 1
Cells0 = 1
Cells1 = 1
Refinement = 4
WriteSubGridToVtk = true
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