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

Externalise evolution printing

parent 4c7fee77
No related branches found
No related tags found
No related merge requests found
...@@ -23,7 +23,7 @@ run-one-body-sample-gdb: one-body-sample ...@@ -23,7 +23,7 @@ run-one-body-sample-gdb: one-body-sample
libtool --mode execute gdb ./one-body-sample libtool --mode execute gdb ./one-body-sample
one_body_sample_SOURCES = \ one_body_sample_SOURCES = \
one-body-sample.cc assemblers.cc LambertW.cc compute_state.cc vtk.cc one-body-sample.cc assemblers.cc LambertW.cc compute_state.cc vtk.cc print.cc
one_body_sample_CPPFLAGS = \ one_body_sample_CPPFLAGS = \
$(AM_CPPFLAGS) -Dsrcdir=\"$(srcdir)\" $(AM_CPPFLAGS) -Dsrcdir=\"$(srcdir)\"
one_body_sample_LDADD = \ one_body_sample_LDADD = \
......
...@@ -64,6 +64,7 @@ ...@@ -64,6 +64,7 @@
#include "assemblers.hh" #include "assemblers.hh"
#include "compute_state.hh" #include "compute_state.hh"
#include "print.hh"
#include "vtk.hh" #include "vtk.hh"
int const dim = 2; int const dim = 2;
...@@ -310,23 +311,10 @@ int main(int argc, char *argv[]) { ...@@ -310,23 +311,10 @@ int main(int argc, char *argv[]) {
} }
} }
if (parset.get<bool>("printEvolution")) { if (parset.get<bool>("printEvolution"))
// Find the first node that belongs to the frictional boundary print_evolution<SingletonVectorType, VectorType>(
for (size_t i = 0; i < frictionalNodes.size(); ++i) { frictionalNodes, *s4_new, u4, functions.get("sampleFunction"),
if (!frictionalNodes[i][0]) run, h * run, octave_writer);
continue;
boost::format const formatter("s[%03d] = %+3e, "
"%|40t|u[%03d] = %+3e");
std::cout << boost::format(formatter) % run % (*s4_new)[i] % run %
u4[i] << std::endl;
double out;
functions.get("sampleFunction").evaluate(h * run, out);
octave_writer << (*s4_new)[i] << " " << u4[i][0] * 1e6 << " " << out
<< std::endl;
break;
}
}
} }
u4 += u4_diff; u4 += u4_diff;
......
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
#include <boost/format.hpp>
#include "print.hh"
template <class SingletonVectorType, class VectorType>
void print_evolution(Dune::BitSetVector<1> const &frictionalNodes,
SingletonVectorType const &state,
VectorType const &displacement,
Dune::VirtualFunction<double, double> const &function,
int run, double time, std::ostream &octave_writer) {
// Find the first node that belongs to the frictional boundary
for (size_t i = 0; i < frictionalNodes.size(); ++i) {
if (!frictionalNodes[i][0])
continue;
double out;
function.evaluate(time, out);
std::cout << boost::format("s[%03d] = %+3e, %|20t|u[%03d] = %+3e") % run %
state[i] % run % displacement[i] << std::endl;
octave_writer << state[i][0] << " " << displacement[i][0] * 1e6 << " "
<< out << std::endl;
break;
}
}
#include "print_tmpl.cc"
#ifndef PRINT_HH
#define PRINT_HH
#include <dune/common/function.hh>
#include <dune/common/bitsetvector.hh>
template <class SingletonVectorType, class VectorType>
void print_evolution(Dune::BitSetVector<1> const &frictionalNodes,
SingletonVectorType const &state,
VectorType const &displacement,
Dune::VirtualFunction<double, double> const &function,
int run, double time, std::ostream &octave_writer);
#endif
#include <dune/common/fvector.hh>
#include <dune/istl/bvector.hh>
typedef Dune::BlockVector<Dune::FieldVector<double, 1>> SingletonVectorType;
// {{{ 2D
typedef Dune::FieldVector<double, 2> SmallVector2;
typedef Dune::BlockVector<SmallVector2> VectorType2;
template void print_evolution<SingletonVectorType, VectorType2>(
Dune::BitSetVector<1> const &frictionalNodes,
SingletonVectorType const &state, VectorType2 const &displacement,
Dune::VirtualFunction<double, double> const &function, int run, double time,
std::ostream &octave_writer);
// }}}
// {{{ 3D
typedef Dune::FieldVector<double, 3> SmallVector3;
typedef Dune::BlockVector<SmallVector3> VectorType3;
template void print_evolution<SingletonVectorType, VectorType3>(
Dune::BitSetVector<1> const &frictionalNodes,
SingletonVectorType const &state, VectorType3 const &displacement,
Dune::VirtualFunction<double, double> const &function, int run, double time,
std::ostream &octave_writer);
// }}}
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