1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH 2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH 26 template<
typename TrialConstra
intsContainer,
typename TestConstra
intsContainer>
36 typedef typename LA::LocalOperator
LOP;
39 typedef typename LA::Traits::Residual
Residual;
43 typedef typename LA::Traits::Solution
Solution;
47 typedef typename LA::LFSU
LFSU;
49 typedef typename LFSU::Traits::GridFunctionSpace
GFSU;
50 typedef typename LA::LFSV
LFSV;
52 typedef typename LFSV::Traits::GridFunctionSpace
GFSV;
54 typedef typename Solution::template ConstLocalView<LFSUCache>
SolutionView;
55 typedef typename Residual::template LocalView<LFSVCache>
ResidualView;
64 : local_assembler(local_assembler_), lop(local_assembler_.localOperator()),
72 {
return local_assembler.doAlphaSkeleton(); }
74 {
return local_assembler.doSkeletonTwoSided(); }
76 {
return local_assembler.doAlphaVolume(); }
78 {
return local_assembler.doAlphaSkeleton(); }
80 {
return local_assembler.doAlphaBoundary(); }
82 {
return local_assembler.doAlphaVolumePostSkeleton(); }
88 return local_assembler;
92 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints&
trialConstraints()
const 98 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints&
testConstraints()
const 107 global_rl_view.attach(residual_);
108 global_rn_view.attach(residual_);
115 global_sl_view.attach(solution_);
116 global_sn_view.attach(solution_);
123 global_zl_view.attach(update_);
124 global_zn_view.attach(update_);
130 template<
typename EG,
typename LFSUC,
typename LFSVC>
131 void onBindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
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());
139 template<
typename EG,
typename LFSVC>
142 global_rl_view.bind(lfsv_cache);
143 rl.
assign(lfsv_cache.size(),0.0);
146 template<
typename IG,
typename LFSUC,
typename LFSVC>
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());
155 template<
typename IG,
typename LFSUC,
typename LFSVC>
157 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
158 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
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());
166 template<
typename IG,
typename LFSVC>
169 global_rl_view.bind(lfsv_cache);
170 rl.
assign(lfsv_cache.size(),0.0);
173 template<
typename IG,
typename LFSVC>
175 const LFSVC & lfsv_s_cache,
176 const LFSVC & lfsv_n_cache)
178 global_rn_view.bind(lfsv_n_cache);
179 rn.
assign(lfsv_n_cache.size(),0.0);
187 template<
typename EG,
typename LFSVC>
190 global_rl_view.add(rl);
191 global_rl_view.commit();
194 template<
typename IG,
typename LFSVC>
197 global_rl_view.add(rl);
198 global_rl_view.commit();
201 template<
typename IG,
typename LFSVC>
203 const LFSVC & lfsv_s_cache,
204 const LFSVC & lfsv_n_cache)
206 global_rn_view.add(rn);
207 global_rn_view.commit();
213 template<
typename LFSUC>
216 global_sl_view.read(xl);
217 global_zl_view.read(zl);
219 template<
typename LFSUC>
222 global_sn_view.read(xn);
223 global_zn_view.read(zn);
225 template<
typename LFSUC>
228 DUNE_THROW(Dune::NotImplemented,
"No coupling lfsu available for ");
237 if(local_assembler.doPostProcessing())
252 template<
typename EG>
255 return LocalAssembler::isNonOverlapping && eg.entity().partitionType() != Dune::InteriorEntity;
258 template<
typename EG,
typename LFSUC,
typename LFSVC>
261 rl_view.
setWeight(local_assembler.weight());
266 template<
typename IG,
typename LFSUC,
typename LFSVC>
268 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
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(),
279 template<
typename IG,
typename LFSUC,
typename LFSVC>
282 rl_view.
setWeight(local_assembler.weight());
287 template<
typename IG,
typename LFSUC,
typename LFSVC>
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)
293 DUNE_THROW(Dune::NotImplemented,
"Assembling of coupling spaces is not implemented for ");
296 template<
typename EG,
typename LFSUC,
typename LFSVC>
299 rl_view.
setWeight(local_assembler.weight());
309 const LocalAssembler & local_assembler;
315 ResidualView global_rl_view;
316 ResidualView global_rn_view;
319 SolutionView global_sl_view;
320 SolutionView global_sn_view;
323 SolutionView global_zl_view;
324 SolutionView global_zn_view;
357 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_NONLINEARJACOBIANAPPLYENGINE_HH An accumulate-only view on a local vector that automatically takes into account an accumulation weigh...
Definition: localvector.hh:26
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: nonlinearjacobianapplyengine.hh:86
const IG & ig
Definition: constraints.hh:148
bool requireUVVolume() const
Definition: nonlinearjacobianapplyengine.hh:75
void setUpdate(const Solution &update_)
Definition: nonlinearjacobianapplyengine.hh:121
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: nonlinearjacobianapplyengine.hh:92
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:131
Definition: localfunctionspacetags.hh:54
LA::LFSUCache LFSUCache
Definition: nonlinearjacobianapplyengine.hh:48
LA LocalAssembler
The type of the wrapping local assembler.
Definition: nonlinearjacobianapplyengine.hh:33
void loadCoefficientsLFSUCoupling(const LFSUC &lfsu_c_cache)
Definition: nonlinearjacobianapplyengine.hh:226
void assign(size_type size, const T &value)
Resize the container to size and assign the passed value to all entries.
Definition: localvector.hh:257
bool requireUVVolumePostSkeleton() const
Definition: nonlinearjacobianapplyengine.hh:81
DefaultLocalNonlinearJacobianApplyAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: nonlinearjacobianapplyengine.hh:63
LFSV::Traits::GridFunctionSpace GFSV
Definition: nonlinearjacobianapplyengine.hh:52
void onBindLFSUVInside(const IG &ig, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:147
void onUnbindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: nonlinearjacobianapplyengine.hh:202
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: nonlinearjacobianapplyengine.hh:267
bool requireSkeleton() const
Definition: nonlinearjacobianapplyengine.hh:71
LFSU::Traits::GridFunctionSpace GFSU
Definition: nonlinearjacobianapplyengine.hh:49
void setResidual(Residual &residual_)
Definition: nonlinearjacobianapplyengine.hh:105
Residual::template LocalView< LFSVCache > ResidualView
Definition: nonlinearjacobianapplyengine.hh:55
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
void onUnbindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:195
void assembleUVVolume(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:259
bool requireUVSkeleton() const
Definition: nonlinearjacobianapplyengine.hh:77
void resize(size_type size)
Resize the container.
Definition: localvector.hh:251
bool assembleCell(const EG &eg)
Definition: nonlinearjacobianapplyengine.hh:253
Residual::ElementType ResidualElement
Definition: nonlinearjacobianapplyengine.hh:40
LA::LFSV LFSV
Definition: nonlinearjacobianapplyengine.hh:50
LA::Traits::Solution Solution
The type of the solution vector.
Definition: nonlinearjacobianapplyengine.hh:43
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: nonlinearjacobianapplyengine.hh:156
void setSolution(const Solution &solution_)
Definition: nonlinearjacobianapplyengine.hh:113
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
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
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:906
bool requireSkeletonTwoSided() const
Definition: nonlinearjacobianapplyengine.hh:73
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
LA::LFSVCache LFSVCache
Definition: nonlinearjacobianapplyengine.hh:51
LA::LFSU LFSU
The local function spaces.
Definition: nonlinearjacobianapplyengine.hh:47
void assembleUVBoundary(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache)
Definition: nonlinearjacobianapplyengine.hh:280
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:188
LA::Traits::Residual Residual
The type of the residual vector.
Definition: nonlinearjacobianapplyengine.hh:39
The local assembler engine for DUNE grids which assembles the local application of the Jacobian...
Definition: nonlinearjacobianapplyengine.hh:21
void onBindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:140
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: nonlinearjacobianapplyengine.hh:235
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: nonlinearjacobianapplyengine.hh:214
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: nonlinearjacobianapplyengine.hh:98
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
Solution::template ConstLocalView< LFSUCache > SolutionView
Definition: nonlinearjacobianapplyengine.hh:54
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:21
Solution::ElementType SolutionElement
Definition: nonlinearjacobianapplyengine.hh:44
void onBindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:167
void setWeight(weight_type weight)
Resets the weighting coefficient of the view.
Definition: localvector.hh:72
void assembleUVVolumePostSkeleton(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: nonlinearjacobianapplyengine.hh:297
bool needsConstraintsCaching(const TrialConstraintsContainer &cu, const TestConstraintsContainer &cv) const
Definition: nonlinearjacobianapplyengine.hh:27
Definition: localfunctionspacetags.hh:48
LA::LocalOperator LOP
The type of the local operator.
Definition: nonlinearjacobianapplyengine.hh:36
void loadCoefficientsLFSUOutside(const LFSUC &lfsu_n_cache)
Definition: nonlinearjacobianapplyengine.hh:220
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: nonlinearjacobianapplyengine.hh:288
bool requireUVBoundary() const
Definition: nonlinearjacobianapplyengine.hh:79
void onBindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: nonlinearjacobianapplyengine.hh:174