DUNE PDELab (git)

assemblerutilities.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLERUTILITIES_HH
5#define DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLERUTILITIES_HH
6
7#include <algorithm>
8#include <tuple>
9
10#include <dune/pdelab/constraints/common/constraintstransformation.hh>
11#include <dune/pdelab/gridoperator/common/localmatrix.hh>
12
13namespace Dune{
14 namespace PDELab{
15
16#ifndef DOXYGEN
17
18 namespace Impl {
19
20 // ********************************************************************************
21 // concept checks that test whether a local operator provides a given flag
22 // ********************************************************************************
23
24 struct HasDoSkipEntity
25 {
26 template<typename LO>
27 auto require(LO&& lo) -> decltype(
28 Concept::requireConvertible<bool>(LO::doSkipEntity)
29 );
30 };
31
32 struct HasDoSkipIntersection
33 {
34 template<typename LO>
35 auto require(LO&& lo) -> decltype(
36 Concept::requireConvertible<bool>(LO::doSkipIntersection)
37 );
38 };
39
40 } // namespace impl
41
42#endif // DOXYGEN
43
49 template<typename GO>
51 {
52
54 typedef typename GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace;
55
57 typedef typename GO::Traits::TestGridFunctionSpace TestGridFunctionSpace;
58
59
61 typedef typename GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints;
62
64 typedef typename GO::Traits::TestGridFunctionSpaceConstraints TestGridFunctionSpaceConstraints;
65
66
68 typedef typename GO::Traits::MatrixBackend MatrixBackend;
69
70
72 typedef typename GO::Traits::DomainField DomainField;
73
75 typedef typename GO::Traits::Domain Domain;
76
78 typedef typename GO::Traits::Domain Solution;
79
80
82 typedef typename GO::Traits::RangeField RangeField;
83
85 typedef typename GO::Traits::Range Range;
86
88 typedef typename GO::Traits::Range Residual;
89
90
92 typedef typename GO::Traits::JacobianField JacobianField;
93
95 typedef typename GO::Traits::Jacobian Jacobian;
96
98 typedef typename MatrixBackend::template Pattern<
103
105 typedef typename GO::BorderDOFExchanger BorderDOFExchanger;
106
107 };
108
109 // ********************************************************************************
110 // default local pattern implementation
111 // ********************************************************************************
112
125 class SparsityLink : public std::tuple<int,int>
126 {
127 public:
130 {}
131
133 SparsityLink (int i, int j)
134 : std::tuple<int,int>(i,j)
135 {}
136
138 inline int i () const
139 {
140 return std::get<0>(*this);
141 }
142
144 inline int j () const
145 {
146 return std::get<1>(*this);
147 }
148
150 void set (int i, int j)
151 {
152 std::get<0>(*this) = i;
153 std::get<1>(*this) = j;
154 }
155 };
156
164 : public std::vector<SparsityLink>
165 {
166
167 // make push_back() inaccessible
168 using std::vector<SparsityLink>::push_back;
169
170 public:
171
173
182 template<typename LFSV, typename LFSU>
183 void addLink(const LFSV& lfsv, std::size_t i, const LFSU& lfsu, std::size_t j)
184 {
187 lfsv.localIndex(i),
188 lfsu.localIndex(j)
189 )
190 );
191 }
192
193 };
194
195
196
197 // ********************************************************************************
198 // Assembler base class
199 // ********************************************************************************
200
213 template<typename B,
214 typename CU=EmptyTransformation,
215 typename CV=EmptyTransformation>
217 {
218 public:
219
220 typedef typename B::size_type SizeType;
221
224 : pconstraintsu(&emptyconstraintsu), pconstraintsv(&emptyconstraintsv)
225 {}
226
228 LocalAssemblerBase (const CU& cu, const CV& cv)
229 :pconstraintsu(&cu), pconstraintsv(&cv)
230 {}
231
233 const CU& trialConstraints() const
234 {
235 return *pconstraintsu;
236 }
237
239 const CV& testConstraints() const
240 {
241 return *pconstraintsv;
242 }
243
244
250 template<typename X>
251 typename std::enable_if<
252 AlwaysTrue<X>::value && !std::is_same<
253 CV,
254 EmptyTransformation
255 >::value
256 >::type
257 forwardtransform(X & x, const bool postrestrict = false) const
258 {
259 typedef typename CV::const_iterator global_col_iterator;
260 for (global_col_iterator cit=pconstraintsv->begin(); cit!=pconstraintsv->end(); ++cit){
261 typedef typename global_col_iterator::value_type::first_type GlobalIndex;
262 const GlobalIndex & contributor = cit->first;
263
264 typedef typename global_col_iterator::value_type::second_type ContributedMap;
265 typedef typename ContributedMap::const_iterator global_row_iterator;
266 const ContributedMap & contributed = cit->second;
267 global_row_iterator it = contributed.begin();
268 global_row_iterator eit = contributed.end();
269
270 for(;it!=eit;++it)
271 {
272 // typename X::block_type block(x[contributor]);
273 // block *= it->second;
274 // x[it->first] += block;
275 x[it->first] += it->second * x[contributor];
276 }
277 }
278
279 if(postrestrict)
280 for (global_col_iterator cit=pconstraintsv->begin(); cit!=pconstraintsv->end(); ++cit)
281 x[cit->first]=0.;
282 }
283
284
285 // Disable forwardtransform for EmptyTransformation
286 template<typename X>
287 typename std::enable_if<
288 AlwaysTrue<X>::value && std::is_same<
289 CV,
290 EmptyTransformation
291 >::value
292 >::type
293 forwardtransform(X & x, const bool postrestrict = false) const
294 {}
295
296
301 template<typename X>
302 typename std::enable_if<
303 AlwaysTrue<X>::value && !std::is_same<
304 CV,
305 EmptyTransformation
306 >::value
307 >::type
308 backtransform(X & x, const bool prerestrict = false) const
309 {
310 typedef typename CV::const_iterator global_col_iterator;
311 for (global_col_iterator cit=pconstraintsv->begin(); cit!=pconstraintsv->end(); ++cit){
312 typedef typename global_col_iterator::value_type::first_type GlobalIndex;
313 const GlobalIndex & contributor = cit->first;
314
315 typedef typename global_col_iterator::value_type::second_type ContributedMap;
316 typedef typename ContributedMap::const_iterator global_row_iterator;
317 const ContributedMap & contributed = cit->second;
318 global_row_iterator it = contributed.begin();
319 global_row_iterator eit = contributed.end();
320
321 if(prerestrict)
322 x[contributor] = 0.;
323
324 for(;it!=eit;++it)
325 {
326 // typename X::block_type block(x[it->first]);
327 // block *= it->second;
328 // x[contributor] += block;
329 x[contributor] += it->second * x[it->first]; // PB: 27 Sep 12 this was the old version
330 }
331 }
332 }
333
334 // disable backtransform for empty transformation
335 template<typename X>
336 typename std::enable_if<
337 AlwaysTrue<X>::value && std::is_same<
338 CV,
339 EmptyTransformation
340 >::value
341 >::type
342 backtransform(X & x, const bool prerestrict = false) const
343 {}
344
345
346 protected:
347
349 template<typename GCView, typename T>
350 void eread (const GCView& globalcontainer_view, LocalMatrix<T>& localcontainer) const
351 {
352 for (int i = 0; i < localcontainer.N(); ++i)
353 for (int j = 0; j < localcontainer.M(); ++j)
354 localcontainer(i,j) = globalcontainer_view(i,j);
355 }
356
358 template<typename T, typename GCView>
359 void ewrite (const LocalMatrix<T>& localcontainer, GCView& globalcontainer_view) const
360 {
361 for (int i = 0; i < localcontainer.N(); ++i)
362 for (int j = 0; j < localcontainer.M(); ++j)
363 globalcontainer_view(i,j) = localcontainer(i,j);
364 }
365
367 template<typename T, typename GCView>
368 void eadd (const LocalMatrix<T>& localcontainer, GCView& globalcontainer_view) const
369 {
370 for (int i = 0; i < localcontainer.N(); ++i)
371 for (int j = 0; j < localcontainer.M(); ++j)
372 globalcontainer_view.add(i,j,localcontainer(i,j));
373 }
374
376 template<typename M, typename GCView>
377 typename std::enable_if<
378 AlwaysTrue<M>::value && !std::is_same<
379 CV,
380 EmptyTransformation
381 >::value
382 >::type
383 scatter_jacobian(M& local_container, GCView& global_container_view, bool symmetric_mode) const
384 {
385 typedef typename GCView::RowIndexCache LFSVIndexCache;
386 typedef typename GCView::ColIndexCache LFSUIndexCache;
387
388 const LFSVIndexCache& lfsv_indices = global_container_view.rowIndexCache();
389 const LFSUIndexCache& lfsu_indices = global_container_view.colIndexCache();
390
391 if (lfsv_indices.constraintsCachingEnabled() && lfsu_indices.constraintsCachingEnabled())
392 if (symmetric_mode)
393 etadd_symmetric(local_container,global_container_view);
394 else
395 etadd(local_container,global_container_view);
396 else
397 {
398
399 typedef typename LFSVIndexCache::LocalFunctionSpace LFSV;
400 const LFSV& lfsv = lfsv_indices.localFunctionSpace();
401
402 typedef typename LFSUIndexCache::LocalFunctionSpace LFSU;
403 const LFSU& lfsu = lfsu_indices.localFunctionSpace();
404
405 // optionally clear out columns that belong to Dirichlet-constrained DOFs to keep matrix symmetric
406 if (symmetric_mode)
407 {
408 typedef typename LFSUIndexCache::ContainerIndex CI;
409
410 for (size_t j = 0; j < lfsu_indices.size(); ++j)
411 {
412 const CI& container_index = lfsu_indices.containerIndex(j);
413 const typename CU::const_iterator cit = pconstraintsu->find(container_index);
414 if (cit != pconstraintsu->end())
415 {
416 // make sure we only have Dirichlet constraints
417 assert(cit->second.empty());
418 // clear out the current column
419 for (size_t i = 0; i < lfsv_indices.size(); ++i)
420 {
421 // we do not need to update the residual, since the solution
422 // (i.e. the correction) for the Dirichlet DOF is 0 by definition
423 local_container(lfsv,i,lfsu,j) = 0.0;
424 }
425 }
426 }
427 }
428
429 // write entries without considering constraints.
430 // Dirichlet-constrained rows will be fixed in a postprocessing step.
431 for (auto it = local_container.begin(); it != local_container.end(); ++it)
432 {
433 // skip 0 entries because they might not be present in the pattern
434 if (*it == 0.0)
435 continue;
436 global_container_view.add(it.row(),it.col(),*it);
437 }
438 }
439 }
440
441 // specialization for empty constraints container
442 template<typename M, typename GCView>
443 typename std::enable_if<
444 AlwaysTrue<M>::value && std::is_same<
445 CV,
446 EmptyTransformation
447 >::value
448 >::type
449 scatter_jacobian(M& local_container, GCView& global_container_view, bool symmetric_mode) const
450 {
451 // write entries without considering constraints.
452 // Dirichlet-constrained rows will be fixed in a postprocessing step.
453 for (auto it = local_container.begin(); it != local_container.end(); ++it)
454 {
455 // skip 0 entries because they might not be present in the pattern
456 if (*it == 0.0)
457 continue;
458 global_container_view.add(it.row(),it.col(),*it);
459 }
460 }
461
465 template<typename M, typename GCView>
466 void etadd_symmetric (M& localcontainer, GCView& globalcontainer_view) const
467 {
468
469 typedef typename GCView::RowIndexCache LFSVIndexCache;
470 typedef typename GCView::ColIndexCache LFSUIndexCache;
471
472 const LFSVIndexCache& lfsv_indices = globalcontainer_view.rowIndexCache();
473 const LFSUIndexCache& lfsu_indices = globalcontainer_view.colIndexCache();
474
475 typedef typename LFSVIndexCache::LocalFunctionSpace LFSV;
476 const LFSV& lfsv = lfsv_indices.localFunctionSpace();
477
478 typedef typename LFSUIndexCache::LocalFunctionSpace LFSU;
479 const LFSU& lfsu = lfsu_indices.localFunctionSpace();
480
481 for (size_t j = 0; j < lfsu_indices.size(); ++j)
482 {
483 if (lfsu_indices.isConstrained(j) && lfsu_indices.isDirichletConstraint(j))
484 {
485 // clear out the current column
486 for (size_t i = 0; i < lfsv_indices.size(); ++i)
487 {
488 // we do not need to update the residual, since the solution
489 // (i.e. the correction) for the Dirichlet DOF is 0 by definition
490 localcontainer(lfsv,i,lfsu,j) = 0.0;
491 }
492 }
493 }
494
495 // hand off to standard etadd() method
496 etadd(localcontainer,globalcontainer_view);
497 }
498
499
500 template<typename M, typename GCView>
501 void etadd (const M& localcontainer, GCView& globalcontainer_view) const
502 {
503
504 typedef typename M::value_type T;
505
506 typedef typename GCView::RowIndexCache LFSVIndexCache;
507 typedef typename GCView::ColIndexCache LFSUIndexCache;
508
509 const LFSVIndexCache& lfsv_indices = globalcontainer_view.rowIndexCache();
510 const LFSUIndexCache& lfsu_indices = globalcontainer_view.colIndexCache();
511
512 typedef typename LFSVIndexCache::LocalFunctionSpace LFSV;
513 const LFSV& lfsv = lfsv_indices.localFunctionSpace();
514
515 typedef typename LFSUIndexCache::LocalFunctionSpace LFSU;
516 const LFSU& lfsu = lfsu_indices.localFunctionSpace();
517
518 for (size_t i = 0; i<lfsv_indices.size(); ++i)
519 for (size_t j = 0; j<lfsu_indices.size(); ++j)
520 {
521
522 if (localcontainer(lfsv,i,lfsu,j) == 0.0)
523 continue;
524
525 const bool constrained_v = lfsv_indices.isConstrained(i);
526 const bool constrained_u = lfsu_indices.isConstrained(j);
527
528 typedef typename LFSVIndexCache::ConstraintsIterator VConstraintsIterator;
529 typedef typename LFSUIndexCache::ConstraintsIterator UConstraintsIterator;
530
531 if (constrained_v)
532 {
533 if (lfsv_indices.isDirichletConstraint(i))
534 continue;
535
536 for (VConstraintsIterator vcit = lfsv_indices.constraintsBegin(i); vcit != lfsv_indices.constraintsEnd(i); ++vcit)
537 if (constrained_u)
538 {
539 if (lfsu_indices.isDirichletConstraint(j))
540 {
541 T value = localcontainer(lfsv,i,lfsu,j) * vcit->weight();
542 if (value != 0.0)
543 globalcontainer_view.add(vcit->containerIndex(),j,value);
544 }
545 else
546 {
547 for (UConstraintsIterator ucit = lfsu_indices.constraintsBegin(j); ucit != lfsu_indices.constraintsEnd(j); ++ucit)
548 {
549 T value = localcontainer(lfsv,i,lfsu,j) * vcit->weight() * ucit->weight();
550 if (value != 0.0)
551 globalcontainer_view.add(vcit->containerIndex(),ucit->containerIndex(),value);
552 }
553 }
554 }
555 else
556 {
557 T value = localcontainer(lfsv,i,lfsu,j) * vcit->weight();
558 if (value != 0.0)
559 globalcontainer_view.add(vcit->containerIndex(),j,value);
560 }
561 }
562 else
563 {
564 if (constrained_u)
565 {
566 if (lfsu_indices.isDirichletConstraint(j))
567 {
568 T value = localcontainer(lfsv,i,lfsu,j);
569 if (value != 0.0)
570 globalcontainer_view.add(i,j,value);
571 }
572 else
573 {
574 for (UConstraintsIterator ucit = lfsu_indices.constraintsBegin(j); ucit != lfsu_indices.constraintsEnd(j); ++ucit)
575 {
576 T value = localcontainer(lfsv,i,lfsu,j) * ucit->weight();
577 if (value != 0.0)
578 globalcontainer_view.add(i,ucit->containerIndex(),value);
579 }
580 }
581 }
582 else
583 globalcontainer_view.add(i,j,localcontainer(lfsv,i,lfsu,j));
584 }
585 }
586 }
587
588
589 template<typename Pattern, typename RI, typename CI>
590 typename std::enable_if<
591 std::is_same<RI,CI>::value
592 >::type
593 add_diagonal_entry(Pattern& pattern, const RI& ri, const CI& ci) const
594 {
595 if (ri == ci)
596 pattern.add_link(ri,ci);
597 }
598
599 template<typename Pattern, typename RI, typename CI>
600 typename std::enable_if<
601 !std::is_same<RI,CI>::value
602 >::type
603 add_diagonal_entry(Pattern& pattern, const RI& ri, const CI& ci) const
604 {}
605
610 template<typename P, typename LFSVIndices, typename LFSUIndices, typename Index>
611 void add_entry(P & globalpattern,
612 const LFSVIndices& lfsv_indices, Index i,
613 const LFSUIndices& lfsu_indices, Index j) const
614 {
615 typedef typename LFSVIndices::ConstraintsIterator VConstraintsIterator;
616 typedef typename LFSUIndices::ConstraintsIterator UConstraintsIterator;
617
618 const bool constrained_v = lfsv_indices.isConstrained(i);
619 const bool constrained_u = lfsu_indices.isConstrained(j);
620
621 add_diagonal_entry(globalpattern,
622 lfsv_indices.containerIndex(i),
623 lfsu_indices.containerIndex(j)
624 );
625
626 if(!constrained_v)
627 {
628 if (!constrained_u || lfsu_indices.isDirichletConstraint(j))
629 {
630 globalpattern.add_link(lfsv_indices.containerIndex(i),lfsu_indices.containerIndex(j));
631 }
632 else
633 {
634 for (UConstraintsIterator gurit = lfsu_indices.constraintsBegin(j); gurit != lfsu_indices.constraintsEnd(j); ++gurit)
635 globalpattern.add_link(lfsv_indices.containerIndex(i),gurit->containerIndex());
636 }
637 }
638 else
639 {
640 if (lfsv_indices.isDirichletConstraint(i))
641 {
642 globalpattern.add_link(lfsv_indices.containerIndex(i),lfsu_indices.containerIndex(j));
643 }
644 else
645 {
646 for(VConstraintsIterator gvrit = lfsv_indices.constraintsBegin(i); gvrit != lfsv_indices.constraintsEnd(i); ++gvrit)
647 {
648 if (!constrained_u || lfsu_indices.isDirichletConstraint(j))
649 {
650 globalpattern.add_link(gvrit->containerIndex(),lfsu_indices.containerIndex(j));
651 }
652 else
653 {
654 for (UConstraintsIterator gurit = lfsu_indices.constraintsBegin(j); gurit != lfsu_indices.constraintsEnd(j); ++gurit)
655 globalpattern.add_link(gvrit->containerIndex(),gurit->containerIndex());
656 }
657 }
658 }
659 }
660 }
661
665 template<typename GFSV, typename GC, typename C>
666 void set_trivial_rows(const GFSV& gfsv, GC& globalcontainer, const C& c) const
667 {
668 typedef typename C::const_iterator global_row_iterator;
669 for (global_row_iterator cit = c.begin(); cit != c.end(); ++cit)
670 globalcontainer.clear_row(cit->first,1);
671 }
672
673 template<typename GFSV, typename GC>
674 void set_trivial_rows(const GFSV& gfsv, GC& globalcontainer, const EmptyTransformation& c) const
675 {
676 }
677
678 template<typename GFSV, typename GC>
679 void handle_dirichlet_constraints(const GFSV& gfsv, GC& globalcontainer) const
680 {
681 globalcontainer.flush();
682 set_trivial_rows(gfsv,globalcontainer,*pconstraintsv);
683 globalcontainer.finalize();
684 }
685
686 /* constraints */
687 const CU* pconstraintsu;
688 const CV* pconstraintsv;
689 static CU emptyconstraintsu;
690 static CV emptyconstraintsv;
691 };
692
693 template<typename B, typename CU, typename CV>
694 CU LocalAssemblerBase<B,CU,CV>::emptyconstraintsu;
695 template<typename B, typename CU, typename CV>
696 CV LocalAssemblerBase<B,CU,CV>::emptyconstraintsv;
697
698 } // namespace PDELab
699} // namespace Dune
700
701#endif // DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLERUTILITIES_HH
Base class for local assembler.
Definition: assemblerutilities.hh:217
std::enable_if< AlwaysTrue< M >::value &&!std::is_same< CV, EmptyTransformation >::value >::type scatter_jacobian(M &local_container, GCView &global_container_view, bool symmetric_mode) const
Scatter local jacobian to global container.
Definition: assemblerutilities.hh:383
void set_trivial_rows(const GFSV &gfsv, GC &globalcontainer, const C &c) const
insert dirichlet constraints for row and assemble T^T_U in constrained rows
Definition: assemblerutilities.hh:666
void add_entry(P &globalpattern, const LFSVIndices &lfsv_indices, Index i, const LFSUIndices &lfsu_indices, Index j) const
Adding matrix entry to pattern with respect to the constraints contributions. This assembles the entr...
Definition: assemblerutilities.hh:611
void etadd_symmetric(M &localcontainer, GCView &globalcontainer_view) const
Add local matrix to global matrix, and apply Dirichlet constraints in a symmetric fashion....
Definition: assemblerutilities.hh:466
std::enable_if< AlwaysTrue< X >::value &&!std::is_same< CV, EmptyTransformation >::value >::type backtransform(X &x, const bool prerestrict=false) const
Transforms a vector from to . If prerestrict == true then is applied instead of the full transform...
Definition: assemblerutilities.hh:308
LocalAssemblerBase(const CU &cu, const CV &cv)
construct AssemblerSpace, with constraints
Definition: assemblerutilities.hh:228
void eread(const GCView &globalcontainer_view, LocalMatrix< T > &localcontainer) const
read local stiffness matrix for entity
Definition: assemblerutilities.hh:350
std::enable_if< AlwaysTrue< X >::value &&!std::is_same< CV, EmptyTransformation >::value >::type forwardtransform(X &x, const bool postrestrict=false) const
Transforms a vector from to . If postrestrict == true then is applied instead of the full transfor...
Definition: assemblerutilities.hh:257
const CU & trialConstraints() const
get the constraints on the trial grid function space
Definition: assemblerutilities.hh:233
void eadd(const LocalMatrix< T > &localcontainer, GCView &globalcontainer_view) const
write local stiffness matrix for entity
Definition: assemblerutilities.hh:368
void ewrite(const LocalMatrix< T > &localcontainer, GCView &globalcontainer_view) const
write local stiffness matrix for entity
Definition: assemblerutilities.hh:359
LocalAssemblerBase()
construct AssemblerSpace
Definition: assemblerutilities.hh:223
const CV & testConstraints() const
get the constraints on the test grid function space
Definition: assemblerutilities.hh:239
A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpac...
Definition: localmatrix.hh:184
Layout description for a sparse linear operator.
Definition: assemblerutilities.hh:165
void addLink(const LFSV &lfsv, std::size_t i, const LFSU &lfsu, std::size_t j)
Adds a link between DOF i of lfsv and DOF j of lfsu.
Definition: assemblerutilities.hh:183
constexpr HybridTreePath< T..., std::size_t > push_back(const HybridTreePath< T... > &tp, std::size_t i)
Appends a run time index to a HybridTreePath.
Definition: treepath.hh:416
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integer_sequence< T, II..., T(IN)> push_back(std::integer_sequence< T, II... >, std::integral_constant< T, IN >={})
Append an index IN to the back of the sequence.
Definition: integersequence.hh:69
STL namespace.
template which always yields a true value
Definition: typetraits.hh:134
Definition: assemblerutilities.hh:51
GO::Traits::Range Residual
The type of the range (residual).
Definition: assemblerutilities.hh:88
MatrixBackend::template Pattern< Jacobian, TestGridFunctionSpace, TrialGridFunctionSpace > MatrixPattern
The matrix pattern.
Definition: assemblerutilities.hh:102
GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace
The trial grid function space.
Definition: assemblerutilities.hh:54
GO::Traits::Jacobian Jacobian
The type of the jacobian.
Definition: assemblerutilities.hh:95
GO::Traits::DomainField DomainField
The field type of the domain (solution).
Definition: assemblerutilities.hh:72
GO::Traits::JacobianField JacobianField
The field type of the jacobian.
Definition: assemblerutilities.hh:92
GO::Traits::TestGridFunctionSpace TestGridFunctionSpace
The test grid function space.
Definition: assemblerutilities.hh:57
GO::BorderDOFExchanger BorderDOFExchanger
The helper class to exchange data on the processor boundary.
Definition: assemblerutilities.hh:105
GO::Traits::MatrixBackend MatrixBackend
The matrix backend of the grid operator.
Definition: assemblerutilities.hh:68
GO::Traits::Domain Solution
The type of the domain (solution).
Definition: assemblerutilities.hh:78
GO::Traits::RangeField RangeField
The field type of the range (residual).
Definition: assemblerutilities.hh:82
GO::Traits::Range Range
The type of the range (residual).
Definition: assemblerutilities.hh:85
GO::Traits::Domain Domain
The type of the domain (solution).
Definition: assemblerutilities.hh:75
GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints
The type of the trial grid function space constraints.
Definition: assemblerutilities.hh:61
GO::Traits::TestGridFunctionSpaceConstraints TestGridFunctionSpaceConstraints
The type of the test grid function space constraints.
Definition: assemblerutilities.hh:64
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)