1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_JACOBIANAPPLYENGINE_HH 2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_JACOBIANAPPLYENGINE_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_),
65 lop(local_assembler_.localOperator()),
73 {
return local_assembler.doAlphaSkeleton(); }
75 {
return local_assembler.doSkeletonTwoSided(); }
77 {
return local_assembler.doAlphaVolume(); }
79 {
return local_assembler.doAlphaSkeleton(); }
81 {
return local_assembler.doAlphaBoundary(); }
83 {
return local_assembler.doAlphaVolumePostSkeleton(); }
89 return local_assembler;
93 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints&
trialConstraints()
const 99 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints&
testConstraints()
const 108 global_rl_view.attach(residual_);
109 global_rn_view.attach(residual_);
116 global_sl_view.attach(solution_);
117 global_sn_view.attach(solution_);
123 template<
typename EG,
typename LFSUC,
typename LFSVC>
124 void onBindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
126 global_sl_view.bind(lfsu_cache);
127 xl.
resize(lfsu_cache.size());
130 template<
typename EG,
typename LFSVC>
133 global_rl_view.bind(lfsv_cache);
134 rl.
assign(lfsv_cache.size(),0.0);
137 template<
typename IG,
typename LFSUC,
typename LFSVC>
140 global_sl_view.bind(lfsu_cache);
141 xl.
resize(lfsu_cache.size());
144 template<
typename IG,
typename LFSUC,
typename LFSVC>
146 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
147 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
149 global_sn_view.bind(lfsu_n_cache);
150 xn.
resize(lfsu_n_cache.size());
153 template<
typename IG,
typename LFSVC>
156 global_rl_view.bind(lfsv_cache);
157 rl.
assign(lfsv_cache.size(),0.0);
160 template<
typename IG,
typename LFSVC>
162 const LFSVC & lfsv_s_cache,
163 const LFSVC & lfsv_n_cache)
165 global_rn_view.bind(lfsv_n_cache);
166 rn.
assign(lfsv_n_cache.size(),0.0);
174 template<
typename EG,
typename LFSVC>
177 global_rl_view.add(rl);
178 global_rl_view.commit();
181 template<
typename IG,
typename LFSVC>
184 global_rl_view.add(rl);
185 global_rl_view.commit();
188 template<
typename IG,
typename LFSVC>
190 const LFSVC & lfsv_s_cache,
191 const LFSVC & lfsv_n_cache)
193 global_rn_view.add(rn);
194 global_rn_view.commit();
200 template<
typename LFSUC>
203 global_sl_view.read(xl);
205 template<
typename LFSUC>
208 global_sn_view.read(xn);
210 template<
typename LFSUC>
213 DUNE_THROW(Dune::NotImplemented,
"No coupling lfsu available for ");
222 if(local_assembler.doPostProcessing())
224 global_rl_view.container());
239 template<
typename EG>
242 return LocalAssembler::isNonOverlapping && eg.entity().partitionType() != Dune::InteriorEntity;
245 template<
typename EG,
typename LFSUC,
typename LFSVC>
248 rl_view.
setWeight(local_assembler.weight());
250 jacobian_apply_volume(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),rl_view);
253 template<
typename IG,
typename LFSUC,
typename LFSVC>
255 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
257 rl_view.
setWeight(local_assembler.weight());
258 rn_view.
setWeight(local_assembler.weight());
261 lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),
262 lfsu_n_cache.localFunctionSpace(),xn,lfsv_n_cache.localFunctionSpace(),
266 template<
typename IG,
typename LFSUC,
typename LFSVC>
269 rl_view.
setWeight(local_assembler.weight());
271 jacobian_apply_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),rl_view);
274 template<
typename IG,
typename LFSUC,
typename LFSVC>
276 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
277 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache,
278 const LFSUC & lfsu_coupling_cache,
const LFSVC & lfsv_coupling_cache)
280 DUNE_THROW(Dune::NotImplemented,
"Assembling of coupling spaces is not implemented for ");
283 template<
typename EG,
typename LFSUC,
typename LFSVC>
286 rl_view.
setWeight(local_assembler.weight());
296 const LocalAssembler & local_assembler;
302 ResidualView global_rl_view;
303 ResidualView global_rn_view;
306 SolutionView global_sl_view;
307 SolutionView global_sn_view;
335 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_JACOBIANAPPLYENGINE_HH bool assembleCell(const EG &eg)
Definition: jacobianapplyengine.hh:240
An accumulate-only view on a local vector that automatically takes into account an accumulation weigh...
Definition: localvector.hh:26
const IG & ig
Definition: constraints.hh:148
void setSolution(const Solution &solution_)
Definition: jacobianapplyengine.hh:114
void setResidual(Residual &residual_)
Definition: jacobianapplyengine.hh:106
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: jacobianapplyengine.hh:145
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: jacobianapplyengine.hh:254
bool requireSkeletonTwoSided() const
Definition: jacobianapplyengine.hh:74
void onBindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:154
void onBindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: jacobianapplyengine.hh:161
Definition: localfunctionspacetags.hh:54
bool requireUVBoundary() const
Definition: jacobianapplyengine.hh:80
bool requireUVVolume() const
Definition: jacobianapplyengine.hh:76
void loadCoefficientsLFSUCoupling(const LFSUC &lfsu_c_cache)
Definition: jacobianapplyengine.hh:211
static void jacobian_apply_volume(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, Y &y)
Definition: callswitch.hh:132
static void jacobian_apply_volume_post_skeleton(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, Y &y)
Definition: callswitch.hh:136
LA::LocalOperator LOP
The type of the local operator.
Definition: jacobianapplyengine.hh:36
Solution::ElementType SolutionElement
Definition: jacobianapplyengine.hh:44
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
Residual::template LocalView< LFSVCache > ResidualView
Definition: jacobianapplyengine.hh:55
void onUnbindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: jacobianapplyengine.hh:189
bool needsConstraintsCaching(const TrialConstraintsContainer &cu, const TestConstraintsContainer &cv) const
Definition: jacobianapplyengine.hh:27
LA::LFSVCache LFSVCache
Definition: jacobianapplyengine.hh:51
void assembleUVBoundary(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache)
Definition: jacobianapplyengine.hh:267
LA::LFSU LFSU
The local function spaces.
Definition: jacobianapplyengine.hh:47
static void jacobian_apply_skeleton(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const LFSV &lfsv_n, Y &y_s, Y &y_n)
Definition: callswitch.hh:140
void assembleUVVolumePostSkeleton(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:284
LFSU::Traits::GridFunctionSpace GFSU
Definition: jacobianapplyengine.hh:49
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
void loadCoefficientsLFSUOutside(const LFSUC &lfsu_n_cache)
Definition: jacobianapplyengine.hh:206
LA LocalAssembler
The type of the wrapping local assembler.
Definition: jacobianapplyengine.hh:33
void resize(size_type size)
Resize the container.
Definition: localvector.hh:251
Solution::template ConstLocalView< LFSUCache > SolutionView
Definition: jacobianapplyengine.hh:54
void onBindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:131
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: jacobianapplyengine.hh:99
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:124
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:906
DefaultLocalJacobianApplyAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: jacobianapplyengine.hh:63
LA::Traits::Residual Residual
The type of the residual vector.
Definition: jacobianapplyengine.hh:39
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:175
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: jacobianapplyengine.hh:220
void assembleUVVolume(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:246
The local assembler engine for DUNE grids which assembles the local application of the Jacobian...
Definition: jacobianapplyengine.hh:21
void onBindLFSUVInside(const IG &ig, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:138
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:21
bool requireSkeleton() const
Definition: jacobianapplyengine.hh:72
LA::Traits::Solution Solution
The type of the solution vector.
Definition: jacobianapplyengine.hh:43
void setWeight(weight_type weight)
Resets the weighting coefficient of the view.
Definition: localvector.hh:72
void onUnbindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: jacobianapplyengine.hh:182
LFSV::Traits::GridFunctionSpace GFSV
Definition: jacobianapplyengine.hh:52
Residual::ElementType ResidualElement
Definition: jacobianapplyengine.hh:40
LA::LFSV LFSV
Definition: jacobianapplyengine.hh:50
Definition: localfunctionspacetags.hh:48
LA::LFSUCache LFSUCache
Definition: jacobianapplyengine.hh:48
static void jacobian_apply_boundary(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, Y &y_s)
Definition: callswitch.hh:147
bool requireUVSkeleton() const
Definition: jacobianapplyengine.hh:78
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: jacobianapplyengine.hh:87
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: jacobianapplyengine.hh:93
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: jacobianapplyengine.hh:201
bool requireUVVolumePostSkeleton() const
Definition: jacobianapplyengine.hh:82
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: jacobianapplyengine.hh:275