version 3.8
staggerednewtonconvergencewriter.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_STAGGERED_NEWTON_CONVERGENCE_WRITER_HH
15#define DUMUX_STAGGERED_NEWTON_CONVERGENCE_WRITER_HH
16
17#include <string>
18
19#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
23
24namespace Dumux {
25
34template <class GridGeometry, class SolutionVector, class ResidualVector>
35class StaggeredNewtonConvergenceWriter : public ConvergenceWriterInterface<SolutionVector, ResidualVector>
36{
37 using GridView = typename GridGeometry::GridView;
38
39 using CellCenterSolutionVector = typename std::decay_t<decltype(std::declval<SolutionVector>()[GridGeometry::cellCenterIdx()])>;
40 using FaceSolutionVector = typename std::decay_t<decltype(std::declval<SolutionVector>()[GridGeometry::faceIdx()])>;
41
42 using Scalar = typename CellCenterSolutionVector::block_type::value_type;
43
44 static constexpr auto numEqCellCenter = CellCenterSolutionVector::block_type::dimension;
45 static constexpr auto numEqFace = FaceSolutionVector::block_type::dimension;
46
47 using Element = typename GridView::template Codim<0>::Entity;
48 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
49
50 static_assert(GridGeometry::discMethod == DiscretizationMethods::staggered,
51 "This convergence writer does only work for the staggered method, use the NewtonConvergenceWriter instead");
52public:
58 StaggeredNewtonConvergenceWriter(const GridGeometry& gridGeometry,
59 const std::string& name = "newton_convergence")
60 : gridGeometry_(gridGeometry)
61 , ccWriter_(gridGeometry.gridView(), name, "", "")
62 , faceWriter_(std::make_shared<PointCloudVtkWriter<Scalar, GlobalPosition>>(coordinates_))
63 , faceSequenceWriter_(faceWriter_, name + "-face", "","",
64 gridGeometry.gridView().comm().rank(),
65 gridGeometry.gridView().comm().size())
66 {
67 resize();
68
69 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
70 {
71 ccWriter_.addCellData(xCellCenter_[eqIdx], "x_" + std::to_string(eqIdx));
72 ccWriter_.addCellData(deltaCellCenter_[eqIdx], "delta_" + std::to_string(eqIdx));
73 ccWriter_.addCellData(defCellCenter_[eqIdx], "defect_" + std::to_string(eqIdx));
74 }
75 }
76
78 void resize()
79 {
80 const auto numCellCenterDofs = gridGeometry_.numCellCenterDofs();
81 const auto numFaceDofs = gridGeometry_.numFaceDofs();
82
83 // resize the cell center output fields
84 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
85 {
86 defCellCenter_[eqIdx].resize(numCellCenterDofs);
87 deltaCellCenter_[eqIdx].resize(numCellCenterDofs);
88 xCellCenter_[eqIdx].resize(numCellCenterDofs);
89 }
90
91 // resize the face output fields
92 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
93 {
94 defFace_[eqIdx].resize(numFaceDofs);
95 deltaFace_[eqIdx].resize(numFaceDofs);
96 xFace_[eqIdx].resize(numFaceDofs);
97 }
98
99 coordinates_.resize(numFaceDofs);
100 for (auto&& facet : facets(gridGeometry_.gridView()))
101 {
102 const auto dofIdxGlobal = gridGeometry_.gridView().indexSet().index(facet);
103 coordinates_[dofIdxGlobal] = facet.geometry().center();
104 }
105 }
106
109 void reset(std::size_t newId = 0UL)
110 { id_ = newId; iteration_ = 0UL; }
111
112 void write(const SolutionVector& uLastIter,
113 const ResidualVector& deltaU,
114 const ResidualVector& residual) override
115 {
116 assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size());
117
118 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU[GridGeometry::cellCenterIdx()].size(); ++dofIdxGlobal)
119 {
120 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
121 {
122 xCellCenter_[eqIdx][dofIdxGlobal] = uLastIter[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
123 deltaCellCenter_[eqIdx][dofIdxGlobal] = - deltaU[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
124 defCellCenter_[eqIdx][dofIdxGlobal] = residual[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
125 }
126 }
127
128 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
129 {
130 faceWriter_->addPointData(xFace_[eqIdx], "x_" + std::to_string(eqIdx));
131 faceWriter_->addPointData(deltaFace_[eqIdx], "delta_" + std::to_string(eqIdx));
132 faceWriter_->addPointData(defFace_[eqIdx], "defect_" + std::to_string(eqIdx));
133 }
134
135 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU[GridGeometry::faceIdx()].size(); ++dofIdxGlobal)
136 {
137 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
138 {
139 xFace_[eqIdx][dofIdxGlobal] = uLastIter[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
140 deltaFace_[eqIdx][dofIdxGlobal] = - deltaU[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
141 defFace_[eqIdx][dofIdxGlobal] = residual[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
142 }
143 }
144
145 ccWriter_.write(static_cast<double>(id_) + static_cast<double>(iteration_)/1000);
146 faceSequenceWriter_.write(static_cast<double>(id_) + static_cast<double>(iteration_)/1000);
147 ++iteration_;
148 }
149
150private:
151 std::size_t id_ = 0UL;
152 std::size_t iteration_ = 0UL;
153
154 const GridGeometry& gridGeometry_;
155
156 Dune::VTKSequenceWriter<GridView> ccWriter_;
157
158 std::vector<GlobalPosition> coordinates_;
159 std::shared_ptr<PointCloudVtkWriter<Scalar, GlobalPosition>> faceWriter_;
161
162 std::array<std::vector<Scalar>, numEqCellCenter> defCellCenter_;
163 std::array<std::vector<Scalar>, numEqCellCenter> deltaCellCenter_;
164 std::array<std::vector<Scalar>, numEqCellCenter> xCellCenter_;
165
166 std::array<std::vector<Scalar>, numEqFace> defFace_;
167 std::array<std::vector<Scalar>, numEqFace> deltaFace_;
168 std::array<std::vector<Scalar>, numEqFace> xFace_;
169};
170
171} // end namespace Dumux
172
173#endif
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:39
Writes the intermediate solutions for every Newton iteration (for staggered grid scheme)
Definition: staggerednewtonconvergencewriter.hh:36
StaggeredNewtonConvergenceWriter(const GridGeometry &gridGeometry, const std::string &name="newton_convergence")
Constructor.
Definition: staggerednewtonconvergencewriter.hh:58
void reset(std::size_t newId=0UL)
Definition: staggerednewtonconvergencewriter.hh:109
void resize()
Resizes the output fields. This has to be called whenever the grid changes.
Definition: staggerednewtonconvergencewriter.hh:78
void write(const SolutionVector &uLastIter, const ResidualVector &deltaU, const ResidualVector &residual) override
Definition: staggerednewtonconvergencewriter.hh:112
Base class to write pvd-files which contains a list of all collected vtk-files. This is a modified ve...
Definition: vtksequencewriter.hh:44
constexpr Staggered staggered
Definition: method.hh:149
Definition: adapt.hh:17
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
A VTK writer specialized for staggered grid implementations with dofs on the faces.
Definition: nonlinear/newtonconvergencewriter.hh:27
Base class to write pvd-files which contains a list of all collected vtk-files. This is a modified ve...