dune-pdelab  2.7-git
fastdg/residualengine.hh
Go to the documentation of this file.
1 #ifndef DUNE_PDELAB_GRIDOPERATOR_FASTDG_RESIDUALENGINE_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_FASTDG_RESIDUALENGINE_HH
3 
10 
11 namespace Dune{
12  namespace PDELab{
13 
21  template<typename LA>
24  {
25  public:
26 
27  template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
28  bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv) const
29  {
30  return false;
31  }
32 
34  typedef LA LocalAssembler;
35 
37  typedef typename LA::LocalOperator LOP;
38 
40  typedef typename LA::Traits::Residual Residual;
41  typedef typename Residual::ElementType ResidualElement;
42 
44  typedef typename LA::Traits::Solution Solution;
45  typedef typename Solution::ElementType SolutionElement;
46 
48  typedef typename LA::LFSU LFSU;
49  typedef typename LA::LFSUCache LFSUCache;
50  typedef typename LFSU::Traits::GridFunctionSpace GFSU;
51  typedef typename LA::LFSV LFSV;
52  typedef typename LA::LFSVCache LFSVCache;
53  typedef typename LFSV::Traits::GridFunctionSpace GFSV;
54 
55  typedef typename Solution::template ConstAliasedLocalView<LFSUCache> SolutionView;
56  typedef typename Residual::template AliasedLocalView<LFSVCache> ResidualView;
57 
65  : local_assembler(local_assembler_),
66  lop(local_assembler_.localOperator())
67  //rl_view(rl,1.0),
68  //rn_view(rn,1.0)
69  {}
70 
72 
87  local_assembler(other.local_assembler), lop(other.lop),
88  global_rl_view(other.global_rl_view),
89  global_rn_view(other.global_rn_view),
90  global_sl_view(other.global_sl_view),
91  global_sn_view(other.global_sn_view)
92  //rl_view(rl,1.0),
93  //rn_view(rn,1.0)
94  { }
95 
98  bool requireSkeleton() const
99  { return ( local_assembler.doAlphaSkeleton() || local_assembler.doLambdaSkeleton() ); }
101  { return local_assembler.doSkeletonTwoSided(); }
102  bool requireUVVolume() const
103  { return local_assembler.doAlphaVolume(); }
104  bool requireVVolume() const
105  { return local_assembler.doLambdaVolume(); }
106  bool requireUVSkeleton() const
107  { return local_assembler.doAlphaSkeleton(); }
108  bool requireVSkeleton() const
109  { return local_assembler.doLambdaSkeleton(); }
110  bool requireUVBoundary() const
111  { return local_assembler.doAlphaBoundary(); }
112  bool requireVBoundary() const
113  { return local_assembler.doLambdaBoundary(); }
115  { return local_assembler.doAlphaVolumePostSkeleton(); }
117  { return local_assembler.doLambdaVolumePostSkeleton(); }
119 
122  {
123  return local_assembler;
124  }
125 
127  const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
128  {
129  return localAssembler().trialConstraints();
130  }
131 
133  const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
134  {
135  return localAssembler().testConstraints();
136  }
137 
140  void setResidual(Residual & residual_)
141  {
142  global_rl_view.attach(residual_);
143  global_rn_view.attach(residual_);
144  }
145 
148  void setSolution(const Solution & solution_)
149  {
150  global_sl_view.attach(solution_);
151  global_sn_view.attach(solution_);
152  }
153 
157  template<typename EG, typename LFSUC, typename LFSVC>
158  void onBindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
159  {
160  global_sl_view.bind(lfsu_cache);
161  //xl.resize(lfsu_cache.size());
162  }
163 
164  template<typename EG, typename LFSVC>
165  void onBindLFSV(const EG & eg, const LFSVC & lfsv_cache)
166  {
167  global_rl_view.bind(lfsv_cache);
168  //rl.assign(lfsv_cache.size(),0.0);
169  }
170 
171  template<typename IG, typename LFSUC, typename LFSVC>
172  void onBindLFSUVInside(const IG & ig, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
173  {
174  global_sl_view.bind(lfsu_cache);
175  //xl.resize(lfsu_cache.size());
176  }
177 
178  template<typename IG, typename LFSUC, typename LFSVC>
179  void onBindLFSUVOutside(const IG & ig,
180  const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
181  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
182  {
183  global_sn_view.bind(lfsu_n_cache);
184  //xn.resize(lfsu_n_cache.size());
185  }
186 
187  template<typename IG, typename LFSVC>
188  void onBindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
189  {
190  global_rl_view.bind(lfsv_cache);
191  //rl.assign(lfsv_cache.size(),0.0);
192  }
193 
194  template<typename IG, typename LFSVC>
195  void onBindLFSVOutside(const IG & ig,
196  const LFSVC & lfsv_s_cache,
197  const LFSVC & lfsv_n_cache)
198  {
199  global_rn_view.bind(lfsv_n_cache);
200  //rn.assign(lfsv_n_cache.size(),0.0);
201  }
202 
204 
208  template<typename EG, typename LFSVC>
209  void onUnbindLFSV(const EG & eg, const LFSVC & lfsv_cache)
210  {
211  //global_rl_view.add(rl);
212  global_rl_view.unbind();
213  global_rl_view.commit();
214  }
215 
216  template<typename IG, typename LFSVC>
217  void onUnbindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
218  {
219  //global_rl_view.add(rl);
220  global_rl_view.commit();
221  global_rl_view.unnbind();
222  }
223 
224  template<typename IG, typename LFSVC>
225  void onUnbindLFSVOutside(const IG & ig,
226  const LFSVC & lfsv_s_cache,
227  const LFSVC & lfsv_n_cache)
228  {
229  //global_rn_view.add(rn);
230  global_rn_view.commit();
231  global_rn_view.unbind();
232  }
234 
237  template<typename LFSUC>
238  void loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)
239  {
240  //global_sl_view.read(xl);
241  }
242  template<typename LFSUC>
243  void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
244  {
245  //global_sn_view.read(xn);
246  }
247  template<typename LFSUC>
248  void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
249  {
250  DUNE_THROW(Dune::NotImplemented,"No coupling lfsu available for ");
251  }
253 
256 
257  void postAssembly(const GFSU& gfsu, const GFSV& gfsv)
258  {
259  if(local_assembler.doPostProcessing())
260  Dune::PDELab::constrain_residual(local_assembler.testConstraints(),
261  global_rl_view.container());
262 
263  }
264 
266 
269 
276  template<typename EG>
277  bool assembleCell(const EG & eg)
278  {
279  return LocalAssembler::isNonOverlapping && eg.entity().partitionType() != Dune::InteriorEntity;
280  }
281 
282  template<typename EG, typename LFSUC, typename LFSVC>
283  void assembleUVVolume(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
284  {
285  //rl_view.setWeight(local_assembler.weight());
286  global_rl_view.setWeight(local_assembler.weight());
288  alpha_volume(lop,eg,lfsu_cache.localFunctionSpace(),global_sl_view,lfsv_cache.localFunctionSpace(),global_rl_view);
289  }
290 
291  template<typename EG, typename LFSVC>
292  void assembleVVolume(const EG & eg, const LFSVC & lfsv_cache)
293  {
294  //rl_view.setWeight(local_assembler.weight());
295  global_rl_view.setWeight(local_assembler.weight());
297  lambda_volume(lop,eg,lfsv_cache.localFunctionSpace(),global_rl_view);
298  }
299 
300  template<typename IG, typename LFSUC, typename LFSVC>
301  void assembleUVSkeleton(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
302  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
303  {
304  //rl_view.setWeight(local_assembler.weight());
305  //rn_view.setWeight(local_assembler.weight());
306  global_rl_view.setWeight(local_assembler.weight());
307  global_rn_view.setWeight(local_assembler.weight());
309  alpha_skeleton(lop,ig,
310  lfsu_s_cache.localFunctionSpace(),global_sl_view,lfsv_s_cache.localFunctionSpace(),
311  lfsu_n_cache.localFunctionSpace(),global_sn_view,lfsv_n_cache.localFunctionSpace(),
312  global_rl_view,global_rn_view);
313  }
314 
315  template<typename IG, typename LFSVC>
316  void assembleVSkeleton(const IG & ig, const LFSVC & lfsv_s_cache, const LFSVC & lfsv_n_cache)
317  {
318  //rl_view.setWeight(local_assembler.weight());
319  //rn_view.setWeight(local_assembler.weight());
320  global_rl_view.setWeight(local_assembler.weight());
321  global_rn_view.setWeight(local_assembler.weight());
323  lambda_skeleton(lop, ig, lfsv_s_cache.localFunctionSpace(), lfsv_n_cache.localFunctionSpace(), global_rl_view, global_rn_view);
324  }
325 
326  template<typename IG, typename LFSUC, typename LFSVC>
327  void assembleUVBoundary(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache)
328  {
329  //rl_view.setWeight(local_assembler.weight());
330  global_rl_view.setWeight(local_assembler.weight());
332  alpha_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),global_sl_view,lfsv_s_cache.localFunctionSpace(),global_rl_view);
333  }
334 
335  template<typename IG, typename LFSVC>
336  void assembleVBoundary(const IG & ig, const LFSVC & lfsv_s_cache)
337  {
338  //rl_view.setWeight(local_assembler.weight());
339  global_rl_view.setWeight(local_assembler.weight());
341  lambda_boundary(lop,ig,lfsv_s_cache.localFunctionSpace(),global_rl_view);
342  }
343 
344  template<typename IG, typename LFSUC, typename LFSVC>
345  static void assembleUVEnrichedCoupling(const IG & ig,
346  const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
347  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache,
348  const LFSUC & lfsu_coupling_cache, const LFSVC & lfsv_coupling_cache)
349  {
350  DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
351  }
352 
353  template<typename IG, typename LFSVC>
354  static void assembleVEnrichedCoupling(const IG & ig,
355  const LFSVC & lfsv_s_cache,
356  const LFSVC & lfsv_n_cache,
357  const LFSVC & lfsv_coupling_cache)
358  {
359  DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
360  }
361 
362  template<typename EG, typename LFSUC, typename LFSVC>
363  void assembleUVVolumePostSkeleton(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
364  {
365  //rl_view.setWeight(local_assembler.weight());
366  global_rl_view.setWeight(local_assembler.weight());
368  alpha_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),global_sl_view,lfsv_cache.localFunctionSpace(),global_rl_view);
369  }
370 
371  template<typename EG, typename LFSVC>
372  void assembleVVolumePostSkeleton(const EG & eg, const LFSVC & lfsv_cache)
373  {
374  //rl_view.setWeight(local_assembler.weight());
375  global_rl_view.setWeight(local_assembler.weight());
377  lambda_volume_post_skeleton(lop,eg,lfsv_cache.localFunctionSpace(),global_rl_view);
378  }
379 
381 
382  private:
385  const LocalAssembler & local_assembler;
386 
388  const LOP & lop;
389 
391  ResidualView global_rl_view;
392  ResidualView global_rn_view;
393 
395  SolutionView global_sl_view;
396  SolutionView global_sn_view;
397 
400  typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
401  typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
402 
405 
407  //typename ResidualVector::WeightedAccumulationView rl_view;
409  //typename ResidualVector::WeightedAccumulationView rn_view;
411 
412  }; // End of class FastDGLocalResidualAssemblerEngine
413 
414  }
415 }
416 #endif
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireSkeletonTwoSided
bool requireSkeletonTwoSided() const
Definition: fastdg/residualengine.hh:100
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleUVSkeleton
void assembleUVSkeleton(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache, const LFSUC &lfsu_n_cache, const LFSVC &lfsv_n_cache)
Definition: fastdg/residualengine.hh:301
gridoperatorutilities.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine::needsConstraintsCaching
bool needsConstraintsCaching(const TrialConstraintsContainer &cu, const TestConstraintsContainer &cv) const
Definition: fastdg/residualengine.hh:28
Dune::PDELab::TrialSpaceTag
Definition: localfunctionspacetags.hh:48
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSVOutside
void onBindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: fastdg/residualengine.hh:195
Dune::PDELab::FastDGLocalResidualAssemblerEngine::SolutionElement
Solution::ElementType SolutionElement
Definition: fastdg/residualengine.hh:45
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LFSU
LA::LFSU LFSU
The local function spaces.
Definition: fastdg/residualengine.hh:48
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleUVBoundary
void assembleUVBoundary(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache)
Definition: fastdg/residualengine.hh:327
Dune::PDELab::FastDGLocalResidualAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: fastdg/residualengine.hh:148
Dune::PDELab::FastDGLocalResidualAssemblerEngine::postAssembly
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: fastdg/residualengine.hh:257
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleUVVolumePostSkeleton
void assembleUVVolumePostSkeleton(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:363
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireVVolumePostSkeleton
bool requireVVolumePostSkeleton() const
Definition: fastdg/residualengine.hh:116
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LFSUCache
LA::LFSUCache LFSUCache
Definition: fastdg/residualengine.hh:49
Dune::PDELab::FastDGLocalResidualAssemblerEngine::setResidual
void setResidual(Residual &residual_)
Definition: fastdg/residualengine.hh:140
Dune
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Dune::PDELab::FastDGLocalResidualAssemblerEngine::Solution
LA::Traits::Solution Solution
The type of the solution vector.
Definition: fastdg/residualengine.hh:44
localassemblerenginebase.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSVInside
void onBindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:188
Dune::PDELab::FastDGLocalResidualAssemblerEngine::loadCoefficientsLFSUInside
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: fastdg/residualengine.hh:238
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireUVVolume
bool requireUVVolume() const
Definition: fastdg/residualengine.hh:102
Dune::PDELab::FastDGLocalResidualAssemblerEngine::ResidualView
Residual::template AliasedLocalView< LFSVCache > ResidualView
Definition: fastdg/residualengine.hh:56
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LOP
LA::LocalOperator LOP
The type of the local operator.
Definition: fastdg/residualengine.hh:37
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireVSkeleton
bool requireVSkeleton() const
Definition: fastdg/residualengine.hh:108
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSUVInside
void onBindLFSUVInside(const IG &ig, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:172
Dune::PDELab::constrain_residual
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:904
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleVVolumePostSkeleton
void assembleVVolumePostSkeleton(const EG &eg, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:372
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LocalAssembler
LA LocalAssembler
The type of the wrapping local assembler.
Definition: fastdg/residualengine.hh:34
Dune::PDELab::FastDGLocalResidualAssemblerEngine::FastDGLocalResidualAssemblerEngine
FastDGLocalResidualAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: fastdg/residualengine.hh:64
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleVVolume
void assembleVVolume(const EG &eg, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:292
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireUVSkeleton
bool requireUVSkeleton() const
Definition: fastdg/residualengine.hh:106
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleCell
bool assembleCell(const EG &eg)
Definition: fastdg/residualengine.hh:277
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleVEnrichedCoupling
static void assembleVEnrichedCoupling(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache, const LFSVC &lfsv_coupling_cache)
Definition: fastdg/residualengine.hh:354
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSUV
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:158
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireVBoundary
bool requireVBoundary() const
Definition: fastdg/residualengine.hh:112
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LFSVCache
LA::LFSVCache LFSVCache
Definition: fastdg/residualengine.hh:52
Dune::PDELab::FastDGLocalResidualAssemblerEngine::SolutionView
Solution::template ConstAliasedLocalView< LFSUCache > SolutionView
Definition: fastdg/residualengine.hh:55
localvector.hh
Dune::PDELab::TestSpaceTag
Definition: localfunctionspacetags.hh:54
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSV
void onBindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:165
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireUVBoundary
bool requireUVBoundary() const
Definition: fastdg/residualengine.hh:110
callswitch.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine::ResidualElement
Residual::ElementType ResidualElement
Definition: fastdg/residualengine.hh:41
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onUnbindLFSVInside
void onUnbindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:217
Dune::PDELab::LocalAssemblerCallSwitch
Impl::LocalAssemblerCallSwitchHelper< LOP, doIt > LocalAssemblerCallSwitch
Definition: callswitch.hh:344
Dune::PDELab::LocalAssemblerEngineBase
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:21
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onUnbindLFSV
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:209
Dune::PDELab::FastDGLocalAssembler
The local assembler for DUNE grids.
Definition: fastdg/localassembler.hh:37
Dune::PDELab::FastDGLocalResidualAssemblerEngine::trialConstraints
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: fastdg/residualengine.hh:127
Dune::PDELab::LocalVector< SolutionElement, LocalTrialSpaceTag >
Dune::PDELab::FastDGLocalResidualAssemblerEngine::testConstraints
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: fastdg/residualengine.hh:133
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleVBoundary
void assembleVBoundary(const IG &ig, const LFSVC &lfsv_s_cache)
Definition: fastdg/residualengine.hh:336
Dune::PDELab::FastDGLocalResidualAssemblerEngine::GFSV
LFSV::Traits::GridFunctionSpace GFSV
Definition: fastdg/residualengine.hh:53
Dune::PDELab::FastDGLocalResidualAssemblerEngine::loadCoefficientsLFSUCoupling
void loadCoefficientsLFSUCoupling(const LFSUC &lfsu_c_cache)
Definition: fastdg/residualengine.hh:248
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleUVVolume
void assembleUVVolume(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: fastdg/residualengine.hh:283
constraints.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine::loadCoefficientsLFSUOutside
void loadCoefficientsLFSUOutside(const LFSUC &lfsu_n_cache)
Definition: fastdg/residualengine.hh:243
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireVVolume
bool requireVVolume() const
Definition: fastdg/residualengine.hh:104
Dune::PDELab::FastDGLocalResidualAssemblerEngine::localAssembler
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: fastdg/residualengine.hh:121
Dune::PDELab::FastDGLocalResidualAssemblerEngine
The fast DG local assembler engine for DUNE grids which assembles the residual vector.
Definition: fastdg/residualengine.hh:22
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleUVEnrichedCoupling
static void assembleUVEnrichedCoupling(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache, const LFSUC &lfsu_n_cache, const LFSVC &lfsv_n_cache, const LFSUC &lfsu_coupling_cache, const LFSVC &lfsv_coupling_cache)
Definition: fastdg/residualengine.hh:345
ig
const IG & ig
Definition: constraints.hh:149
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireSkeleton
bool requireSkeleton() const
Definition: fastdg/residualengine.hh:98
Dune::PDELab::FastDGLocalResidualAssemblerEngine::GFSU
LFSU::Traits::GridFunctionSpace GFSU
Definition: fastdg/residualengine.hh:50
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onBindLFSUVOutside
void onBindLFSUVOutside(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache, const LFSUC &lfsu_n_cache, const LFSVC &lfsv_n_cache)
Definition: fastdg/residualengine.hh:179
Dune::PDELab::FastDGLocalResidualAssemblerEngine::LFSV
LA::LFSV LFSV
Definition: fastdg/residualengine.hh:51
Dune::PDELab::FastDGLocalResidualAssemblerEngine::assembleVSkeleton
void assembleVSkeleton(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: fastdg/residualengine.hh:316
Dune::PDELab::FastDGLocalResidualAssemblerEngine::onUnbindLFSVOutside
void onUnbindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: fastdg/residualengine.hh:225
assemblerutilities.hh
Dune::PDELab::FastDGLocalResidualAssemblerEngine::requireUVVolumePostSkeleton
bool requireUVVolumePostSkeleton() const
Definition: fastdg/residualengine.hh:114
Dune::PDELab::FastDGLocalResidualAssemblerEngine::Residual
LA::Traits::Residual Residual
The type of the residual vector.
Definition: fastdg/residualengine.hh:40