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

[Cleanup] Allow vtk writer to be const

parent 3428c9a1
No related branches found
No related tags found
No related merge requests found
......@@ -287,7 +287,7 @@ int main(int argc, char *argv[]) {
writer.writeInfo(alpha_initial, u_initial, v_initial);
MyVTKWriter<typename MyAssembler::VertexBasis,
typename MyAssembler::CellBasis>
typename MyAssembler::CellBasis> const
vtkWriter(myAssembler.cellBasis, myAssembler.vertexBasis, "obs");
// Set up TNNMG solver
......@@ -442,7 +442,7 @@ int main(int argc, char *argv[]) {
ScalarVector stress;
myAssembler.assembleVonMisesStress(body.getYoungModulus(),
body.getPoissonRatio(), u, stress);
vtkWriter.write(u, v, alpha, stress);
vtkWriter.write(timeStep - 1, u, v, alpha, stress);
}
}
iterationWriter.close();
......
......@@ -12,17 +12,13 @@ template <class VertexBasis, class CellBasis>
MyVTKWriter<VertexBasis, CellBasis>::MyVTKWriter(
CellBasis const &_cellBasis, VertexBasis const &_vertexBasis,
std::string _prefix)
: cellBasis(_cellBasis),
vertexBasis(_vertexBasis),
prefix(_prefix),
counter(0) {}
: cellBasis(_cellBasis), vertexBasis(_vertexBasis), prefix(_prefix) {}
template <class VertexBasis, class CellBasis>
template <class Vector, class ScalarVector>
void MyVTKWriter<VertexBasis, CellBasis>::write(Vector const &u,
Vector const &v,
ScalarVector const &alpha,
ScalarVector const &stress) {
void MyVTKWriter<VertexBasis, CellBasis>::write(
size_t record, Vector const &u, Vector const &v, ScalarVector const &alpha,
ScalarVector const &stress) const {
Dune::VTKWriter<typename VertexBasis::GridView> writer(
vertexBasis.getGridView());
......@@ -46,7 +42,7 @@ void MyVTKWriter<VertexBasis, CellBasis>::write(Vector const &u,
cellBasis, stress, "stress");
writer.addCellData(stressPointer);
std::string const filename = prefix + std::to_string(counter++);
std::string const filename = prefix + std::to_string(record);
writer.write(filename.c_str());
}
......
......@@ -8,15 +8,13 @@ template <class VertexBasis, class CellBasis> class MyVTKWriter {
VertexBasis const &vertexBasis;
std::string const prefix;
size_t counter;
public:
MyVTKWriter(CellBasis const &cellBasis, VertexBasis const &vertexBasis,
std::string prefix);
template <class Vector, class ScalarVector>
void write(Vector const &u, Vector const &v, ScalarVector const &alpha,
ScalarVector const &stress);
void write(size_t record, Vector const &u, Vector const &v,
ScalarVector const &alpha, ScalarVector const &stress) const;
};
#endif
......@@ -14,5 +14,5 @@ using P1Basis = P1NodalBasis<GridView, double>;
template class MyVTKWriter<P1Basis, MyP0Basis>;
template void MyVTKWriter<P1Basis, MyP0Basis>::write<Vector, ScalarVector>(
Vector const &u, Vector const &v, ScalarVector const &alpha,
ScalarVector const &stress);
size_t record, Vector const &u, Vector const &v, ScalarVector const &alpha,
ScalarVector const &stress) const;
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