Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Experimental | ||
10 | * \ingroup Assembly | ||
11 | * \brief A linear system assembler (residual and Jacobian) for finite volume schemes | ||
12 | */ | ||
13 | #ifndef DUMUX_EXPERIMENTAL_MULTISTAGE_FV_ASSEMBLER_HH | ||
14 | #define DUMUX_EXPERIMENTAL_MULTISTAGE_FV_ASSEMBLER_HH | ||
15 | |||
16 | #include <vector> | ||
17 | #include <deque> | ||
18 | #include <type_traits> | ||
19 | #include <memory> | ||
20 | |||
21 | #include <dune/istl/matrixindexset.hh> | ||
22 | |||
23 | #include <dumux/common/properties.hh> | ||
24 | #include <dumux/common/gridcapabilities.hh> | ||
25 | #include <dumux/common/typetraits/vector.hh> | ||
26 | |||
27 | #include <dumux/discretization/method.hh> | ||
28 | #include <dumux/linear/parallelhelpers.hh> | ||
29 | #include <dumux/linear/dunevectors.hh> | ||
30 | |||
31 | #include <dumux/assembly/coloring.hh> | ||
32 | #include <dumux/assembly/jacobianpattern.hh> | ||
33 | #include <dumux/assembly/diffmethod.hh> | ||
34 | |||
35 | #include <dumux/parallel/multithreading.hh> | ||
36 | #include <dumux/parallel/parallel_for.hh> | ||
37 | |||
38 | #include <dumux/experimental/assembly/cvfelocalassembler.hh> | ||
39 | #include <dumux/experimental/assembly/cclocalassembler.hh> | ||
40 | #include <dumux/experimental/assembly/multistagefvlocaloperator.hh> | ||
41 | |||
42 | #include <dumux/experimental/timestepping/multistagemethods.hh> | ||
43 | #include <dumux/experimental/timestepping/multistagetimestepper.hh> | ||
44 | |||
45 | namespace Dumux::Experimental::Detail { | ||
46 | |||
47 | template<class DiscretizationMethod> | ||
48 | struct LocalAssemblerChooser; | ||
49 | |||
50 | template<class DM> | ||
51 | struct LocalAssemblerChooser<DiscretizationMethods::CVFE<DM>> | ||
52 | { | ||
53 | template<class TypeTag, class Impl, DiffMethod diffMethod> | ||
54 | using type = Experimental::CVFELocalAssembler<TypeTag, Impl, diffMethod>; | ||
55 | }; | ||
56 | |||
57 | template<> | ||
58 | struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa> | ||
59 | { | ||
60 | template<class TypeTag, class Impl, DiffMethod diffMethod> | ||
61 | using type = Experimental::CCLocalAssembler<TypeTag, Impl, diffMethod>; | ||
62 | }; | ||
63 | |||
64 | template<> | ||
65 | struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa> | ||
66 | { | ||
67 | template<class TypeTag, class Impl, DiffMethod diffMethod> | ||
68 | using type = Experimental::CCLocalAssembler<TypeTag, Impl, diffMethod>; | ||
69 | }; | ||
70 | |||
71 | template<class TypeTag, class Impl, DiffMethod diffMethod> | ||
72 | using LocalAssemblerChooser_t = typename LocalAssemblerChooser< | ||
73 | typename GetPropType<TypeTag, Properties::GridGeometry>::DiscretizationMethod | ||
74 | >::template type<TypeTag, Impl, diffMethod>; | ||
75 | |||
76 | } // end namespace Dumux::Detail | ||
77 | |||
78 | namespace Dumux::Experimental { | ||
79 | |||
80 | /*! | ||
81 | * \ingroup Experimental | ||
82 | * \ingroup Assembly | ||
83 | * \brief A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa, mpfa, ...) | ||
84 | * \tparam TypeTag The TypeTag | ||
85 | * \tparam diffMethod The differentiation method to residual compute derivatives | ||
86 | */ | ||
87 | template<class TypeTag, DiffMethod diffMethod> | ||
88 | class MultiStageFVAssembler | ||
89 | { | ||
90 | using GridGeo = GetPropType<TypeTag, Properties::GridGeometry>; | ||
91 | using GridView = typename GridGeo::GridView; | ||
92 | using LocalResidual = GetPropType<TypeTag, Properties::LocalResidual>; | ||
93 | using Element = typename GridView::template Codim<0>::Entity; | ||
94 | using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed; | ||
95 | |||
96 | static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box; | ||
97 | |||
98 | using ThisType = MultiStageFVAssembler<TypeTag, diffMethod>; | ||
99 | using LocalAssembler = typename Detail::LocalAssemblerChooser_t<TypeTag, ThisType, diffMethod>; | ||
100 | |||
101 | public: | ||
102 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
103 | using StageParams = Experimental::MultiStageParams<Scalar>; | ||
104 | using JacobianMatrix = GetPropType<TypeTag, Properties::JacobianMatrix>; | ||
105 | using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>; | ||
106 | using ResidualType = typename Dumux::Detail::NativeDuneVectorType<SolutionVector>::type; | ||
107 | |||
108 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
109 | |||
110 | using GridGeometry = GridGeo; | ||
111 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
112 | |||
113 | /*! | ||
114 | * \brief The constructor for instationary problems | ||
115 | * \note the grid variables might be temporarily changed during assembly (if caching is enabled) | ||
116 | * it is however guaranteed that the state after assembly will be the same as before | ||
117 | */ | ||
118 | 8 | MultiStageFVAssembler(std::shared_ptr<const Problem> problem, | |
119 | std::shared_ptr<const GridGeometry> gridGeometry, | ||
120 | std::shared_ptr<GridVariables> gridVariables, | ||
121 | std::shared_ptr<const Experimental::MultiStageMethod<Scalar>> msMethod, | ||
122 | const SolutionVector& prevSol) | ||
123 | : timeSteppingMethod_(msMethod) | ||
124 | , problem_(problem) | ||
125 | , gridGeometry_(gridGeometry) | ||
126 | , gridVariables_(gridVariables) | ||
127 |
2/20✓ Branch 9 taken 8 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 8 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
8 | , prevSol_(&prevSol) |
128 | { | ||
129 | 8 | enableMultithreading_ = SupportsColoring<typename GridGeometry::DiscretizationMethod>::value | |
130 |
3/6✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
|
24 | && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView()) |
131 | && !Multithreading::isSerial() | ||
132 |
2/4✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
|
8 | && getParam<bool>("Assembly.Multithreading", true); |
133 | |||
134 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | maybeComputeColors_(); |
135 | 8 | } | |
136 | |||
137 | /*! | ||
138 | * \brief Assembles the global Jacobian of the residual | ||
139 | * and the residual for the current solution. | ||
140 | */ | ||
141 | template<class PartialReassembler = DefaultPartialReassembler> | ||
142 | 2158 | void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr) | |
143 | { | ||
144 | 2158 | resetJacobian_(partialReassembler); | |
145 | |||
146 | 2158 | resetResidual_(); | |
147 | 4316 | spatialOperatorEvaluations_.back() = 0.0; | |
148 | 4316 | temporalOperatorEvaluations_.back() = 0.0; | |
149 | |||
150 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 2158 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2158 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2158 times.
|
6474 | if (stageParams_->size() != spatialOperatorEvaluations_.size()) |
151 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Wrong number of residuals"); | |
152 | |||
153 | 2158 | assemble_([&](const Element& element) | |
154 | { | ||
155 |
1/4✓ Branch 1 taken 10368 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
12641696 | LocalAssembler localAssembler(*this, element, curSol); |
156 |
1/2✓ Branch 1 taken 5060672 times.
✗ Branch 2 not taken.
|
5060672 | localAssembler.assembleJacobianAndResidual( |
157 |
3/6✓ Branch 1 taken 5060672 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5060672 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 5060672 times.
✗ Branch 8 not taken.
|
15182016 | *jacobian_, *residual_, *gridVariables_, |
158 |
2/4✓ Branch 1 taken 5060672 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5060672 times.
✗ Branch 5 not taken.
|
10121344 | *stageParams_, |
159 | temporalOperatorEvaluations_.back(), | ||
160 | spatialOperatorEvaluations_.back(), | ||
161 | 5060672 | constrainedDofs_, | |
162 |
1/2✓ Branch 1 taken 5060672 times.
✗ Branch 2 not taken.
|
5060672 | partialReassembler |
163 | ); | ||
164 | }); | ||
165 | |||
166 | // assemble the full residual for the time integration stage | ||
167 |
1/2✓ Branch 2 taken 2158 times.
✗ Branch 3 not taken.
|
6474 | auto constantResidualComponent = (*residual_); |
168 | constantResidualComponent = 0.0; | ||
169 |
4/4✓ Branch 0 taken 2158 times.
✓ Branch 1 taken 2158 times.
✓ Branch 2 taken 2158 times.
✓ Branch 3 taken 2158 times.
|
6474 | for (std::size_t k = 0; k < stageParams_->size()-1; ++k) |
170 | { | ||
171 |
3/6✓ Branch 0 taken 2158 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2158 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2158 times.
✗ Branch 5 not taken.
|
6474 | if (!stageParams_->skipTemporal(k)) |
172 | 8632 | constantResidualComponent.axpy(stageParams_->temporalWeight(k), temporalOperatorEvaluations_[k]); | |
173 |
6/6✓ Branch 0 taken 2000 times.
✓ Branch 1 taken 158 times.
✓ Branch 2 taken 2000 times.
✓ Branch 3 taken 158 times.
✓ Branch 4 taken 2000 times.
✓ Branch 5 taken 158 times.
|
6474 | if (!stageParams_->skipSpatial(k)) |
174 | 8000 | constantResidualComponent.axpy(stageParams_->spatialWeight(k), spatialOperatorEvaluations_[k]); | |
175 | } | ||
176 | |||
177 | // masked summation of constant residual component onto this stage's resiudal component | ||
178 |
2/2✓ Branch 0 taken 5062886 times.
✓ Branch 1 taken 2158 times.
|
5065044 | for (std::size_t i = 0; i < constantResidualComponent.size(); ++i) |
179 |
4/4✓ Branch 0 taken 10125772 times.
✓ Branch 1 taken 5062886 times.
✓ Branch 2 taken 10125772 times.
✓ Branch 3 taken 5062886 times.
|
30377316 | for (std::size_t ii = 0; ii < constantResidualComponent[i].size(); ++ii) |
180 |
6/6✓ Branch 0 taken 10123936 times.
✓ Branch 1 taken 1836 times.
✓ Branch 2 taken 10123936 times.
✓ Branch 3 taken 1836 times.
✓ Branch 4 taken 10123936 times.
✓ Branch 5 taken 1836 times.
|
30377316 | (*residual_)[i][ii] += constrainedDofs_[i][ii] > 0.5 ? 0.0 : constantResidualComponent[i][ii]; |
181 | 2158 | } | |
182 | |||
183 | //! compute the residuals using the internal residual | ||
184 | ✗ | void assembleResidual(const SolutionVector& curSol) | |
185 | ✗ | { DUNE_THROW(Dune::NotImplemented, "residual"); } | |
186 | |||
187 | /*! | ||
188 | * \brief The version without arguments uses the default constructor to create | ||
189 | * the jacobian and residual objects in this assembler if you don't need them outside this class | ||
190 | */ | ||
191 | 8 | void setLinearSystem() | |
192 | { | ||
193 |
2/4✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
|
8 | jacobian_ = std::make_shared<JacobianMatrix>(); |
194 | 16 | jacobian_->setBuildMode(JacobianMatrix::random); | |
195 |
2/4✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
|
8 | residual_ = std::make_shared<ResidualType>(); |
196 | |||
197 | 14 | setResidualSize_(*residual_); | |
198 | 8 | setJacobianPattern_(); | |
199 | 8 | } | |
200 | |||
201 | /*! | ||
202 | * \brief Resizes jacobian and residual and recomputes colors | ||
203 | */ | ||
204 | void updateAfterGridAdaption() | ||
205 | { | ||
206 | setResidualSize_(*residual_); | ||
207 | setJacobianPattern_(); | ||
208 | maybeComputeColors_(); | ||
209 | } | ||
210 | |||
211 | //! Returns the number of degrees of freedom | ||
212 | std::size_t numDofs() const | ||
213 | 6270 | { return gridGeometry_->numDofs(); } | |
214 | |||
215 | //! The problem | ||
216 | const Problem& problem() const | ||
217 |
14/33✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 4 times.
✗ Branch 18 not taken.
✓ Branch 20 taken 4 times.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 4 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✓ Branch 26 taken 4 times.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✓ Branch 33 taken 2 times.
✗ Branch 34 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✗ Branch 46 not taken.
✗ Branch 47 not taken.
|
23688164 | { return *problem_; } |
218 | |||
219 | //! The global finite volume geometry | ||
220 | const GridGeometry& gridGeometry() const | ||
221 |
6/17✓ Branch 1 taken 10000000 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 10000000 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 52224 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 52224 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 13440 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 13440 times.
✗ Branch 19 not taken.
|
20131328 | { return *gridGeometry_; } |
222 | |||
223 | //! The gridview | ||
224 | const GridView& gridView() const | ||
225 |
4/16✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2158 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 5060672 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 14 taken 2048 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 5018432 times.
|
10083310 | { return gridGeometry().gridView(); } |
226 | |||
227 | //! The global grid variables | ||
228 | GridVariables& gridVariables() | ||
229 | 4412 | { return *gridVariables_; } | |
230 | |||
231 | //! The global grid variables | ||
232 | const GridVariables& gridVariables() const | ||
233 |
4/8✓ Branch 1 taken 10065664 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 10065664 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 10065664 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 10065664 times.
✗ Branch 11 not taken.
|
40316416 | { return *gridVariables_; } |
234 | |||
235 | //! The jacobian matrix | ||
236 | JacobianMatrix& jacobian() | ||
237 | 8632 | { return *jacobian_; } | |
238 | |||
239 | //! The residual vector (rhs) | ||
240 | ResidualType& residual() | ||
241 |
2/4✓ Branch 3 taken 48 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 48 times.
✗ Branch 7 not taken.
|
4412 | { return *residual_; } |
242 | |||
243 | //! The solution of the previous time step | ||
244 | ✗ | const SolutionVector& prevSol() const | |
245 | ✗ | { return *prevSol_; } | |
246 | |||
247 | /*! | ||
248 | * \brief Create a local residual object (used by the local assembler) | ||
249 | */ | ||
250 | MultiStageFVLocalOperator<LocalResidual> localResidual() const | ||
251 |
3/6✓ Branch 1 taken 10065664 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 10065664 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 5052608 times.
✗ Branch 8 not taken.
|
25210816 | { return { LocalResidual{problem_.get(), nullptr} }; } |
252 | |||
253 | /*! | ||
254 | * \brief Update the grid variables | ||
255 | */ | ||
256 | ✗ | void updateGridVariables(const SolutionVector &cursol) | |
257 |
2/4✓ Branch 1 taken 2000 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2000 times.
✗ Branch 5 not taken.
|
4158 | { gridVariables().update(cursol); } |
258 | |||
259 | /*! | ||
260 | * \brief Reset the gridVariables | ||
261 | */ | ||
262 | ✗ | void resetTimeStep(const SolutionVector &cursol) | |
263 | { | ||
264 | ✗ | gridVariables().resetTimeStep(cursol); | |
265 | ✗ | this->clearStages(); | |
266 | ✗ | } | |
267 | |||
268 | 4096 | void clearStages() | |
269 | { | ||
270 | 4096 | spatialOperatorEvaluations_.clear(); | |
271 | 4096 | temporalOperatorEvaluations_.clear(); | |
272 |
2/2✓ Branch 0 taken 2048 times.
✓ Branch 1 taken 2048 times.
|
4096 | stageParams_.reset(); |
273 | 4096 | } | |
274 | |||
275 | template<class StageParams> | ||
276 | 2048 | void prepareStage(SolutionVector& x, StageParams params) | |
277 | { | ||
278 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2048 times.
|
2048 | stageParams_ = std::move(params); |
279 |
1/2✓ Branch 0 taken 2048 times.
✗ Branch 1 not taken.
|
2048 | const auto curStage = stageParams_->size() - 1; |
280 | |||
281 | // in the first stage, also assemble the residual | ||
282 | // at the previous time level (stage 0 residual) | ||
283 |
1/2✓ Branch 0 taken 2048 times.
✗ Branch 1 not taken.
|
2048 | if (curStage == 1) |
284 | { | ||
285 | // update time in variables? | ||
286 | 6144 | setProblemTime_(*problem_, stageParams_->timeAtStage(0)); | |
287 | |||
288 | 2048 | resetResidual_(); // residual resized and zero | |
289 | |||
290 | 2048 | assert(spatialOperatorEvaluations_.size() >= 0); | |
291 |
1/2✓ Branch 0 taken 2048 times.
✗ Branch 1 not taken.
|
2048 | if (spatialOperatorEvaluations_.size() == 0) |
292 | { | ||
293 | 4096 | spatialOperatorEvaluations_.push_back(*residual_); | |
294 | 4096 | temporalOperatorEvaluations_.push_back(*residual_); | |
295 | |||
296 | // assemble stage 0 residuals | ||
297 | 20100352 | assemble_([&](const auto& element) | |
298 | { | ||
299 |
1/5✓ Branch 1 taken 3072 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
12543008 | LocalAssembler localAssembler(*this, element, *prevSol_); |
300 |
1/4✓ Branch 1 taken 5018432 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
5018432 | localAssembler.localResidual().spatialWeight(1.0); |
301 |
1/4✓ Branch 1 taken 5018432 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
5018432 | localAssembler.localResidual().temporalWeight(1.0); |
302 |
3/12✓ Branch 1 taken 5018432 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 5018432 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 5018432 times.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
15055296 | localAssembler.assembleCurrentResidual(spatialOperatorEvaluations_.back(), temporalOperatorEvaluations_.back()); |
303 | }); | ||
304 | } | ||
305 | |||
306 | // we don't delete the first stage so it can be reused in a restarted | ||
307 | // time integration step. The evaluations are only deleted | ||
308 | // when explicitly requested by calling clearStages(). | ||
309 | // So if here the vector is non-empty, we don't need to evaluate again | ||
310 | // (this should only occur if we are restarting time integration, e.g. | ||
311 | // with a different time step size) | ||
312 | ✗ | else if (spatialOperatorEvaluations_.size() > 0) | |
313 | { | ||
314 | ✗ | updateGridVariables(x); | |
315 | ✗ | spatialOperatorEvaluations_.resize(1); | |
316 | ✗ | temporalOperatorEvaluations_.resize(1); | |
317 | } | ||
318 | } | ||
319 | |||
320 | // update time in variables? | ||
321 | 6144 | setProblemTime_(*problem_, stageParams_->timeAtStage(curStage)); | |
322 | |||
323 | 2048 | resetResidual_(); // residual resized and zero | |
324 | |||
325 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 2048 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2048 times.
|
4096 | if (spatialOperatorEvaluations_.size() != curStage) |
326 | ✗ | DUNE_THROW(Dune::InvalidStateException, | |
327 | "Invalid state. Maybe you forgot to call clearStages()"); | ||
328 | |||
329 | // allocate memory for this stage | ||
330 | 4096 | spatialOperatorEvaluations_.push_back(*residual_); | |
331 | 4096 | temporalOperatorEvaluations_.push_back(*residual_); | |
332 | 2048 | } | |
333 | |||
334 | //! TODO get rid of this (called by Newton but shouldn't be necessary) | ||
335 | ✗ | bool isStationaryProblem() const | |
336 | ✗ | { return false; } | |
337 | |||
338 | bool isImplicit() const | ||
339 | 20158208 | { return timeSteppingMethod_->implicit(); } | |
340 | |||
341 | private: | ||
342 | /*! | ||
343 | * \brief Resizes the jacobian and sets the jacobian' sparsity pattern. | ||
344 | */ | ||
345 | 8 | void setJacobianPattern_() | |
346 | { | ||
347 | // resize the jacobian and the residual | ||
348 | 8 | const auto numDofs = this->numDofs(); | |
349 | 16 | jacobian_->setSize(numDofs, numDofs); | |
350 | |||
351 | // create occupation pattern of the jacobian | ||
352 |
2/2✓ Branch 2 taken 6 times.
✓ Branch 3 taken 2 times.
|
16 | if (timeSteppingMethod_->implicit()) |
353 |
2/4✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 6 times.
✗ Branch 7 not taken.
|
12 | getJacobianPattern<true>(gridGeometry()).exportIdx(*jacobian_); |
354 | else | ||
355 |
2/4✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
|
4 | getJacobianPattern<false>(gridGeometry()).exportIdx(*jacobian_); |
356 | 8 | } | |
357 | |||
358 | //! Resizes the residual | ||
359 | 6174 | void setResidualSize_(ResidualType& res) | |
360 | 6264 | { res.resize(numDofs()); } | |
361 | |||
362 | //! Computes the colors | ||
363 | 8 | void maybeComputeColors_() | |
364 | { | ||
365 |
1/2✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
|
8 | if (enableMultithreading_) |
366 | 16 | elementSets_ = computeColoring(gridGeometry()).sets; | |
367 | 8 | } | |
368 | |||
369 | // reset the residual vector to 0.0 | ||
370 | 6254 | void resetResidual_() | |
371 | { | ||
372 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6254 times.
|
6254 | if(!residual_) |
373 | { | ||
374 | ✗ | residual_ = std::make_shared<ResidualType>(); | |
375 | ✗ | setResidualSize_(*residual_); | |
376 | } | ||
377 | |||
378 | 6254 | setResidualSize_(constrainedDofs_); | |
379 | |||
380 | 12508 | (*residual_) = 0.0; | |
381 | 6254 | constrainedDofs_ = 0.0; | |
382 | 6254 | } | |
383 | |||
384 | // reset the Jacobian matrix to 0.0 | ||
385 | template <class PartialReassembler = DefaultPartialReassembler> | ||
386 | 2158 | void resetJacobian_(const PartialReassembler *partialReassembler = nullptr) | |
387 | { | ||
388 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2158 times.
|
2158 | if(!jacobian_) |
389 | { | ||
390 | ✗ | jacobian_ = std::make_shared<JacobianMatrix>(); | |
391 | ✗ | jacobian_->setBuildMode(JacobianMatrix::random); | |
392 | ✗ | setJacobianPattern_(); | |
393 | } | ||
394 | |||
395 |
2/2✓ Branch 0 taken 2000 times.
✓ Branch 1 taken 158 times.
|
2158 | if (partialReassembler) |
396 | ✗ | partialReassembler->resetJacobian(*this); | |
397 | else | ||
398 | 4316 | *jacobian_ = 0.0; | |
399 | 2158 | } | |
400 | |||
401 | /*! | ||
402 | * \brief A method assembling something per element | ||
403 | * \note Handles exceptions for parallel runs | ||
404 | * \throws NumericalProblem on all processes if an exception is thrown during assembly | ||
405 | */ | ||
406 | template<typename AssembleElementFunc> | ||
407 | 8412 | void assemble_(AssembleElementFunc&& assembleElement) const | |
408 | { | ||
409 | // a state that will be checked on all processes | ||
410 | 8412 | bool succeeded = false; | |
411 | |||
412 | // try assembling using the local assembly function | ||
413 | try | ||
414 | { | ||
415 |
1/2✓ Branch 0 taken 4206 times.
✗ Branch 1 not taken.
|
8412 | if (enableMultithreading_) |
416 | { | ||
417 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 4206 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4206 times.
|
16824 | assert(elementSets_.size() > 0); |
418 | |||
419 | // make this element loop run in parallel | ||
420 | // for this we have to color the elements so that we don't get | ||
421 | // race conditions when writing into the global matrix | ||
422 | // each color can be assembled using multiple threads | ||
423 |
4/4✓ Branch 0 taken 37436 times.
✓ Branch 1 taken 4206 times.
✓ Branch 2 taken 37436 times.
✓ Branch 3 taken 4206 times.
|
166568 | for (const auto& elements : elementSets_) |
424 | { | ||
425 |
1/2✓ Branch 2 taken 37436 times.
✗ Branch 3 not taken.
|
20307952 | Dumux::parallelFor(elements.size(), [&](const std::size_t i) |
426 | { | ||
427 |
4/8✗ Branch 0 not taken.
✓ Branch 1 taken 5060672 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5060672 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 5018432 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 5018432 times.
|
20158208 | const auto element = gridView().grid().entity(elements[i]); |
428 | 10079104 | assembleElement(element); | |
429 | }); | ||
430 | } | ||
431 | } | ||
432 | else | ||
433 | ✗ | for (const auto& element : elements(gridView())) | |
434 | ✗ | assembleElement(element); | |
435 | |||
436 | // if we get here, everything worked well on this process | ||
437 | 8412 | succeeded = true; | |
438 | } | ||
439 | // throw exception if a problem occurred | ||
440 | ✗ | catch (NumericalProblem &e) | |
441 | { | ||
442 | ✗ | std::cout << "rank " << gridView().comm().rank() | |
443 | ✗ | << " caught an exception while assembling:" << e.what() | |
444 | ✗ | << "\n"; | |
445 | ✗ | succeeded = false; | |
446 | } | ||
447 | |||
448 | // make sure everything worked well on all processes | ||
449 |
4/8✗ Branch 0 not taken.
✓ Branch 1 taken 4206 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4206 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 4206 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 4206 times.
|
16824 | if (gridView().comm().size() > 1) |
450 | ✗ | succeeded = gridView().comm().min(succeeded); | |
451 | |||
452 | // if not succeeded rethrow the error on all processes | ||
453 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4206 times.
|
8412 | if (!succeeded) |
454 | ✗ | DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system"); | |
455 | 8412 | } | |
456 | |||
457 | // TODO make this nicer with a is_detected trait in a common location | ||
458 | template<class P> | ||
459 | ✗ | void setProblemTime_(const P& p, const Scalar t) | |
460 | 4096 | { setProblemTimeImpl_(p, t, 0); } | |
461 | |||
462 | template<class P> | ||
463 | auto setProblemTimeImpl_(const P& p, const Scalar t, int) -> decltype(p.setTime(0)) | ||
464 | { p.setTime(t); } | ||
465 | |||
466 | template<class P> | ||
467 | ✗ | void setProblemTimeImpl_(const P& p, const Scalar t, long) | |
468 | ✗ | {} | |
469 | |||
470 | std::shared_ptr<const Experimental::MultiStageMethod<Scalar>> timeSteppingMethod_; | ||
471 | std::vector<ResidualType> spatialOperatorEvaluations_; | ||
472 | std::vector<ResidualType> temporalOperatorEvaluations_; | ||
473 | ResidualType constrainedDofs_; | ||
474 | std::shared_ptr<const StageParams> stageParams_; | ||
475 | |||
476 | //! pointer to the problem to be solved | ||
477 | std::shared_ptr<const Problem> problem_; | ||
478 | |||
479 | //! the finite volume geometry of the grid | ||
480 | std::shared_ptr<const GridGeometry> gridGeometry_; | ||
481 | |||
482 | //! the variables container for the grid | ||
483 | std::shared_ptr<GridVariables> gridVariables_; | ||
484 | |||
485 | //! an observing pointer to the previous solution for instationary problems | ||
486 | const SolutionVector* prevSol_ = nullptr; | ||
487 | |||
488 | //! shared pointers to the jacobian matrix and residual | ||
489 | std::shared_ptr<JacobianMatrix> jacobian_; | ||
490 | std::shared_ptr<ResidualType> residual_; | ||
491 | |||
492 | //! element sets for parallel assembly | ||
493 | bool enableMultithreading_ = false; | ||
494 | std::deque<std::vector<ElementSeed>> elementSets_; | ||
495 | }; | ||
496 | |||
497 | } // namespace Dumux | ||
498 | |||
499 | #endif | ||
500 |