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

[Cleanup] Friction writing: Move half of friction writing into new class

parent 09b97b02
No related branches found
No related tags found
No related merge requests found
......@@ -3,6 +3,7 @@ bin_PROGRAMS = \
SOURCES = \
assemblers.cc \
boundary_writer.cc \
friction_writer.cc \
state/compute_state_dieterich_euler.cc \
state/compute_state_ruina.cc \
......
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
#include "boundary_writer.hh"
template <class ScalarVector, class Vector>
BoundaryWriter<ScalarVector, Vector>::BoundaryWriter(
Vector const &vertexCoordinates,
Dune::BitSetVector<1> const &_boundaryNodes)
: displacementWriter("displacements", std::fstream::out),
velocityWriter("velocities", std::fstream::out),
boundaryNodes(_boundaryNodes) {
std::fstream vertexCoordinateWriter("coordinates", std::fstream::out);
for (size_t i = 0; i < boundaryNodes.size(); ++i)
if (boundaryNodes[i][0])
vertexCoordinateWriter << vertexCoordinates[i] << std::endl;
vertexCoordinateWriter.close();
}
template <class ScalarVector, class Vector>
BoundaryWriter<ScalarVector, Vector>::~BoundaryWriter() {
displacementWriter.close();
velocityWriter.close();
}
template <class ScalarVector, class Vector>
void BoundaryWriter<ScalarVector, Vector>::writeKinetics(Vector const &u,
Vector const &v) {
for (size_t i = 0; i < boundaryNodes.size(); ++i) {
if (!boundaryNodes[i][0])
continue;
displacementWriter << u[i][0] << " ";
velocityWriter << v[i][0] << " ";
}
displacementWriter << std::endl;
velocityWriter << std::endl;
}
#include "boundary_writer_tmpl.cc"
#ifndef BOUNDARY_WRITER_HH
#define BOUNDARY_WRITER_HH
#include <fstream>
#include <string>
#include <dune/common/bitsetvector.hh>
template <class ScalarVector, class Vector> class BoundaryWriter {
public:
BoundaryWriter(Vector const &vertexCoordinates,
Dune::BitSetVector<1> const &_boundaryNodes);
virtual ~BoundaryWriter();
void writeKinetics(Vector const &u, Vector const &v);
protected:
std::fstream displacementWriter;
std::fstream velocityWriter;
Dune::BitSetVector<1> const &boundaryNodes;
};
#endif
#ifndef DIM
#error DIM unset
#endif
#include "explicitvectors.hh"
template class BoundaryWriter<ScalarVector, Vector>;
......@@ -8,41 +8,16 @@ template <class ScalarVector, class Vector>
FrictionWriter<ScalarVector, Vector>::FrictionWriter(
Vector const &vertexCoordinates,
Dune::BitSetVector<1> const &_boundaryNodes)
: coefficientWriter("coefficients", std::fstream::out),
displacementWriter("displacements", std::fstream::out),
stateWriter("states", std::fstream::out),
velocityWriter("velocities", std::fstream::out),
boundaryNodes(_boundaryNodes) {
std::fstream vertexCoordinateWriter("coordinates", std::fstream::out);
for (size_t i = 0; i < boundaryNodes.size(); ++i)
if (boundaryNodes[i][0])
vertexCoordinateWriter << vertexCoordinates[i] << std::endl;
vertexCoordinateWriter.close();
}
: BW(vertexCoordinates, _boundaryNodes),
coefficientWriter("coefficients", std::fstream::out),
stateWriter("states", std::fstream::out) {}
template <class ScalarVector, class Vector>
FrictionWriter<ScalarVector, Vector>::~FrictionWriter() {
stateWriter.close();
displacementWriter.close();
velocityWriter.close();
coefficientWriter.close();
}
template <class ScalarVector, class Vector>
void FrictionWriter<ScalarVector, Vector>::writeKinetics(Vector const &u,
Vector const &v) {
for (size_t i = 0; i < boundaryNodes.size(); ++i) {
if (!boundaryNodes[i][0])
continue;
displacementWriter << u[i][0] << " ";
velocityWriter << v[i][0] << " ";
}
displacementWriter << std::endl;
velocityWriter << std::endl;
}
template <class ScalarVector, class Vector>
void FrictionWriter<ScalarVector, Vector>::writeOther(
ScalarVector const &coefficient, ScalarVector const &alpha) {
......@@ -51,7 +26,7 @@ void FrictionWriter<ScalarVector, Vector>::writeOther(
continue;
coefficientWriter << coefficient[i] << " ";
stateWriter << alpha[i][0] << " ";
stateWriter << alpha[i] << " ";
}
stateWriter << std::endl;
......
......@@ -2,24 +2,27 @@
#define FRICTION_WRITER_HH
#include <fstream>
#include <string>
#include <dune/common/bitsetvector.hh>
template <class ScalarVector, class Vector> class FrictionWriter {
#include "boundary_writer.hh"
template <class ScalarVector, class Vector>
class FrictionWriter : public BoundaryWriter<ScalarVector, Vector> {
using BW = BoundaryWriter<ScalarVector, Vector>;
public:
FrictionWriter(Vector const &vertexCoordinates,
Dune::BitSetVector<1> const &_boundaryNodes);
~FrictionWriter();
void writeKinetics(Vector const &u, Vector const &v);
void writeOther(ScalarVector const &coefficient, ScalarVector const &alpha);
private:
std::fstream coefficientWriter;
std::fstream displacementWriter;
std::fstream stateWriter;
std::fstream velocityWriter;
Dune::BitSetVector<1> const &boundaryNodes;
using BW::boundaryNodes;
};
#endif
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