From f8522b998f9902fa2b0da5517f36a0757b8f43b1 Mon Sep 17 00:00:00 2001 From: Kilian Weishaupt <kilian.weishaupt@iws.uni-stuttgart.de> Date: Tue, 26 Mar 2019 12:27:34 +0100 Subject: [PATCH] [newtonconvergencewriter] Use FVGridGeometry instead of GridView * non-backward compatible change * store const ref to fvGridGeometry * call resize() without argument * adapt test --- dumux/nonlinear/newtonconvergencewriter.hh | 46 ++++++++++--------- .../2p/implicit/nonisothermal/main.cc | 5 +- 2 files changed, 27 insertions(+), 24 deletions(-) diff --git a/dumux/nonlinear/newtonconvergencewriter.hh b/dumux/nonlinear/newtonconvergencewriter.hh index 0ebfcb6577..6d06350ed5 100644 --- a/dumux/nonlinear/newtonconvergencewriter.hh +++ b/dumux/nonlinear/newtonconvergencewriter.hh @@ -25,9 +25,10 @@ #ifndef 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 <dumux/common/parameters.hh> +#include <dumux/discretization/method.hh> namespace Dumux { @@ -49,12 +50,15 @@ struct ConvergenceWriterInterface * 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. */ -// template <class Scalar, class GridView, int numEq> -template <class GridView, class SolutionVector> +template <class FVGridGeometry, class SolutionVector> class NewtonConvergenceWriter : public ConvergenceWriterInterface<SolutionVector> { + using GridView = typename FVGridGeometry::GridView; static constexpr auto numEq = SolutionVector::block_type::dimension; 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: /*! * \brief Constructor @@ -62,14 +66,14 @@ public: * \param size the size of the solution data * \param name base name of the vtk output */ - NewtonConvergenceWriter(const GridView& gridView, - std::size_t size, + NewtonConvergenceWriter(const FVGridGeometry& fvGridGeometry, 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) { @@ -78,7 +82,7 @@ public: writer_.addVertexData(def_[eqIdx], "defect_" + std::to_string(eqIdx)); } } - else if (size == gridView.size(0)) + else { for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) { @@ -87,21 +91,19 @@ public: 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 - void resize(std::size_t size) + void resize() { + const auto numDofs = fvGridGeometry_.numDofs(); + // resize the output fields for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) { - def_[eqIdx].resize(size); - delta_[eqIdx].resize(size); - x_[eqIdx].resize(size); + def_[eqIdx].resize(numDofs); + delta_[eqIdx].resize(numDofs); + x_[eqIdx].resize(numDofs); } } @@ -110,9 +112,9 @@ public: void reset(std::size_t newId = 0UL) { id_ = newId; iteration_ = 0UL; } - void write(const SolutionVector &uLastIter, - const SolutionVector &deltaU, - const SolutionVector &residual) override + void write(const SolutionVector& uLastIter, + const SolutionVector& deltaU, + const SolutionVector& residual) override { assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size()); @@ -134,6 +136,8 @@ private: std::size_t id_ = 0UL; std::size_t iteration_ = 0UL; + const FVGridGeometry& fvGridGeometry_; + Dune::VTKSequenceWriter<GridView> writer_; std::array<std::vector<Scalar>, numEq> def_; diff --git a/test/porousmediumflow/2p/implicit/nonisothermal/main.cc b/test/porousmediumflow/2p/implicit/nonisothermal/main.cc index 21e0366c7e..3c5aa3face 100644 --- a/test/porousmediumflow/2p/implicit/nonisothermal/main.cc +++ b/test/porousmediumflow/2p/implicit/nonisothermal/main.cc @@ -155,9 +155,8 @@ int main(int argc, char** argv) try NewtonSolver nonLinearSolver(assembler, linearSolver); //the convergence writer - using GridView = GetPropType<TypeTag, Properties::GridView>; - using NewtonConvergenceWriter = Dumux::NewtonConvergenceWriter<GridView, SolutionVector>; - auto convergenceWriter = std::make_shared<NewtonConvergenceWriter>(leafGridView, fvGridGeometry->numDofs()); + using NewtonConvergenceWriter = Dumux::NewtonConvergenceWriter<FVGridGeometry, SolutionVector>; + auto convergenceWriter = std::make_shared<NewtonConvergenceWriter>(*fvGridGeometry); nonLinearSolver.attachConvergenceWriter(convergenceWriter); // time loop -- GitLab