25#ifndef DUMUX_NEWTON_CONVERGENCE_WRITER_HH
26#define DUMUX_NEWTON_CONVERGENCE_WRITER_HH
30#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
37template <
class SolutionVector>
42 virtual void write(
const SolutionVector &uLastIter,
const SolutionVector &deltaU,
const SolutionVector &residual) {}
53template <
class Gr
idGeometry,
class SolutionVector>
56 using GridView =
typename GridGeometry::GridView;
57 static constexpr auto numEq = SolutionVector::block_type::dimension;
58 using Scalar =
typename SolutionVector::block_type::value_type;
61 "This convergence writer does not work for the staggered method, use the StaggeredNewtonConvergenceWriter instead");
69 const std::string& name =
"newton_convergence")
70 : gridGeometry_(gridGeometry)
71 , writer_(gridGeometry.gridView(), name,
"",
"")
77 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
79 writer_.addVertexData(x_[eqIdx],
"x_" + std::to_string(eqIdx));
80 writer_.addVertexData(delta_[eqIdx],
"delta_" + std::to_string(eqIdx));
81 writer_.addVertexData(def_[eqIdx],
"defect_" + std::to_string(eqIdx));
86 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
88 writer_.addCellData(x_[eqIdx],
"x_" + std::to_string(eqIdx));
89 writer_.addCellData(delta_[eqIdx],
"delta_" + std::to_string(eqIdx));
90 writer_.addCellData(def_[eqIdx],
"defect_" + std::to_string(eqIdx));
98 const auto numDofs = gridGeometry_.numDofs();
101 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
103 def_[eqIdx].resize(numDofs);
104 delta_[eqIdx].resize(numDofs);
105 x_[eqIdx].resize(numDofs);
112 { id_ = newId; iteration_ = 0UL; }
114 void write(
const SolutionVector& uLastIter,
115 const SolutionVector& deltaU,
116 const SolutionVector& residual)
override
118 assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size());
120 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU.size(); ++dofIdxGlobal)
122 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
124 x_[eqIdx][dofIdxGlobal] = uLastIter[dofIdxGlobal][eqIdx];
125 delta_[eqIdx][dofIdxGlobal] = - deltaU[dofIdxGlobal][eqIdx];
126 def_[eqIdx][dofIdxGlobal] = residual[dofIdxGlobal][eqIdx];
130 writer_.write(
static_cast<double>(id_) +
static_cast<double>(iteration_)/1000);
135 std::size_t id_ = 0UL;
136 std::size_t iteration_ = 0UL;
138 const GridGeometry& gridGeometry_;
140 Dune::VTKSequenceWriter<GridView> writer_;
142 std::array<std::vector<Scalar>, numEq> def_;
143 std::array<std::vector<Scalar>, numEq> delta_;
144 std::array<std::vector<Scalar>, numEq> x_;
The available discretization methods in Dumux.
Definition: newtonconvergencewriter.hh:39
virtual void write(const SolutionVector &uLastIter, const SolutionVector &deltaU, const SolutionVector &residual)
Definition: newtonconvergencewriter.hh:42
virtual ~ConvergenceWriterInterface()=default
Writes the intermediate solutions for every Newton iteration.
Definition: newtonconvergencewriter.hh:55
void reset(std::size_t newId=0UL)
Definition: newtonconvergencewriter.hh:111
void resize()
Resizes the output fields. This has to be called whenever the grid changes.
Definition: newtonconvergencewriter.hh:96
void write(const SolutionVector &uLastIter, const SolutionVector &deltaU, const SolutionVector &residual) override
Definition: newtonconvergencewriter.hh:114
NewtonConvergenceWriter(const GridGeometry &gridGeometry, const std::string &name="newton_convergence")
Constructor.
Definition: newtonconvergencewriter.hh:68