12#ifndef DUMUX_FV_ASSEMBLER_HH
13#define DUMUX_FV_ASSEMBLER_HH
20#include <dune/istl/matrixindexset.hh>
45template<
class DiscretizationMethod>
51 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
58 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
65 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
72 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
76template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
79>::template type<TypeTag, Impl, diffMethod, isImplicit>;
92template<
class TypeTag, DiffMethod diffMethod,
bool isImplicit = true>
96 using GridView =
typename GridGeo::GridView;
98 using Element =
typename GridView::template Codim<0>::Entity;
99 using ElementSeed =
typename GridView::Grid::template Codim<0>::EntitySeed;
130 , isStationaryProblem_(true)
132 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
138 maybeComputeColors_();
149 std::shared_ptr<const TimeLoop> timeLoop,
154 , timeLoop_(timeLoop)
156 , isStationaryProblem_(!timeLoop)
163 maybeComputeColors_();
170 template<
class PartialReassembler = DefaultPartialReassembler>
173 checkAssemblerState_();
174 resetJacobian_(partialReassembler);
177 assemble_([&](
const Element& element)
179 LocalAssembler localAssembler(*
this, element, curSol);
180 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
183 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
191 checkAssemblerState_();
194 assemble_([&](
const Element& element)
196 LocalAssembler localAssembler(*
this, element, curSol);
197 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
211 checkAssemblerState_();
213 assemble_([&](
const Element& element)
215 LocalAssembler localAssembler(*
this, element, curSol);
216 localAssembler.assembleResidual(r);
226 std::shared_ptr<ResidualType> r)
232 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
233 jacobian_->setBuildMode(JacobianMatrix::random);
234 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
235 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
238 setJacobianPattern_();
247 jacobian_ = std::make_shared<JacobianMatrix>();
248 jacobian_->setBuildMode(JacobianMatrix::random);
249 residual_ = std::make_shared<ResidualType>();
252 setJacobianPattern_();
261 setJacobianPattern_();
262 maybeComputeColors_();
267 {
return gridGeometry_->numDofs(); }
271 {
return *problem_; }
275 {
return *gridGeometry_; }
283 {
return *gridVariables_; }
287 {
return *gridVariables_; }
291 {
return *jacobian_; }
295 {
return *residual_; }
299 {
return *prevSol_; }
306 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
319 {
return isStationaryProblem_; }
325 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
347 void setJacobianPattern_()
357 occupationPattern.exportIdx(*jacobian_);
361 void setResidualSize_()
362 { residual_->resize(
numDofs()); }
365 void maybeComputeColors_()
367 if (enableMultithreading_)
372 void resetResidual_()
376 residual_ = std::make_shared<ResidualType>();
384 template <
class PartialReassembler = DefaultPartialReassembler>
385 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
389 jacobian_ = std::make_shared<JacobianMatrix>();
390 jacobian_->setBuildMode(JacobianMatrix::random);
391 setJacobianPattern_();
394 if (partialReassembler)
395 partialReassembler->resetJacobian(*
this);
401 void checkAssemblerState_()
const
403 if (!isStationaryProblem_ && !prevSol_)
404 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
412 template<
typename AssembleElementFunc>
413 void assemble_(AssembleElementFunc&& assembleElement)
const
416 bool succeeded =
false;
421 if (enableMultithreading_)
423 assert(elementSets_.size() > 0);
429 for (
const auto& elements : elementSets_)
433 const auto element = gridView().grid().entity(elements[i]);
434 assembleElement(element);
439 for (
const auto& element : elements(
gridView()))
446 catch (NumericalProblem &e)
448 std::cout <<
"rank " <<
gridView().comm().rank()
449 <<
" caught an exception while assembling:" << e.what()
456 succeeded =
gridView().comm().min(succeeded);
460 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
470 if (m.first < m.second)
473 res[m.first] += res[m.second];
474 const auto end = jac[m.second].end();
475 for (
auto it = jac[m.second].begin(); it != end; ++it)
476 jac[m.first][it.index()] += (*it);
479 res[m.second] = curSol[m.second] - curSol[m.first];
482 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
484 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
485 matrixBlock[eIdx][eIdx] = diagValue;
488 for (
auto it = jac[m.second].begin(); it != end; ++it)
490 auto& matrixBlock = *it;
493 assert(matrixBlock.N() == matrixBlock.M());
494 if(it.index() == m.second)
495 setMatrixBlock(matrixBlock, 1.0);
497 if(it.index() == m.first)
498 setMatrixBlock(matrixBlock, -1.0);
507 std::shared_ptr<const Problem> problem_;
510 std::shared_ptr<const GridGeometry> gridGeometry_;
513 std::shared_ptr<GridVariables> gridVariables_;
516 std::shared_ptr<const TimeLoop> timeLoop_;
522 bool isStationaryProblem_;
525 std::shared_ptr<JacobianMatrix> jacobian_;
526 std::shared_ptr<ResidualType> residual_;
529 bool enableMultithreading_ =
false;
530 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:302
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition assembly/fvassembler.hh:94
ResidualType & residual()
The residual vector (rhs)
Definition assembly/fvassembler.hh:294
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition assembly/fvassembler.hh:108
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition assembly/fvassembler.hh:189
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition assembly/fvassembler.hh:258
const GridView & gridView() const
The gridview.
Definition assembly/fvassembler.hh:278
const Problem & problem() const
The problem.
Definition assembly/fvassembler.hh:270
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition assembly/fvassembler.hh:111
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition assembly/fvassembler.hh:245
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition assembly/fvassembler.hh:318
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition assembly/fvassembler.hh:109
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition assembly/fvassembler.hh:324
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition assembly/fvassembler.hh:110
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:146
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition assembly/fvassembler.hh:266
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:225
GridVariables & gridVariables()
The global grid variables.
Definition assembly/fvassembler.hh:282
JacobianMatrix & jacobian()
The jacobian matrix.
Definition assembly/fvassembler.hh:290
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition assembly/fvassembler.hh:209
GetPropType< TypeTag, Properties::Problem > Problem
Definition assembly/fvassembler.hh:116
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition assembly/fvassembler.hh:305
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition assembly/fvassembler.hh:202
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition assembly/fvassembler.hh:113
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:312
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition assembly/fvassembler.hh:274
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition assembly/fvassembler.hh:338
const GridVariables & gridVariables() const
The global grid variables.
Definition assembly/fvassembler.hh:286
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:171
GridGeo GridGeometry
Definition assembly/fvassembler.hh:115
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:123
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition assembly/fvassembler.hh:298
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition assembly/fvassembler.hh:330
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
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for cell-centered methods.
Definition jacobianpattern.hh:28
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
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
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:77
constexpr bool hasPeriodicDofMap()
Definition periodic.hh:24
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.
Type traits to detect periodicity support.
Definition assembly/fvassembler.hh:46
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.