Skip to content
Snippets Groups Projects
Commit 2ab057f2 authored by Elias Pipping's avatar Elias Pipping Committed by Elias Pipping
Browse files

Take the gridview from the basis

parent 9b43f618
No related branches found
No related tags found
No related merge requests found
......@@ -506,9 +506,8 @@ int main(int argc, char *argv[]) {
E, nu, gridDisplacement);
p0Assembler.assembleFunctional(localStressAssembler, vonMisesStress);
writeVtk<P1Basis, P0Basis, VectorType, SingletonVectorType, GridView>(
p1Basis, u, alpha, p0Basis, vonMisesStress, leafView,
(boost::format("obs%d") % run).str());
writeVtk(p1Basis, u, alpha, p0Basis, vonMisesStress,
(boost::format("obs%d") % run).str());
}
}
if (parset.get<bool>("enableTimer"))
......
......@@ -5,12 +5,12 @@
#include "vtk.hh"
template <class VertexBasis, class CellBasis, class VectorType,
class SingletonVectorType, class GridView>
class SingletonVectorType>
void writeVtk(VertexBasis const &vertexBasis, VectorType const &displacement,
SingletonVectorType const &state, CellBasis const &cellBasis,
SingletonVectorType const &stress, GridView const &gridView,
std::string const &filename) {
Dune::VTKWriter<GridView> writer(gridView);
SingletonVectorType const &stress, std::string const &filename) {
Dune::VTKWriter<typename VertexBasis::GridView> writer(
vertexBasis.getGridView());
auto const displacement_ptr =
Dune::make_shared<VTKBasisGridFunction<VertexBasis, VectorType> const>(
......
......@@ -7,10 +7,9 @@
#include <dune/fufem/functions/vtkbasisgridfunction.hh>
template <class VertexBasis, class CellBasis, class VectorType,
class SingletonVectorType, class GridView>
class SingletonVectorType>
void writeVtk(VertexBasis const &vertexBasis, VectorType const &displacement,
SingletonVectorType const &state, CellBasis const &cellBasis,
SingletonVectorType const &stress, GridView const &gridView,
std::string const &filename);
SingletonVectorType const &stress, std::string const &filename);
#endif
......@@ -21,11 +21,7 @@ using GridView = GridType::LeafGridView;
using P1Basis = P1NodalBasis<GridView, double>;
using MyP0Basis = P0Basis<GridView, double>;
template void writeVtk<P1Basis, MyP0Basis, VectorType, SingletonVectorType,
GridView>(P1Basis const &vertexBasis,
VectorType const &displacement,
SingletonVectorType const &state,
MyP0Basis const &cellBasis,
SingletonVectorType const &stress,
GridView const &gridView,
std::string const &filename);
template void writeVtk<P1Basis, MyP0Basis, VectorType, SingletonVectorType>(
P1Basis const &vertexBasis, VectorType const &displacement,
SingletonVectorType const &state, MyP0Basis const &cellBasis,
SingletonVectorType const &stress, std::string const &filename);
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment