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);
220 [[deprecated(
"Use the linear solver's norm. Will be deleted after 3.7")]]
224 static bool warningIssued =
false;
230 using DM =
typename GridGeometry::VertexMapper;
233 PVHelper vectorHelper(
gridView(), gridGeometry_->vertexMapper());
235 vectorHelper.makeNonOverlappingConsistent(
residual);
238 else if (!warningIssued)
241 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
242 <<
"overlapping entities multiple times. Please use the norm\n"
243 <<
"function provided by a linear solver instead." << std::endl;
245 warningIssued =
true;
251 result2 =
gridView().comm().sum(result2);
254 return sqrt(result2);
258 [[deprecated(
"Use assembleResidual and the linear solver's norm. Will be deleted after 3.7")]]
272 std::shared_ptr<ResidualType> r)
278 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
279 jacobian_->setBuildMode(JacobianMatrix::random);
280 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
281 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
284 setJacobianPattern_();
293 jacobian_ = std::make_shared<JacobianMatrix>();
294 jacobian_->setBuildMode(JacobianMatrix::random);
295 residual_ = std::make_shared<ResidualType>();
298 setJacobianPattern_();
307 setJacobianPattern_();
308 maybeComputeColors_();
313 {
return gridGeometry_->numDofs(); }
317 {
return *problem_; }
321 {
return *gridGeometry_; }
329 {
return *gridVariables_; }
333 {
return *gridVariables_; }
337 {
return *jacobian_; }
341 {
return *residual_; }
345 {
return *prevSol_; }
352 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
365 {
return isStationaryProblem_; }
371 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
393 void setJacobianPattern_()
400 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
403 occupationPattern.exportIdx(*jacobian_);
407 void setResidualSize_()
408 { residual_->resize(
numDofs()); }
411 void maybeComputeColors_()
413 if (enableMultithreading_)
418 void resetResidual_()
422 residual_ = std::make_shared<ResidualType>();
430 template <
class PartialReassembler = DefaultPartialReassembler>
431 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
435 jacobian_ = std::make_shared<JacobianMatrix>();
436 jacobian_->setBuildMode(JacobianMatrix::random);
437 setJacobianPattern_();
440 if (partialReassembler)
441 partialReassembler->resetJacobian(*
this);
447 void checkAssemblerState_()
const
449 if (!isStationaryProblem_ && !prevSol_)
450 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
458 template<
typename AssembleElementFunc>
459 void assemble_(AssembleElementFunc&& assembleElement)
const
462 bool succeeded =
false;
467 if (enableMultithreading_)
469 assert(elementSets_.size() > 0);
475 for (
const auto& elements : elementSets_)
479 const auto element = gridView().grid().entity(elements[i]);
480 assembleElement(element);
485 for (
const auto& element : elements(
gridView()))
486 assembleElement(element);
492 catch (NumericalProblem &e)
494 std::cout <<
"rank " <<
gridView().comm().rank()
495 <<
" caught an exception while assembling:" << e.what()
502 succeeded =
gridView().comm().min(succeeded);
506 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
509 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethods::box, void>
514 if (m.first < m.second)
517 res[m.first] += res[m.second];
518 const auto end = jac[m.second].end();
519 for (
auto it = jac[m.second].begin(); it != end; ++it)
520 jac[m.first][it.index()] += (*it);
523 res[m.second] = curSol[m.second] - curSol[m.first];
524 for (
auto it = jac[m.second].begin(); it != end; ++it)
525 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
530 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethods::box, void>
534 std::shared_ptr<const Problem> problem_;
537 std::shared_ptr<const GridGeometry> gridGeometry_;
540 std::shared_ptr<GridVariables> gridVariables_;
543 std::shared_ptr<const TimeLoop> timeLoop_;
549 bool isStationaryProblem_;
552 std::shared_ptr<JacobianMatrix> jacobian_;
553 std::shared_ptr<ResidualType> residual_;
556 bool enableMultithreading_ =
false;
557 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 (cell-centered methods)
Definition: cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: 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:340
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:304
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:324
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:316
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:291
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:364
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:370
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
Scalar normOfResidual(ResidualType &residual) const
compute a residual's vector norm (this is a temporary interface introduced during the deprecation per...
Definition: assembly/fvassembler.hh:221
Scalar residualNorm(const SolutionVector &curSol) const
compute the residual and return it's vector norm
Definition: assembly/fvassembler.hh:259
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:312
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:271
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:328
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:336
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:351
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:358
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:320
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:384
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:332
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:344
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:376
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:284
Definition: parallelhelpers.hh:476
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:56
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
An assembler for Jacobian and residual contribution per element (CVFE methods)
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:267
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:78
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:59
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:292
Type traits to be used with vector types.