version 3.10-dev
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
32template <class GridGeometry, class SolutionVector, class ResidualVector>
33class StaggeredNewtonConvergenceWriter : public ConvergenceWriterInterface<SolutionVector, ResidualVector>
34{
35 using GridView = typename GridGeometry::GridView;
36
37 using CellCenterSolutionVector = typename std::decay_t<decltype(std::declval<SolutionVector>()[GridGeometry::cellCenterIdx()])>;
38 using FaceSolutionVector = typename std::decay_t<decltype(std::declval<SolutionVector>()[GridGeometry::faceIdx()])>;
39
40 using Scalar = typename CellCenterSolutionVector::block_type::value_type;
41
42 static constexpr auto numEqCellCenter = CellCenterSolutionVector::block_type::dimension;
43 static constexpr auto numEqFace = FaceSolutionVector::block_type::dimension;
44
45 using Element = typename GridView::template Codim<0>::Entity;
46 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
47
48 static_assert(GridGeometry::discMethod == DiscretizationMethods::staggered,
49 "This convergence writer does only work for the staggered method, use the NewtonConvergenceWriter instead");
50public:
56 StaggeredNewtonConvergenceWriter(const GridGeometry& gridGeometry,
57 const std::string& name = "newton_convergence")
58 : gridGeometry_(gridGeometry)
59 , ccWriter_(gridGeometry.gridView(), name, "", "")
60 , faceWriter_(std::make_shared<PointCloudVtkWriter<Scalar, GlobalPosition>>(coordinates_))
61 , faceSequenceWriter_(faceWriter_, name + "-face", "","",
62 gridGeometry.gridView().comm().rank(),
63 gridGeometry.gridView().comm().size())
64 {
65 resize();
66
67 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
68 {
69 ccWriter_.addCellData(xCellCenter_[eqIdx], "x_" + std::to_string(eqIdx));
70 ccWriter_.addCellData(deltaCellCenter_[eqIdx], "delta_" + std::to_string(eqIdx));
71 ccWriter_.addCellData(defCellCenter_[eqIdx], "defect_" + std::to_string(eqIdx));
72 }
73 }
74
76 void resize()
77 {
78 const auto numCellCenterDofs = gridGeometry_.numCellCenterDofs();
79 const auto numFaceDofs = gridGeometry_.numFaceDofs();
80
81 // resize the cell center output fields
82 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
83 {
84 defCellCenter_[eqIdx].resize(numCellCenterDofs);
85 deltaCellCenter_[eqIdx].resize(numCellCenterDofs);
86 xCellCenter_[eqIdx].resize(numCellCenterDofs);
87 }
88
89 // resize the face output fields
90 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
91 {
92 defFace_[eqIdx].resize(numFaceDofs);
93 deltaFace_[eqIdx].resize(numFaceDofs);
94 xFace_[eqIdx].resize(numFaceDofs);
95 }
96
97 coordinates_.resize(numFaceDofs);
98 for (auto&& facet : facets(gridGeometry_.gridView()))
99 {
100 const auto dofIdxGlobal = gridGeometry_.gridView().indexSet().index(facet);
101 coordinates_[dofIdxGlobal] = facet.geometry().center();
102 }
103 }
104
107 void reset(std::size_t newId = 0UL)
108 { id_ = newId; iteration_ = 0UL; }
109
110 void write(const SolutionVector& uLastIter,
111 const ResidualVector& deltaU,
112 const ResidualVector& residual) override
113 {
114 assert(uLastIter.size() == deltaU.size() && uLastIter.size() == residual.size());
115
116 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU[GridGeometry::cellCenterIdx()].size(); ++dofIdxGlobal)
117 {
118 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
119 {
120 xCellCenter_[eqIdx][dofIdxGlobal] = uLastIter[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
121 deltaCellCenter_[eqIdx][dofIdxGlobal] = - deltaU[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
122 defCellCenter_[eqIdx][dofIdxGlobal] = residual[GridGeometry::cellCenterIdx()][dofIdxGlobal][eqIdx];
123 }
124 }
125
126 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
127 {
128 faceWriter_->addPointData(xFace_[eqIdx], "x_" + std::to_string(eqIdx));
129 faceWriter_->addPointData(deltaFace_[eqIdx], "delta_" + std::to_string(eqIdx));
130 faceWriter_->addPointData(defFace_[eqIdx], "defect_" + std::to_string(eqIdx));
131 }
132
133 for (std::size_t dofIdxGlobal = 0; dofIdxGlobal < deltaU[GridGeometry::faceIdx()].size(); ++dofIdxGlobal)
134 {
135 for (int eqIdx = 0; eqIdx < numEqFace; ++eqIdx)
136 {
137 xFace_[eqIdx][dofIdxGlobal] = uLastIter[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
138 deltaFace_[eqIdx][dofIdxGlobal] = - deltaU[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
139 defFace_[eqIdx][dofIdxGlobal] = residual[GridGeometry::faceIdx()][dofIdxGlobal][eqIdx];
140 }
141 }
142
143 ccWriter_.write(static_cast<double>(id_) + static_cast<double>(iteration_)/1000);
144 faceSequenceWriter_.write(static_cast<double>(id_) + static_cast<double>(iteration_)/1000);
145 ++iteration_;
146 }
147
148private:
149 std::size_t id_ = 0UL;
150 std::size_t iteration_ = 0UL;
151
152 const GridGeometry& gridGeometry_;
153
154 Dune::VTKSequenceWriter<GridView> ccWriter_;
155
156 std::vector<GlobalPosition> coordinates_;
157 std::shared_ptr<PointCloudVtkWriter<Scalar, GlobalPosition>> faceWriter_;
159
160 std::array<std::vector<Scalar>, numEqCellCenter> defCellCenter_;
161 std::array<std::vector<Scalar>, numEqCellCenter> deltaCellCenter_;
162 std::array<std::vector<Scalar>, numEqCellCenter> xCellCenter_;
163
164 std::array<std::vector<Scalar>, numEqFace> defFace_;
165 std::array<std::vector<Scalar>, numEqFace> deltaFace_;
166 std::array<std::vector<Scalar>, numEqFace> xFace_;
167};
168
169} // end namespace Dumux
170
171#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:34
StaggeredNewtonConvergenceWriter(const GridGeometry &gridGeometry, const std::string &name="newton_convergence")
Constructor.
Definition: staggerednewtonconvergencewriter.hh:56
void reset(std::size_t newId=0UL)
Definition: staggerednewtonconvergencewriter.hh:107
void resize()
Resizes the output fields. This has to be called whenever the grid changes.
Definition: staggerednewtonconvergencewriter.hh:76
void write(const SolutionVector &uLastIter, const ResidualVector &deltaU, const ResidualVector &residual) override
Definition: staggerednewtonconvergencewriter.hh:110
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.
A convergence writer interface Provide an interface that show the minimal requirements a convergence ...
Definition: nonlinear/newtonconvergencewriter.hh:32
Base class to write pvd-files which contains a list of all collected vtk-files. This is a modified ve...