DUNE PDELab (git)

fastdg.hh
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 
7 
8 #include <dune/pdelab/gridfunctionspace/interpolate.hh>
9 #include <dune/pdelab/gridoperator/common/borderdofexchanger.hh>
10 #include <dune/pdelab/gridoperator/common/gridoperatorutilities.hh>
11 #include <dune/pdelab/gridoperator/fastdg/assembler.hh>
12 #include <dune/pdelab/gridoperator/fastdg/localassembler.hh>
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,
31  typename CU=Dune::PDELab::EmptyTransformation,
32  typename CV=Dune::PDELab::EmptyTransformation>
34  {
35  public:
36 
39 
41  using Domain = Dune::PDELab::Backend::Vector<GFSU,DF>;
43  using Range = Dune::PDELab::Backend::Vector<GFSV,RF>;
45  using Jacobian = Dune::PDELab::Backend::Matrix<MB,Domain,Range,JF>;
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,
61  OverlappingBorderDOFExchanger<FastDGGridOperator>
62  >::type BorderDOFExchanger;
63 
66  <GFSU,GFSV,MB,DF,RF,JF,CU,CV,Assembler,LocalAssembler> Traits;
67 
68  template <typename MFT>
69  struct MatrixContainer{
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 
107  CU & trialGridFunctionSpaceConstraints()
108  {
109  return cu;
110  }
111  const CU & trialGridFunctionSpaceConstraints() const
112  {
113  return cu;
114  }
115 
116  CV & testGridFunctionSpaceConstraints()
117  {
118  return cv;
119  }
120  const CV & testGridFunctionSpaceConstraints() const
121  {
122  return cv;
123  }
124 
125  LOP & localOperator()
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,
181  [&](auto &el) { setup_visitor.visit(el); });
182  }
183 
185  template<typename F, typename X>
186  void interpolate (const X& xold, F& f, X& x) const
187  {
188  // Interpolate f into grid function space and set corresponding coefficients
189  Dune::PDELab::interpolate(f,global_assembler.trialGridFunctionSpace(),x);
190 
191  // Copy non-constrained dofs from old time step
193  }
194 
196  void fill_pattern(Pattern & p) const
197  {
198  typedef typename LocalAssembler::LocalPatternAssemblerEngine PatternEngine;
199  PatternEngine & pattern_engine = local_assembler.localPatternAssemblerEngine(p);
200  global_assembler.assemble(pattern_engine);
201  }
202 
204  void residual(const Domain & x, Range & r) const
205  {
206  typedef typename LocalAssembler::LocalResidualAssemblerEngine ResidualEngine;
207  ResidualEngine & residual_engine = local_assembler.localResidualAssemblerEngine(r,x);
208  global_assembler.assemble(residual_engine);
209  }
210 
212  void jacobian(const Domain & x, Jacobian & a) const
213  {
214  typedef typename LocalAssembler::LocalJacobianAssemblerEngine JacobianEngine;
215  JacobianEngine & jacobian_engine = local_assembler.localJacobianAssemblerEngine(a,x);
216  global_assembler.assemble(jacobian_engine);
217  }
218 
220  void jacobian_apply(const Domain & update, Range & result) const
221  {
222  if (not local_assembler.localOperator().isLinear)
223  DUNE_THROW(Dune::Exception, "Your trying to use a linear jacobian apply for a non linear problem.");
224  global_assembler.assemble(local_assembler.localJacobianApplyAssemblerEngine(update, result));
225  }
226 
228  void jacobian_apply(const Domain & solution, const Domain & update, Range & result) const
229  {
230  if (local_assembler.localOperator().isLinear)
231  DUNE_THROW(Dune::Exception, "Your trying to use a non linear jacobian apply for a linear problem.");
232  global_assembler.assemble(local_assembler.localJacobianApplyAssemblerEngine(solution, update, result));
233  }
234 
236  [[deprecated("nonlinear_jacobian_apply(x,z,r) is deprecated. Please use jacobian_apply(solution, update, result) instead!")]]
237  void nonlinear_jacobian_apply(const Domain & solution, const Domain & update, Range & result) const
238  {
239  if (local_assembler.localOperator().isLinear)
240  DUNE_THROW(Dune::Exception, "Your trying to use a non linear jacobian apply for a linear problem.");
241  global_assembler.assemble(local_assembler.localJacobianApplyAssemblerEngine(solution, update, result));
242  }
243 
244  void make_consistent(Jacobian& a) const
245  {
246  dof_exchanger->accumulateBorderEntries(*this,a);
247  }
248 
249  void update()
250  {
251  // the DOF exchanger has matrix information, so we need to update it
252  dof_exchanger->update(*this);
253  }
254 
256  const typename Traits::MatrixBackend& matrixBackend() const
257  {
258  return backend;
259  }
260 
261  private:
262  Assembler global_assembler;
263  std::shared_ptr<BorderDOFExchanger> dof_exchanger;
264 
265  mutable LocalAssembler local_assembler;
266  MB backend;
267 
268  const LOP& lop;
269  const CU& cu;
270  const CV& cv;
271 
272  }; // end class FastDGGridOperator
273  } // end namespace PDELab
274 } // end namespace Dune
275 #endif
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
The fast DG assembler for standard DUNE grid.
Definition: assembler.hh:25
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: assembler.hh:70
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: assembler.hh:76
Definition: fastdg.hh:34
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: fastdg.hh:102
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::GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, Assembler, LocalAssembler > Traits
The grid operator traits.
Definition: fastdg.hh:66
FastDGLocalAssembler< FastDGGridOperator, LOP, GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_Partition > LocalAssembler
The local assembler type.
Definition: fastdg.hh:55
FastDGAssembler< GFSU, GFSV, CU, CV > Assembler
The global assembler type.
Definition: fastdg.hh:38
GFSV::Traits::SizeType globalSizeV() const
Get dimension of space v.
Definition: fastdg.hh:142
Dune::PDELab::Backend::Vector< GFSV, RF > Range
The type of the range (residual).
Definition: fastdg.hh:43
MB::template Pattern< Jacobian, GFSV, GFSU > Pattern
The sparsity pattern container for the jacobian matrix.
Definition: fastdg.hh:48
GFSU::Traits::SizeType globalSizeU() const
Get dimension of space u.
Definition: fastdg.hh:136
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: fastdg.hh:45
void jacobian_apply(const Domain &update, Range &result) const
Apply jacobian matrix to the vector update without explicitly assembling it.
Definition: fastdg.hh:220
FastDGGridOperator(const GFSU &gfsu_, const GFSV &gfsv_, LOP &lop_, const MB &mb_=MB())
Constructor for empty constraints.
Definition: fastdg.hh:85
void jacobian_apply(const Domain &solution, const Domain &update, Range &result) const
Apply jacobian matrix to the vector update without explicitly assembling it.
Definition: fastdg.hh:228
const Traits::MatrixBackend & matrixBackend() const
Get the matrix backend for this grid operator.
Definition: fastdg.hh:256
static void setupGridOperators(GridOperatorTuple tuple)
Definition: fastdg.hh:177
void interpolate(const X &xold, F &f, X &x) const
Interpolate the constrained dofs from given function.
Definition: fastdg.hh:186
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: fastdg.hh:41
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: fastdg.hh:96
void residual(const Domain &x, Range &r) const
Assemble residual.
Definition: fastdg.hh:204
void jacobian(const Domain &x, Jacobian &a) const
Assembler jacobian.
Definition: fastdg.hh:212
void nonlinear_jacobian_apply(const Domain &solution, const Domain &update, Range &result) const
Apply jacobian matrix to the vector update without explicitly assembling it.
Definition: fastdg.hh:237
void fill_pattern(Pattern &p) const
Fill pattern of jacobian matrix.
Definition: fastdg.hh:196
The local assembler for DUNE grids.
Definition: localassembler.hh:41
LOP & localOperator()
get a reference to the local operator
Definition: localassembler.hh:111
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(const typename Traits::Domain &update, typename Traits::Range &result)
Definition: localassembler.hh:226
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: localassembler.hh:206
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: localassembler.hh:216
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: localassembler.hh:197
const CU & trialConstraints() const
get the constraints on the trial grid function space
Definition: assemblerutilities.hh:233
Helper class for adding up matrix entries on border.
Definition: borderdofexchanger.hh:68
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
@ InteriorBorder_Partition
interior and border entities
Definition: gridenums.hh:138
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:987
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:256
Dune namespace.
Definition: alignedallocator.hh:13
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:34
MB MatrixBackend
The matrix backend of the grid operator.
Definition: gridoperatorutilities.hh:51
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
Contains utility classes which can be used with std::tuple.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)