12#ifndef DUMUX_FV_ASSEMBLER_HH
13#define DUMUX_FV_ASSEMBLER_HH
20#include <dune/istl/matrixindexset.hh>
44template<
class DiscretizationMethod>
50 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
57 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
64 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
71 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
75template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
78>::template type<TypeTag, Impl, diffMethod, isImplicit>;
91template<
class TypeTag, DiffMethod diffMethod,
bool isImplicit = true>
95 using GridView =
typename GridGeo::GridView;
97 using Element =
typename GridView::template Codim<0>::Entity;
98 using ElementSeed =
typename GridView::Grid::template Codim<0>::EntitySeed;
129 , isStationaryProblem_(true)
131 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
135 && getParam<bool>(
"Assembly.Multithreading",
true);
137 maybeComputeColors_();
148 std::shared_ptr<const TimeLoop> timeLoop,
153 , timeLoop_(timeLoop)
155 , isStationaryProblem_(!timeLoop)
160 && getParam<bool>(
"Assembly.Multithreading",
true);
162 maybeComputeColors_();
169 template<
class PartialReassembler = DefaultPartialReassembler>
172 checkAssemblerState_();
173 resetJacobian_(partialReassembler);
176 assemble_([&](
const Element& element)
178 LocalAssembler localAssembler(*
this, element, curSol);
179 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
182 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
190 checkAssemblerState_();
193 assemble_([&](
const Element& element)
195 LocalAssembler localAssembler(*
this, element, curSol);
196 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
210 checkAssemblerState_();
212 assemble_([&](
const Element& element)
214 LocalAssembler localAssembler(*
this, element, curSol);
215 localAssembler.assembleResidual(r);
225 std::shared_ptr<ResidualType> r)
231 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
232 jacobian_->setBuildMode(JacobianMatrix::random);
233 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
234 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
237 setJacobianPattern_();
246 jacobian_ = std::make_shared<JacobianMatrix>();
247 jacobian_->setBuildMode(JacobianMatrix::random);
248 residual_ = std::make_shared<ResidualType>();
251 setJacobianPattern_();
260 setJacobianPattern_();
261 maybeComputeColors_();
266 {
return gridGeometry_->numDofs(); }
270 {
return *problem_; }
274 {
return *gridGeometry_; }
282 {
return *gridVariables_; }
286 {
return *gridVariables_; }
290 {
return *jacobian_; }
294 {
return *residual_; }
298 {
return *prevSol_; }
305 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
318 {
return isStationaryProblem_; }
324 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
346 void setJacobianPattern_()
353 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
356 occupationPattern.exportIdx(*jacobian_);
360 void setResidualSize_()
361 { residual_->resize(
numDofs()); }
364 void maybeComputeColors_()
366 if (enableMultithreading_)
371 void resetResidual_()
375 residual_ = std::make_shared<ResidualType>();
383 template <
class PartialReassembler = DefaultPartialReassembler>
384 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
388 jacobian_ = std::make_shared<JacobianMatrix>();
389 jacobian_->setBuildMode(JacobianMatrix::random);
390 setJacobianPattern_();
393 if (partialReassembler)
394 partialReassembler->resetJacobian(*
this);
400 void checkAssemblerState_()
const
402 if (!isStationaryProblem_ && !prevSol_)
403 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
411 template<
typename AssembleElementFunc>
412 void assemble_(AssembleElementFunc&& assembleElement)
const
415 bool succeeded =
false;
420 if (enableMultithreading_)
422 assert(elementSets_.size() > 0);
428 for (
const auto& elements : elementSets_)
432 const auto element = gridView().grid().entity(elements[i]);
433 assembleElement(element);
438 for (
const auto& element : elements(
gridView()))
445 catch (NumericalProblem &e)
447 std::cout <<
"rank " <<
gridView().comm().rank()
448 <<
" caught an exception while assembling:" << e.what()
455 succeeded =
gridView().comm().min(succeeded);
459 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
462 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethods::box, void>
467 if (m.first < m.second)
470 res[m.first] += res[m.second];
471 const auto end = jac[m.second].end();
472 for (
auto it = jac[m.second].begin(); it != end; ++it)
473 jac[m.first][it.index()] += (*it);
476 res[m.second] = curSol[m.second] - curSol[m.first];
477 for (
auto it = jac[m.second].begin(); it != end; ++it)
478 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
483 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethods::box, void>
487 std::shared_ptr<const Problem> problem_;
490 std::shared_ptr<const GridGeometry> gridGeometry_;
493 std::shared_ptr<GridVariables> gridVariables_;
496 std::shared_ptr<const TimeLoop> timeLoop_;
502 bool isStationaryProblem_;
505 std::shared_ptr<JacobianMatrix> jacobian_;
506 std::shared_ptr<ResidualType> residual_;
509 bool enableMultithreading_ =
false;
510 std::deque<std::vector<ElementSeed>> elementSets_;
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An assembler for Jacobian and residual contribution per element (CVFE methods)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition assembly/cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition assembly/cvfelocalassembler.hh:299
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition assembly/fvassembler.hh:93
ResidualType & residual()
The residual vector (rhs)
Definition assembly/fvassembler.hh:293
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition assembly/fvassembler.hh:107
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition assembly/fvassembler.hh:188
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition assembly/fvassembler.hh:257
const GridView & gridView() const
The gridview.
Definition assembly/fvassembler.hh:277
const Problem & problem() const
The problem.
Definition assembly/fvassembler.hh:269
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition assembly/fvassembler.hh:110
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition assembly/fvassembler.hh:244
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition assembly/fvassembler.hh:317
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition assembly/fvassembler.hh:108
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition assembly/fvassembler.hh:323
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition assembly/fvassembler.hh:109
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition assembly/fvassembler.hh:145
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition assembly/fvassembler.hh:265
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< ResidualType > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition assembly/fvassembler.hh:224
GridVariables & gridVariables()
The global grid variables.
Definition assembly/fvassembler.hh:281
JacobianMatrix & jacobian()
The jacobian matrix.
Definition assembly/fvassembler.hh:289
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition assembly/fvassembler.hh:208
GetPropType< TypeTag, Properties::Problem > Problem
Definition assembly/fvassembler.hh:115
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition assembly/fvassembler.hh:304
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition assembly/fvassembler.hh:201
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition assembly/fvassembler.hh:112
void setPreviousSolution(const SolutionVector &u)
Sets the solution from which to start the time integration. Has to be called prior to assembly for ti...
Definition assembly/fvassembler.hh:311
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition assembly/fvassembler.hh:273
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition assembly/fvassembler.hh:337
const GridVariables & gridVariables() const
The global grid variables.
Definition assembly/fvassembler.hh:285
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition assembly/fvassembler.hh:170
GridGeo GridGeometry
Definition assembly/fvassembler.hh:114
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables)
The constructor for stationary problems.
Definition assembly/fvassembler.hh:122
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition assembly/fvassembler.hh:297
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition assembly/fvassembler.hh:329
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition fclocalassembler.hh:284
detects which entries in the Jacobian have to be recomputed
Definition partialreassembler.hh:420
Manages the handling of time dependent problems.
Definition common/timeloop.hh:84
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper to extract native Dune vector types from particular Dumux types.
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition multithreading.hh:45
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
Distance implementation details.
Definition cvfelocalresidual.hh:25
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition assembly/fvassembler.hh:76
constexpr Box box
Definition method.hh:147
bool supportsMultithreading(const GridView &gridView)
Definition gridcapabilities.hh:74
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition coloring.hh:239
Parallel for loop (multithreading)
Provides a helper class for nonoverlapping decomposition.
Definition assembly/fvassembler.hh:45
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition dunevectors.hh:57
Traits specifying if a given discretization tag supports coloring.
Definition coloring.hh:292
Type traits to be used with vector types.