dune-pdelab  2.5-dev
fastdg.hh
Go to the documentation of this file.
1 // -*- tab-width: 2; indent-tabs-mode: nil -*-
2 // vi: set et ts=2 sw=2 sts=2:
3 #ifndef DUNE_PDELAB_GRIDOPERATOR_FASTDG_HH
4 #define DUNE_PDELAB_GRIDOPERATOR_FASTDG_HH
5 
6 #include <dune/common/tupleutility.hh>
7 
13 
14 namespace Dune {
15  namespace PDELab {
16 
29  template<typename GFSU, typename GFSV, typename LOP,
30  typename MB, typename DF, typename RF, typename JF,
34  {
35  public:
36 
39 
46 
48  typedef typename MB::template Pattern<Jacobian,GFSV,GFSU> Pattern;
49 
51  typedef FastDGLocalAssembler<
53  LOP,
54  GFSU::Traits::EntitySet::Partitions::partitionIterator() == InteriorBorder_Partition
56 
57  // Fix this as soon as the default Partitions are constexpr
58  typedef typename std::conditional<
59  GFSU::Traits::EntitySet::Partitions::partitionIterator() == InteriorBorder_Partition,
63 
66  <GFSU,GFSV,MB,DF,RF,JF,CU,CV,Assembler,LocalAssembler> Traits;
67 
68  template <typename MFT>
70  typedef typename Traits::Jacobian Type;
71  };
72 
74  FastDGGridOperator(const GFSU & gfsu_, const CU & cu_, const GFSV & gfsv_, const CV & cv_, LOP & lop_, const MB& mb_ = MB())
75  : global_assembler(gfsu_,gfsv_,cu_,cv_)
76  , dof_exchanger(std::make_shared<BorderDOFExchanger>(*this))
77  , local_assembler(lop_, cu_, cv_,dof_exchanger)
78  , backend(mb_)
79  , lop(lop_)
80  , cu(cu_)
81  , cv(cv_)
82  {}
83 
85  FastDGGridOperator(const GFSU & gfsu_, const GFSV & gfsv_, LOP & lop_, const MB& mb_ = MB())
86  : global_assembler(gfsu_,gfsv_)
87  , dof_exchanger(std::make_shared<BorderDOFExchanger>(*this))
88  , local_assembler(lop_,dof_exchanger)
89  , backend(mb_)
90  , lop(lop_)
91  , cu(*std::make_shared<Dune::PDELab::EmptyTransformation>())
92  , cv(*std::make_shared<Dune::PDELab::EmptyTransformation>())
93  {}
94 
96  const GFSU& trialGridFunctionSpace() const
97  {
98  return global_assembler.trialGridFunctionSpace();
99  }
100 
102  const GFSV& testGridFunctionSpace() const
103  {
104  return global_assembler.testGridFunctionSpace();
105  }
106 
108  {
109  return cu;
110  }
112  {
113  return cu;
114  }
115 
117  {
118  return cv;
119  }
121  {
122  return cv;
123  }
124 
126  {
127  return lop;
128  }
129 
130  const LOP & localOperator() const
131  {
132  return lop;
133  }
134 
136  typename GFSU::Traits::SizeType globalSizeU () const
137  {
138  return trialGridFunctionSpace().globalSize();
139  }
140 
142  typename GFSV::Traits::SizeType globalSizeV () const
143  {
144  return testGridFunctionSpace().globalSize();
145  }
146 
147  Assembler & assembler() { return global_assembler; }
148 
149  const Assembler & assembler() const { return global_assembler; }
150 
151  LocalAssembler & localAssembler() const { return local_assembler; }
152 
153 
156  template <typename GridOperatorTuple>
158  {
160  : index(0), size(std::tuple_size<GridOperatorTuple>::value) {}
161 
162  template <typename T>
163  void visit(T& elem) {
164  elem.localAssembler().preProcessing(index == 0);
165  elem.localAssembler().postProcessing(index == size-1);
166  ++index;
167  }
168 
169  int index;
170  const int size;
171  };
172 
176  template<typename GridOperatorTuple>
177  static void setupGridOperators(GridOperatorTuple tuple)
178  {
180  Hybrid::forEach(tuple, [&](auto &el) { setup_visitor.visit(el); });
181  }
182 
184  template<typename F, typename X>
185  void interpolate (const X& xold, F& f, X& x) const
186  {
187  // Interpolate f into grid function space and set corresponding coefficients
188  Dune::PDELab::interpolate(f,global_assembler.trialGridFunctionSpace(),x);
189 
190  // Copy non-constrained dofs from old time step
192  }
193 
195  void fill_pattern(Pattern & p) const
196  {
197  typedef typename LocalAssembler::LocalPatternAssemblerEngine PatternEngine;
198  global_assembler.assemble
199  ([&p](LocalAssembler &la) -> PatternEngine&
200  {
201  return la.localPatternAssemblerEngine(p);
202  },
203  local_assembler
204  );
205  }
206 
208  void residual(const Domain & x, Range & r) const
209  {
210  typedef typename LocalAssembler::LocalResidualAssemblerEngine ResidualEngine;
211  global_assembler.assemble
212  ([&r,&x](LocalAssembler &la) -> ResidualEngine&
213  {
214  return la.localResidualAssemblerEngine(r,x);
215  },
216  local_assembler
217  );
218  }
219 
221  void jacobian(const Domain & x, Jacobian & a) const
222  {
223  typedef typename LocalAssembler::LocalJacobianAssemblerEngine JacobianEngine;
224  global_assembler.assemble
225  ([&a,&x](LocalAssembler &la) -> JacobianEngine&
226  {
227  return la.localJacobianAssemblerEngine(a,x);
228  },
229  local_assembler
230  );
231  }
232 
234  void jacobian_apply(const Domain & x, Range & r) const
235  {
236  typedef typename LocalAssembler::LocalJacobianApplyAssemblerEngine JacobianApplyEngine;
237  global_assembler.assemble
238  ([&r,&x](LocalAssembler &la) -> JacobianApplyEngine&
239  {
240  return la.localJacobianApplyAssemblerEngine(r,x);
241  },
242  local_assembler
243  );
244  }
245 
247  void nonlinear_jacobian_apply(const Domain & x, const Domain & z, Range & r) const
248  {
249  typedef typename LocalAssembler::LocalNonlinearJacobianApplyAssemblerEngine NonlinearJacobianApplyEngine;
250  global_assembler.assemble
251  ([&r,&x,&z](LocalAssembler &la) -> NonlinearJacobianApplyEngine&
252  {
254  },
255  local_assembler
256  );
257  }
258 
259  void make_consistent(Jacobian& a) const
260  {
261  dof_exchanger->accumulateBorderEntries(*this,a);
262  }
263 
264  void update()
265  {
266  // the DOF exchanger has matrix information, so we need to update it
267  dof_exchanger->update(*this);
268  }
269 
271  const typename Traits::MatrixBackend& matrixBackend() const
272  {
273  return backend;
274  }
275 
276  private:
277  Assembler global_assembler;
278  shared_ptr<BorderDOFExchanger> dof_exchanger;
279 
280  mutable LocalAssembler local_assembler;
281  MB backend;
282 
283  const LOP& lop;
284  const CU& cu;
285  const CV& cv;
286 
287  }; // end class FastDGGridOperator
288  } // end namespace PDELab
289 } // end namespace Dune
290 #endif
Dune::PDELab::FastDGGridOperator::fill_pattern
void fill_pattern(Pattern &p) const
Fill pattern of jacobian matrix.
Definition: fastdg.hh:195
Dune::PDELab::FastDGGridOperator::MatrixContainer
Definition: fastdg.hh:69
Dune::PDELab::FastDGGridOperator::Domain
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: fastdg.hh:41
Dune::PDELab::FastDGGridOperator::Traits
Dune::PDELab::GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, Assembler, LocalAssembler > Traits
The grid operator traits.
Definition: fastdg.hh:66
Dune::PDELab::copy_nonconstrained_dofs
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:989
Dune::PDELab::FastDGGridOperator::testGridFunctionSpaceConstraints
const CV & testGridFunctionSpaceConstraints() const
Definition: fastdg.hh:120
Dune::PDELab::FastDGGridOperator::BorderDOFExchanger
std::conditional< GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_Partition, NonOverlappingBorderDOFExchanger< FastDGGridOperator >, OverlappingBorderDOFExchanger< FastDGGridOperator > >::type BorderDOFExchanger
Definition: fastdg.hh:62
Dune::PDELab::FastDGGridOperator::SetupGridOperator::size
const int size
Definition: fastdg.hh:170
assembler.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine< FastDGLocalAssembler >
Dune::PDELab::FastDGGridOperator::localOperator
LOP & localOperator()
Definition: fastdg.hh:125
Dune::PDELab::FastDGLocalAssembler
The local assembler for DUNE grids.
Definition: fastdg/localassembler.hh:38
Dune::PDELab::FastDGLocalAssembler::localPatternAssemblerEngine
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: fastdg/localassembler.hh:157
Dune::PDELab::FastDGGridOperator::localOperator
const LOP & localOperator() const
Definition: fastdg.hh:130
Dune::PDELab::FastDGGridOperator::trialGridFunctionSpace
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: fastdg.hh:96
Dune::PDELab::FastDGGridOperator::SetupGridOperator::index
int index
Definition: fastdg.hh:169
Dune::PDELab::FastDGLocalJacobianApplyAssemblerEngine< FastDGLocalAssembler >
localassembler.hh
Dune::PDELab::FastDGLocalAssembler::localJacobianAssemblerEngine
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: fastdg/localassembler.hh:176
Dune
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Dune::PDELab::FastDGGridOperator::Range
Dune::PDELab::Backend::Vector< GFSV, RF > Range
The type of the range (residual).
Definition: fastdg.hh:43
Dune::PDELab::FastDGLocalPatternAssemblerEngine< FastDGLocalAssembler >
Dune::PDELab::Backend::Matrix
typename impl::BackendMatrixSelector< Backend, VU, VV, E >::Type Matrix
alias of the return type of BackendMatrixSelector
Definition: backend/interface.hh:127
Dune::PDELab::GridOperatorTraits
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:33
Dune::PDELab::FastDGAssembler
The fast DG assembler for standard DUNE grid.
Definition: fastdg/assembler.hh:25
Dune::PDELab::FastDGGridOperator::Assembler
FastDGAssembler< GFSU, GFSV, CU, CV > Assembler
The global assembler type.
Definition: fastdg.hh:38
Dune::PDELab::FastDGGridOperator::FastDGGridOperator
FastDGGridOperator(const GFSU &gfsu_, const GFSV &gfsv_, LOP &lop_, const MB &mb_=MB())
Constructor for empty constraints.
Definition: fastdg.hh:85
gridoperatorutilities.hh
Dune::PDELab::FastDGGridOperator::jacobian
void jacobian(const Domain &x, Jacobian &a) const
Assembler jacobian.
Definition: fastdg.hh:221
Dune::PDELab::FastDGGridOperator::make_consistent
void make_consistent(Jacobian &a) const
Definition: fastdg.hh:259
Dune::PDELab::FastDGGridOperator::SetupGridOperator::SetupGridOperator
SetupGridOperator()
Definition: fastdg.hh:159
Dune::PDELab::GridOperatorTraits::MatrixBackend
MB MatrixBackend
The matrix backend of the grid operator.
Definition: gridoperatorutilities.hh:51
Dune::PDELab::FastDGAssembler::testGridFunctionSpace
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: fastdg/assembler.hh:76
Dune::PDELab::FastDGGridOperator::SetupGridOperator::visit
void visit(T &elem)
Definition: fastdg.hh:163
Dune::PDELab::FastDGGridOperator::nonlinear_jacobian_apply
void nonlinear_jacobian_apply(const Domain &x, const Domain &z, Range &r) const
Apply jacobian matrix without explicitly assembling it.
Definition: fastdg.hh:247
Dune::PDELab::FastDGGridOperator::assembler
Assembler & assembler()
Definition: fastdg.hh:147
Dune::PDELab::FastDGGridOperator::MatrixContainer::Type
Traits::Jacobian Type
Definition: fastdg.hh:70
Dune::PDELab::FastDGGridOperator::assembler
const Assembler & assembler() const
Definition: fastdg.hh:149
Dune::PDELab::FastDGGridOperator::matrixBackend
const Traits::MatrixBackend & matrixBackend() const
Get the matrix backend for this grid operator.
Definition: fastdg.hh:271
borderdofexchanger.hh
Dune::PDELab::FastDGGridOperator::testGridFunctionSpace
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: fastdg.hh:102
Dune::PDELab::FastDGGridOperator::interpolate
void interpolate(const X &xold, F &f, X &x) const
Interpolate the constrained dofs from given function.
Definition: fastdg.hh:185
Dune::PDELab::EmptyTransformation
Definition: constraintstransformation.hh:111
Dune::PDELab::FastDGGridOperator
Definition: fastdg.hh:33
Dune::PDELab::FastDGGridOperator::testGridFunctionSpaceConstraints
CV & testGridFunctionSpaceConstraints()
Definition: fastdg.hh:116
Dune::PDELab::FastDGGridOperator::globalSizeU
GFSU::Traits::SizeType globalSizeU() const
Get dimension of space u.
Definition: fastdg.hh:136
Dune::PDELab::GridOperatorTraits::Jacobian
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
Dune::PDELab::FastDGLocalAssembler::localJacobianApplyAssemblerEngine
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: fastdg/localassembler.hh:186
interpolate.hh
Dune::PDELab::FastDGLocalAssembler::localResidualAssemblerEngine
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: fastdg/localassembler.hh:166
Dune::PDELab::FastDGGridOperator::trialGridFunctionSpaceConstraints
CU & trialGridFunctionSpaceConstraints()
Definition: fastdg.hh:107
Dune::PDELab::FastDGLocalJacobianAssemblerEngine< FastDGLocalAssembler >
Dune::PDELab::FastDGAssembler::assemble
void assemble(const EngineFactory &engineFactory, LocalAssembler &la) const
do the assembly
Definition: fastdg/assembler.hh:93
Dune::PDELab::FastDGGridOperator::SetupGridOperator
Definition: fastdg.hh:157
Dune::PDELab::OverlappingBorderDOFExchanger
Definition: borderdofexchanger.hh:577
Dune::PDELab::FastDGGridOperator::globalSizeV
GFSV::Traits::SizeType globalSizeV() const
Get dimension of space v.
Definition: fastdg.hh:142
Dune::PDELab::FastDGGridOperator::LocalAssembler
FastDGLocalAssembler< FastDGGridOperator, LOP, GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_Partition > LocalAssembler
The local assembler type.
Definition: fastdg.hh:55
Dune::PDELab::FastDGGridOperator::trialGridFunctionSpaceConstraints
const CU & trialGridFunctionSpaceConstraints() const
Definition: fastdg.hh:111
Dune::PDELab::FastDGLocalAssembler::localNonlinearJacobianApplyAssemblerEngine
LocalNonlinearJacobianApplyAssemblerEngine & localNonlinearJacobianApplyAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x, const typename Traits::Solution &z)
Definition: fastdg/localassembler.hh:196
Dune::PDELab::FastDGAssembler::trialGridFunctionSpace
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: fastdg/assembler.hh:70
Dune::PDELab::FastDGGridOperator::localAssembler
LocalAssembler & localAssembler() const
Definition: fastdg.hh:151
Dune::PDELab::FastDGGridOperator::Pattern
MB::template Pattern< Jacobian, GFSV, GFSU > Pattern
The sparsity pattern container for the jacobian matrix.
Definition: fastdg.hh:48
Dune::PDELab::FastDGLocalNonlinearJacobianApplyAssemblerEngine< FastDGLocalAssembler >
Dune::PDELab::NonOverlappingBorderDOFExchanger
Helper class for adding up matrix entries on border.
Definition: borderdofexchanger.hh:67
Dune::PDELab::LocalAssemblerBase::trialConstraints
const CU & trialConstraints() const
get the constraints on the trial grid function space
Definition: assemblerutilities.hh:199
value
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
Dune::PDELab::Backend::Vector
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:106
Dune::PDELab::FastDGGridOperator::jacobian_apply
void jacobian_apply(const Domain &x, Range &r) const
Apply jacobian matrix without explicitly assembling it.
Definition: fastdg.hh:234
Dune::PDELab::interpolate
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:204
Dune::PDELab::FastDGGridOperator::FastDGGridOperator
FastDGGridOperator(const GFSU &gfsu_, const CU &cu_, const GFSV &gfsv_, const CV &cv_, LOP &lop_, const MB &mb_=MB())
Constructor for non trivial constraints.
Definition: fastdg.hh:74
Dune::PDELab::FastDGGridOperator::setupGridOperators
static void setupGridOperators(GridOperatorTuple tuple)
Definition: fastdg.hh:177
Dune::PDELab::FastDGGridOperator::update
void update()
Definition: fastdg.hh:264
Dune::PDELab::FastDGGridOperator::residual
void residual(const Domain &x, Range &r) const
Assemble residual.
Definition: fastdg.hh:208
Dune::PDELab::FastDGGridOperator::Jacobian
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: fastdg.hh:45
p
const P & p
Definition: constraints.hh:147