dune-pdelab  2.5-dev
default/nonlinearjacobianapplyengine.hh
Go to the documentation of this file.
1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH
3 
9 
10 namespace Dune{
11  namespace PDELab{
12 
20  template<typename LA>
23  {
24  public:
25 
26  template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
27  bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv) const
28  {
29  return false;
30  }
31 
33  typedef LA LocalAssembler;
34 
36  typedef typename LA::LocalOperator LOP;
37 
39  typedef typename LA::Traits::Residual Residual;
40  typedef typename Residual::ElementType ResidualElement;
41 
43  typedef typename LA::Traits::Solution Solution;
44  typedef typename Solution::ElementType SolutionElement;
45 
47  typedef typename LA::LFSU LFSU;
48  typedef typename LA::LFSUCache LFSUCache;
49  typedef typename LFSU::Traits::GridFunctionSpace GFSU;
50  typedef typename LA::LFSV LFSV;
51  typedef typename LA::LFSVCache LFSVCache;
52  typedef typename LFSV::Traits::GridFunctionSpace GFSV;
53 
54  typedef typename Solution::template ConstLocalView<LFSUCache> SolutionView;
55  typedef typename Residual::template LocalView<LFSVCache> ResidualView;
56 
64  : local_assembler(local_assembler_), lop(local_assembler_.localOperator()),
65  rl_view(rl,1.0),
66  rn_view(rn,1.0)
67  {}
68 
71  bool requireSkeleton() const
72  { return local_assembler.doAlphaSkeleton(); }
74  { return local_assembler.doSkeletonTwoSided(); }
75  bool requireUVVolume() const
76  { return local_assembler.doAlphaVolume(); }
77  bool requireUVSkeleton() const
78  { return local_assembler.doAlphaSkeleton(); }
79  bool requireUVBoundary() const
80  { return local_assembler.doAlphaBoundary(); }
82  { return local_assembler.doAlphaVolumePostSkeleton(); }
84 
87  {
88  return local_assembler;
89  }
90 
92  const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
93  {
94  return localAssembler().trialConstraints();
95  }
96 
98  const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
99  {
100  return localAssembler().testConstraints();
101  }
102 
105  void setResidual(Residual & residual_)
106  {
107  global_rl_view.attach(residual_);
108  global_rn_view.attach(residual_);
109  }
110 
113  void setSolution(const Solution & solution_)
114  {
115  global_sl_view.attach(solution_);
116  global_sn_view.attach(solution_);
117  }
118 
121  void setUpdate(const Solution & update_)
122  {
123  global_zl_view.attach(update_);
124  global_zn_view.attach(update_);
125  }
126 
130  template<typename EG, typename LFSUC, typename LFSVC>
131  void onBindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
132  {
133  global_sl_view.bind(lfsu_cache);
134  xl.resize(lfsu_cache.size());
135  global_zl_view.bind(lfsu_cache);
136  zl.resize(lfsu_cache.size());
137  }
138 
139  template<typename EG, typename LFSVC>
140  void onBindLFSV(const EG & eg, const LFSVC & lfsv_cache)
141  {
142  global_rl_view.bind(lfsv_cache);
143  rl.assign(lfsv_cache.size(),0.0);
144  }
145 
146  template<typename IG, typename LFSUC, typename LFSVC>
147  void onBindLFSUVInside(const IG & ig, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
148  {
149  global_sl_view.bind(lfsu_cache);
150  xl.resize(lfsu_cache.size());
151  global_zl_view.bind(lfsu_cache);
152  zl.resize(lfsu_cache.size());
153  }
154 
155  template<typename IG, typename LFSUC, typename LFSVC>
156  void onBindLFSUVOutside(const IG & ig,
157  const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
158  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
159  {
160  global_sn_view.bind(lfsu_n_cache);
161  xn.resize(lfsu_n_cache.size());
162  global_zn_view.bind(lfsu_n_cache);
163  zn.resize(lfsu_n_cache.size());
164  }
165 
166  template<typename IG, typename LFSVC>
167  void onBindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
168  {
169  global_rl_view.bind(lfsv_cache);
170  rl.assign(lfsv_cache.size(),0.0);
171  }
172 
173  template<typename IG, typename LFSVC>
174  void onBindLFSVOutside(const IG & ig,
175  const LFSVC & lfsv_s_cache,
176  const LFSVC & lfsv_n_cache)
177  {
178  global_rn_view.bind(lfsv_n_cache);
179  rn.assign(lfsv_n_cache.size(),0.0);
180  }
181 
183 
187  template<typename EG, typename LFSVC>
188  void onUnbindLFSV(const EG & eg, const LFSVC & lfsv_cache)
189  {
190  global_rl_view.add(rl);
191  global_rl_view.commit();
192  }
193 
194  template<typename IG, typename LFSVC>
195  void onUnbindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
196  {
197  global_rl_view.add(rl);
198  global_rl_view.commit();
199  }
200 
201  template<typename IG, typename LFSVC>
202  void onUnbindLFSVOutside(const IG & ig,
203  const LFSVC & lfsv_s_cache,
204  const LFSVC & lfsv_n_cache)
205  {
206  global_rn_view.add(rn);
207  global_rn_view.commit();
208  }
210 
213  template<typename LFSUC>
214  void loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)
215  {
216  global_sl_view.read(xl);
217  global_zl_view.read(zl);
218  }
219  template<typename LFSUC>
220  void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
221  {
222  global_sn_view.read(xn);
223  global_zn_view.read(zn);
224  }
225  template<typename LFSUC>
226  void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
227  {
228  DUNE_THROW(Dune::NotImplemented,"No coupling lfsu available for ");
229  }
231 
234 
235  void postAssembly(const GFSU& gfsu, const GFSV& gfsv)
236  {
237  if(local_assembler.doPostProcessing())
238  Dune::PDELab::constrain_residual(local_assembler.testConstraints(),global_rl_view.container());
239  }
240 
242 
245 
252  template<typename EG>
253  bool assembleCell(const EG & eg)
254  {
255  return LocalAssembler::isNonOverlapping && eg.entity().partitionType() != Dune::InteriorEntity;
256  }
257 
258  template<typename EG, typename LFSUC, typename LFSVC>
259  void assembleUVVolume(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
260  {
261  rl_view.setWeight(local_assembler.weight());
263  nonlinear_jacobian_apply_volume(lop,eg,lfsu_cache.localFunctionSpace(),xl,zl,lfsv_cache.localFunctionSpace(),rl_view);
264  }
265 
266  template<typename IG, typename LFSUC, typename LFSVC>
267  void assembleUVSkeleton(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
268  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
269  {
270  rl_view.setWeight(local_assembler.weight());
271  rn_view.setWeight(local_assembler.weight());
274  lfsu_s_cache.localFunctionSpace(),xl,zl,lfsv_s_cache.localFunctionSpace(),
275  lfsu_n_cache.localFunctionSpace(),xn,zn,lfsv_n_cache.localFunctionSpace(),
276  rl_view,rn_view);
277  }
278 
279  template<typename IG, typename LFSUC, typename LFSVC>
280  void assembleUVBoundary(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache)
281  {
282  rl_view.setWeight(local_assembler.weight());
284  nonlinear_jacobian_apply_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,zl,lfsv_s_cache.localFunctionSpace(),rl_view);
285  }
286 
287  template<typename IG, typename LFSUC, typename LFSVC>
288  static void assembleUVEnrichedCoupling(const IG & ig,
289  const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
290  const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache,
291  const LFSUC & lfsu_coupling_cache, const LFSVC & lfsv_coupling_cache)
292  {
293  DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
294  }
295 
296  template<typename EG, typename LFSUC, typename LFSVC>
297  void assembleUVVolumePostSkeleton(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
298  {
299  rl_view.setWeight(local_assembler.weight());
301  nonlinear_jacobian_apply_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),xl,zl,lfsv_cache.localFunctionSpace(),rl_view);
302  }
303 
305 
306  private:
309  const LocalAssembler & local_assembler;
310 
312  const LOP & lop;
313 
315  ResidualView global_rl_view;
316  ResidualView global_rn_view;
317 
319  SolutionView global_sl_view;
320  SolutionView global_sn_view;
321 
323  SolutionView global_zl_view;
324  SolutionView global_zn_view;
325 
328  typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
329  typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
330 
333 
335  SolutionVector xl;
337  SolutionVector xn;
339  SolutionVector zl;
341  SolutionVector zn;
343  ResidualVector rl;
345  ResidualVector rn;
351 
352  }; // End of class DefaultLocalJacobianAssemblerEngine
353 
354  }
355 }
356 
357 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH
Dune::PDELab::TestSpaceTag
Definition: localfunctionspacetags.hh:54
assemblerutilities.hh
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::testConstraints
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: default/nonlinearjacobianapplyengine.hh:98
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::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: default/nonlinearjacobianapplyengine.hh:288
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::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: default/nonlinearjacobianapplyengine.hh:156
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::DefaultLocalNonlinearJacobianApplyAssemblerEngine
DefaultLocalNonlinearJacobianApplyAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: default/nonlinearjacobianapplyengine.hh:63
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::assembleUVBoundary
void assembleUVBoundary(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache)
Definition: default/nonlinearjacobianapplyengine.hh:280
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine
The local assembler engine for DUNE grids which assembles the local application of the Jacobian.
Definition: default/nonlinearjacobianapplyengine.hh:21
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LOP
LA::LocalOperator LOP
The type of the local operator.
Definition: default/nonlinearjacobianapplyengine.hh:36
Dune::PDELab::WeightedVectorAccumulationView::setWeight
void setWeight(weight_type weight)
Resets the weighting coefficient of the view.
Definition: localvector.hh:72
Dune::PDELab::LocalVector< SolutionElement, LocalTrialSpaceTag >
callswitch.hh
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::localAssembler
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: default/nonlinearjacobianapplyengine.hh:86
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onBindLFSUVInside
void onBindLFSUVInside(const IG &ig, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:147
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::ResidualView
Residual::template LocalView< LFSVCache > ResidualView
Definition: default/nonlinearjacobianapplyengine.hh:55
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::postAssembly
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: default/nonlinearjacobianapplyengine.hh:235
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireUVVolume
bool requireUVVolume() const
Definition: default/nonlinearjacobianapplyengine.hh:75
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::assembleUVVolume
void assembleUVVolume(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:259
Dune
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::assembleUVVolumePostSkeleton
void assembleUVVolumePostSkeleton(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:297
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LocalAssembler
LA LocalAssembler
The type of the wrapping local assembler.
Definition: default/nonlinearjacobianapplyengine.hh:33
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LFSV
LA::LFSV LFSV
Definition: default/nonlinearjacobianapplyengine.hh:50
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onBindLFSVInside
void onBindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:167
Dune::PDELab::LocalVector::resize
void resize(size_type size)
Resize the container.
Definition: localvector.hh:257
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::Solution
LA::Traits::Solution Solution
The type of the solution vector.
Definition: default/nonlinearjacobianapplyengine.hh:43
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onUnbindLFSV
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:188
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::needsConstraintsCaching
bool needsConstraintsCaching(const TrialConstraintsContainer &cu, const TestConstraintsContainer &cv) const
Definition: default/nonlinearjacobianapplyengine.hh:27
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::loadCoefficientsLFSUOutside
void loadCoefficientsLFSUOutside(const LFSUC &lfsu_n_cache)
Definition: default/nonlinearjacobianapplyengine.hh:220
localvector.hh
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onBindLFSV
void onBindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:140
ig
const IG & ig
Definition: constraints.hh:148
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireUVSkeleton
bool requireUVSkeleton() const
Definition: default/nonlinearjacobianapplyengine.hh:77
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::ResidualElement
Residual::ElementType ResidualElement
Definition: default/nonlinearjacobianapplyengine.hh:40
Dune::PDELab::LocalAssemblerCallSwitch::nonlinear_jacobian_apply_skeleton
static void nonlinear_jacobian_apply_skeleton(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const X &z_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const X &z_n, const LFSV &lfsv_n, Y &y_s, Y &y_n)
Definition: callswitch.hh:163
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::GFSV
LFSV::Traits::GridFunctionSpace GFSV
Definition: default/nonlinearjacobianapplyengine.hh:52
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LFSUCache
LA::LFSUCache LFSUCache
Definition: default/nonlinearjacobianapplyengine.hh:48
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: default/nonlinearjacobianapplyengine.hh:113
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::SolutionView
Solution::template ConstLocalView< LFSUCache > SolutionView
Definition: default/nonlinearjacobianapplyengine.hh:54
Dune::PDELab::constrain_residual
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:906
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onUnbindLFSVOutside
void onUnbindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: default/nonlinearjacobianapplyengine.hh:202
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setUpdate
void setUpdate(const Solution &update_)
Definition: default/nonlinearjacobianapplyengine.hh:121
Dune::PDELab::TrialSpaceTag
Definition: localfunctionspacetags.hh:48
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::SolutionElement
Solution::ElementType SolutionElement
Definition: default/nonlinearjacobianapplyengine.hh:44
localassemblerenginebase.hh
constraints.hh
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LFSU
LA::LFSU LFSU
The local function spaces.
Definition: default/nonlinearjacobianapplyengine.hh:47
Dune::PDELab::LocalAssemblerEngineBase
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:21
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireUVBoundary
bool requireUVBoundary() const
Definition: default/nonlinearjacobianapplyengine.hh:79
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireSkeletonTwoSided
bool requireSkeletonTwoSided() const
Definition: default/nonlinearjacobianapplyengine.hh:73
Dune::PDELab::LocalAssemblerCallSwitch::nonlinear_jacobian_apply_volume
static void nonlinear_jacobian_apply_volume(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const X &z, const LFSV &lfsv, Y &y)
Definition: callswitch.hh:155
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::GFSU
LFSU::Traits::GridFunctionSpace GFSU
Definition: default/nonlinearjacobianapplyengine.hh:49
Dune::PDELab::LocalVector::assign
void assign(size_type size, const T &value)
Resize the container to size and assign the passed value to all entries.
Definition: localvector.hh:263
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::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: default/nonlinearjacobianapplyengine.hh:267
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::Residual
LA::Traits::Residual Residual
The type of the residual vector.
Definition: default/nonlinearjacobianapplyengine.hh:39
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::loadCoefficientsLFSUCoupling
void loadCoefficientsLFSUCoupling(const LFSUC &lfsu_c_cache)
Definition: default/nonlinearjacobianapplyengine.hh:226
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onBindLFSVOutside
void onBindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: default/nonlinearjacobianapplyengine.hh:174
Dune::PDELab::LocalAssemblerCallSwitch::nonlinear_jacobian_apply_boundary
static void nonlinear_jacobian_apply_boundary(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const X &z_s, const LFSV &lfsv_s, Y &y_s)
Definition: callswitch.hh:170
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireUVVolumePostSkeleton
bool requireUVVolumePostSkeleton() const
Definition: default/nonlinearjacobianapplyengine.hh:81
Dune::PDELab::LocalVector< ResidualElement, LocalTestSpaceTag >::WeightedAccumulationView
WeightedVectorAccumulationView< LocalVector > WeightedAccumulationView
An accumulate-only view of this container that automatically applies a weight to all contributions.
Definition: localvector.hh:198
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setResidual
void setResidual(Residual &residual_)
Definition: default/nonlinearjacobianapplyengine.hh:105
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onUnbindLFSVInside
void onUnbindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:195
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::LFSVCache
LA::LFSVCache LFSVCache
Definition: default/nonlinearjacobianapplyengine.hh:51
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::assembleCell
bool assembleCell(const EG &eg)
Definition: default/nonlinearjacobianapplyengine.hh:253
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::requireSkeleton
bool requireSkeleton() const
Definition: default/nonlinearjacobianapplyengine.hh:71
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::loadCoefficientsLFSUInside
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: default/nonlinearjacobianapplyengine.hh:214
Dune::PDELab::LocalAssemblerCallSwitch::nonlinear_jacobian_apply_volume_post_skeleton
static void nonlinear_jacobian_apply_volume_post_skeleton(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const X &z, const LFSV &lfsv, Y &y)
Definition: callswitch.hh:159
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::trialConstraints
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: default/nonlinearjacobianapplyengine.hh:92
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::onBindLFSUV
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/nonlinearjacobianapplyengine.hh:131