Skip to content
Snippets Groups Projects
Commit f8522b99 authored by Kilian Weishaupt's avatar Kilian Weishaupt
Browse files

[newtonconvergencewriter] Use FVGridGeometry instead of GridView

* non-backward compatible change
* store const ref to fvGridGeometry
* call resize() without argument

* adapt test
parent 64746501
No related branches found
No related tags found
1 merge request!1529Improve Newton convergence writer, add implementation for staggered
...@@ -25,9 +25,10 @@ ...@@ -25,9 +25,10 @@
#ifndef DUMUX_NEWTON_CONVERGENCE_WRITER_HH #ifndef DUMUX_NEWTON_CONVERGENCE_WRITER_HH
#define DUMUX_NEWTON_CONVERGENCE_WRITER_HH #define DUMUX_NEWTON_CONVERGENCE_WRITER_HH
#include <dune/common/exceptions.hh> #include <string>
#include <dune/grid/io/file/vtk/vtksequencewriter.hh> #include <dune/grid/io/file/vtk/vtksequencewriter.hh>
#include <dumux/common/parameters.hh> #include <dumux/discretization/method.hh>
namespace Dumux { namespace Dumux {
...@@ -49,12 +50,15 @@ struct ConvergenceWriterInterface ...@@ -49,12 +50,15 @@ struct ConvergenceWriterInterface
* to write out multiple Newton solves with a unique id, if you don't call use all * to write out multiple Newton solves with a unique id, if you don't call use all
* Newton iterations just come after each other in the pvd file. * Newton iterations just come after each other in the pvd file.
*/ */
// template <class Scalar, class GridView, int numEq> template <class FVGridGeometry, class SolutionVector>
template <class GridView, class SolutionVector>
class NewtonConvergenceWriter : public ConvergenceWriterInterface<SolutionVector> class NewtonConvergenceWriter : public ConvergenceWriterInterface<SolutionVector>
{ {
using GridView = typename FVGridGeometry::GridView;
static constexpr auto numEq = SolutionVector::block_type::dimension; static constexpr auto numEq = SolutionVector::block_type::dimension;
using Scalar = typename SolutionVector::block_type::value_type; using Scalar = typename SolutionVector::block_type::value_type;
static_assert(FVGridGeometry::discMethod != DiscretizationMethod::staggered,
"This convergence writer does not work for the staggered method, use the StaggeredNewtonConvergenceWriter instead");
public: public:
/*! /*!
* \brief Constructor * \brief Constructor
...@@ -62,14 +66,14 @@ public: ...@@ -62,14 +66,14 @@ public:
* \param size the size of the solution data * \param size the size of the solution data
* \param name base name of the vtk output * \param name base name of the vtk output
*/ */
NewtonConvergenceWriter(const GridView& gridView, NewtonConvergenceWriter(const FVGridGeometry& fvGridGeometry,
std::size_t size,
const std::string& name = "newton_convergence") const std::string& name = "newton_convergence")
: writer_(gridView, name, "", "") : fvGridGeometry_(fvGridGeometry)
, writer_(fvGridGeometry.gridView(), name, "", "")
{ {
resize(size); resize();
if (size == gridView.size(GridView::dimension)) if (FVGridGeometry::discMethod == DiscretizationMethod::box)
{ {
for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
{ {
...@@ -78,7 +82,7 @@ public: ...@@ -78,7 +82,7 @@ public:
writer_.addVertexData(def_[eqIdx], "defect_" + std::to_string(eqIdx)); writer_.addVertexData(def_[eqIdx], "defect_" + std::to_string(eqIdx));
} }
} }
else if (size == gridView.size(0)) else
{ {
for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
{ {
...@@ -87,21 +91,19 @@ public: ...@@ -87,21 +91,19 @@ public:
writer_.addCellData(def_[eqIdx], "defect_" + std::to_string(eqIdx)); writer_.addCellData(def_[eqIdx], "defect_" + std::to_string(eqIdx));
} }
} }
else
{
DUNE_THROW(Dune::InvalidStateException, "Wrong size of output fields in the Newton convergence writer!");
}
} }
//! Resizes the output fields. This has to be called whenever the grid changes //! Resizes the output fields. This has to be called whenever the grid changes
void resize(std::size_t size) void resize()
{ {
const auto numDofs = fvGridGeometry_.numDofs();
// resize the output fields // resize the output fields
for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
{ {
def_[eqIdx].resize(size); def_[eqIdx].resize(numDofs);
delta_[eqIdx].resize(size); delta_[eqIdx].resize(numDofs);
x_[eqIdx].resize(size); x_[eqIdx].resize(numDofs);
} }
} }
...@@ -110,9 +112,9 @@ public: ...@@ -110,9 +112,9 @@ public:
void reset(std::size_t newId = 0UL) void reset(std::size_t newId = 0UL)
{ id_ = newId; iteration_ = 0UL; } { id_ = newId; iteration_ = 0UL; }
void write(const SolutionVector &uLastIter, void write(const SolutionVector& uLastIter,
const SolutionVector &deltaU, const SolutionVector& deltaU,
const SolutionVector &residual) override const SolutionVector& residual) override
{ {
assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size()); assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size());
...@@ -134,6 +136,8 @@ private: ...@@ -134,6 +136,8 @@ private:
std::size_t id_ = 0UL; std::size_t id_ = 0UL;
std::size_t iteration_ = 0UL; std::size_t iteration_ = 0UL;
const FVGridGeometry& fvGridGeometry_;
Dune::VTKSequenceWriter<GridView> writer_; Dune::VTKSequenceWriter<GridView> writer_;
std::array<std::vector<Scalar>, numEq> def_; std::array<std::vector<Scalar>, numEq> def_;
......
...@@ -155,9 +155,8 @@ int main(int argc, char** argv) try ...@@ -155,9 +155,8 @@ int main(int argc, char** argv) try
NewtonSolver nonLinearSolver(assembler, linearSolver); NewtonSolver nonLinearSolver(assembler, linearSolver);
//the convergence writer //the convergence writer
using GridView = GetPropType<TypeTag, Properties::GridView>; using NewtonConvergenceWriter = Dumux::NewtonConvergenceWriter<FVGridGeometry, SolutionVector>;
using NewtonConvergenceWriter = Dumux::NewtonConvergenceWriter<GridView, SolutionVector>; auto convergenceWriter = std::make_shared<NewtonConvergenceWriter>(*fvGridGeometry);
auto convergenceWriter = std::make_shared<NewtonConvergenceWriter>(leafGridView, fvGridGeometry->numDofs());
nonLinearSolver.attachConvergenceWriter(convergenceWriter); nonLinearSolver.attachConvergenceWriter(convergenceWriter);
// time loop // time loop
......
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