3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
nonlinear/newtonconvergencewriter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_NEWTON_CONVERGENCE_WRITER_HH
26#define DUMUX_NEWTON_CONVERGENCE_WRITER_HH
27
28#include <string>
29
30#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
32
33namespace Dumux {
34
37template <class SolutionVector>
39{
40 virtual ~ConvergenceWriterInterface() = default;
41
42 virtual void write(const SolutionVector &uLastIter, const SolutionVector &deltaU, const SolutionVector &residual) {}
43};
44
53template <class GridGeometry, class SolutionVector>
55{
56 using GridView = typename GridGeometry::GridView;
57 static constexpr auto numEq = SolutionVector::block_type::dimension;
58 using Scalar = typename SolutionVector::block_type::value_type;
59
60 static_assert(GridGeometry::discMethod != DiscretizationMethod::staggered,
61 "This convergence writer does not work for the staggered method, use the StaggeredNewtonConvergenceWriter instead");
62public:
68 NewtonConvergenceWriter(const GridGeometry& gridGeometry,
69 const std::string& name = "newton_convergence")
70 : gridGeometry_(gridGeometry)
71 , writer_(gridGeometry.gridView(), name, "", "")
72 {
73 resize();
74
75 if (GridGeometry::discMethod == DiscretizationMethod::box)
76 {
77 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
78 {
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));
82 }
83 }
84 else
85 {
86 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
87 {
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));
91 }
92 }
93 }
94
96 void resize()
97 {
98 const auto numDofs = gridGeometry_.numDofs();
99
100 // resize the output fields
101 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
102 {
103 def_[eqIdx].resize(numDofs);
104 delta_[eqIdx].resize(numDofs);
105 x_[eqIdx].resize(numDofs);
106 }
107 }
108
111 void reset(std::size_t newId = 0UL)
112 { id_ = newId; iteration_ = 0UL; }
113
114 void write(const SolutionVector& uLastIter,
115 const SolutionVector& deltaU,
116 const SolutionVector& residual) override
117 {
118 assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size());
119
120 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU.size(); ++dofIdxGlobal)
121 {
122 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
123 {
124 x_[eqIdx][dofIdxGlobal] = uLastIter[dofIdxGlobal][eqIdx];
125 delta_[eqIdx][dofIdxGlobal] = - deltaU[dofIdxGlobal][eqIdx];
126 def_[eqIdx][dofIdxGlobal] = residual[dofIdxGlobal][eqIdx];
127 }
128 }
129
130 writer_.write(static_cast<double>(id_) + static_cast<double>(iteration_)/1000);
131 ++iteration_;
132 }
133
134private:
135 std::size_t id_ = 0UL;
136 std::size_t iteration_ = 0UL;
137
138 const GridGeometry& gridGeometry_;
139
140 Dune::VTKSequenceWriter<GridView> writer_;
141
142 std::array<std::vector<Scalar>, numEq> def_;
143 std::array<std::vector<Scalar>, numEq> delta_;
144 std::array<std::vector<Scalar>, numEq> x_;
145};
146
147} // end namespace Dumux
148
149#endif
The available discretization methods in Dumux.
Definition: adapt.hh:29
Definition: nonlinear/newtonconvergencewriter.hh:39
virtual void write(const SolutionVector &uLastIter, const SolutionVector &deltaU, const SolutionVector &residual)
Definition: nonlinear/newtonconvergencewriter.hh:42
virtual ~ConvergenceWriterInterface()=default
Writes the intermediate solutions for every Newton iteration.
Definition: nonlinear/newtonconvergencewriter.hh:55
void reset(std::size_t newId=0UL)
Definition: nonlinear/newtonconvergencewriter.hh:111
void resize()
Resizes the output fields. This has to be called whenever the grid changes.
Definition: nonlinear/newtonconvergencewriter.hh:96
void write(const SolutionVector &uLastIter, const SolutionVector &deltaU, const SolutionVector &residual) override
Definition: nonlinear/newtonconvergencewriter.hh:114
NewtonConvergenceWriter(const GridGeometry &gridGeometry, const std::string &name="newton_convergence")
Constructor.
Definition: nonlinear/newtonconvergencewriter.hh:68