1 #ifndef DUNE_PDELAB_GRIDOPERATOR_GRIDOPERATOR_HH 2 #define DUNE_PDELAB_GRIDOPERATOR_GRIDOPERATOR_HH 6 #include <dune/common/hybridutilities.hh> 30 template<
typename GFSU,
typename GFSV,
typename LOP,
31 typename MB,
typename DF,
typename RF,
typename JF,
56 GFSU::Traits::EntitySet::Partitions::partitionIterator() == InteriorBorder_Partition
60 typedef typename std::conditional<
61 GFSU::Traits::EntitySet::Partitions::partitionIterator() == InteriorBorder_Partition,
70 template <
typename MFT>
76 GridOperator(
const GFSU & gfsu_,
const CU & cu_,
const GFSV & gfsv_,
const CV & cv_, LOP & lop_,
const MB& mb_ = MB())
77 : global_assembler(gfsu_,gfsv_,cu_,cv_)
78 , dof_exchanger(
std::make_shared<BorderDOFExchanger>(*this))
79 , local_assembler(lop_, cu_, cv_,dof_exchanger)
84 GridOperator(
const GFSU & gfsu_,
const GFSV & gfsv_, LOP & lop_,
const MB& mb_ = MB())
85 : global_assembler(gfsu_,gfsv_)
86 , dof_exchanger(
std::make_shared<BorderDOFExchanger>(*this))
87 , local_assembler(lop_,dof_exchanger)
94 return global_assembler.trialGridFunctionSpace();
100 return global_assembler.testGridFunctionSpace();
117 const Assembler &
assembler()
const {
return global_assembler; }
124 template <
typename Gr
idOperatorTuple>
128 :
index(0), size(
std::tuple_size<GridOperatorTuple>::
value) {}
130 template <
typename T>
132 elem.localAssembler().preProcessing(
index == 0);
133 elem.localAssembler().postProcessing(
index == size-1);
144 template<
typename Gr
idOperatorTuple>
147 SetupGridOperator<GridOperatorTuple> setup_visitor;
148 Hybrid::forEach(tuple,
149 [&](
auto &el) { setup_visitor.visit(el); });
153 template<
typename F,
typename X>
167 PatternEngine & pattern_engine = local_assembler.localPatternAssemblerEngine(p);
168 global_assembler.assemble(pattern_engine);
175 ResidualEngine & residual_engine = local_assembler.localResidualAssemblerEngine(r,x);
176 global_assembler.assemble(residual_engine);
183 JacobianEngine & jacobian_engine = local_assembler.localJacobianAssemblerEngine(a,x);
184 global_assembler.assemble(jacobian_engine);
191 JacobianApplyEngine & jacobian_apply_engine = local_assembler.localJacobianApplyAssemblerEngine(r,z);
192 global_assembler.assemble(jacobian_apply_engine);
198 global_assembler.assemble(local_assembler.localNonlinearJacobianApplyAssemblerEngine(r,x,z));
204 dof_exchanger->accumulateBorderEntries(*
this,a);
210 dof_exchanger->update(*
this);
220 Assembler global_assembler;
221 shared_ptr<BorderDOFExchanger> dof_exchanger;
230 #endif // DUNE_PDELAB_GRIDOPERATOR_GRIDOPERATOR_HH The assembler for standard DUNE grid.
Definition: default/assembler.hh:22
void visit(T &elem)
Definition: gridoperator.hh:131
Dune::PDELab::Backend::Matrix< MBE, Domain, Range, field_type > Jacobian
The type of the jacobian.
Definition: gridoperator.hh:47
GridOperator(const GFSU &gfsu_, const CU &cu_, const GFSV &gfsv_, const CV &cv_, LOP &lop_, const MB &mb_=MB())
Constructor for non trivial constraints.
Definition: gridoperator.hh:76
int index
Definition: gridoperator.hh:137
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:106
Definition: constraintstransformation.hh:111
void make_consistent(Jacobian &a) const
Definition: gridoperator.hh:202
void residual(const Domain &x, Range &r) const
Assemble residual.
Definition: gridoperator.hh:172
MBE ::template Pattern< Jacobian, GFS, CGGFS > Pattern
The sparsity pattern container for the jacobian matrix.
Definition: gridoperator.hh:50
DefaultAssembler< GFSU, GFSV, CU, CV > Assembler
The global assembler type.
Definition: gridoperator.hh:40
SetupGridOperator()
Definition: gridoperator.hh:127
const Traits::MatrixBackend & matrixBackend() const
Get the matrix backend for this grid operator.
Definition: gridoperator.hh:214
Dune::PDELab::Backend::Vector< CGGFS, field_type > Domain
The type of the domain (solution).
Definition: gridoperator.hh:43
Definition: gridoperator.hh:71
std::conditional< GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_Partition, NonOverlappingBorderDOFExchanger< GridOperator >, OverlappingBorderDOFExchanger< GridOperator > >::type BorderDOFExchanger
Definition: gridoperator.hh:64
Dune::PDELab::GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, Assembler, LocalAssembler > Traits
The grid operator traits.
Definition: gridoperator.hh:68
typename impl::BackendMatrixSelector< Backend, VU, VV, E >::Type Matrix
alias of the return type of BackendMatrixSelector
Definition: backend/interface.hh:127
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
void update()
Definition: gridoperator.hh:207
const int size
Definition: gridoperator.hh:138
Assembler & assembler()
Definition: gridoperator.hh:115
static void setupGridOperators(GridOperatorTuple tuple)
Definition: gridoperator.hh:145
Standard grid operator implementation.
Definition: gridoperator.hh:35
MB MatrixBackend
The matrix backend of the grid operator.
Definition: gridoperatorutilities.hh:51
void nonlinear_jacobian_apply(const Domain &x, const Domain &z, Range &r) const
Apply jacobian matrix without explicitly assembling it.
Definition: gridoperator.hh:196
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
GFSU::Traits::SizeType globalSizeU() const
Get dimension of space u.
Definition: gridoperator.hh:104
const Assembler & assembler() const
Definition: gridoperator.hh:117
Dune::PDELab::Backend::Vector< GFS, field_type > Range
The type of the range (residual).
Definition: gridoperator.hh:45
void interpolate(const X &xold, F &f, X &x) const
Interpolate the constrained dofs from given function.
Definition: gridoperator.hh:154
void jacobian(const Domain &x, Jacobian &a) const
Assembler jacobian.
Definition: gridoperator.hh:180
void fill_pattern(Pattern &p) const
Fill pattern of jacobian matrix.
Definition: gridoperator.hh:164
Traits::Jacobian Type
Definition: gridoperator.hh:72
GFSV::Traits::SizeType globalSizeV() const
Get dimension of space v.
Definition: gridoperator.hh:110
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: gridoperator.hh:92
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
GridOperator(const GFSU &gfsu_, const GFSV &gfsv_, LOP &lop_, const MB &mb_=MB())
Constructor for empty constraints.
Definition: gridoperator.hh:84
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: gridoperator.hh:98
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:989
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:204
Helper class for adding up matrix entries on border.
Definition: borderdofexchanger.hh:67
std::size_t index
Definition: interpolate.hh:118
DefaultLocalAssembler< GridOperator, LOP, GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_Partition > LocalAssembler
The local assembler type.
Definition: gridoperator.hh:57
const P & p
Definition: constraints.hh:147
LocalAssembler & localAssembler() const
Definition: gridoperator.hh:119
void jacobian_apply(const Domain &z, Range &r) const
Apply jacobian matrix without explicitly assembling it.
Definition: gridoperator.hh:188
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:33
Definition: gridoperator.hh:125
The local assembler for DUNE grids.
Definition: default/localassembler.hh:32
Definition: borderdofexchanger.hh:577