DUNE PDELab (git)

pdelab.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3#ifndef DUNE_PDELAB_BOILERPLATE_PDELAB_HH
4#define DUNE_PDELAB_BOILERPLATE_PDELAB_HH
5
16// first of all we include a lot of dune grids and pdelab files
17#include <array>
18#include <iostream>
19#include <memory>
20
21#include <dune/common/parallel/mpihelper.hh> // include mpi helper class
27
28#include <dune/geometry/type.hh>
30
31#include <dune/grid/onedgrid.hh>
33#include <dune/grid/yaspgrid.hh>
34#if HAVE_DUNE_UGGRID
35#include <dune/grid/uggrid.hh>
36#endif
37#if HAVE_ALBERTA
38#include<dune/grid/albertagrid.hh>
39#include <dune/grid/albertagrid/dgfparser.hh>
40#endif
41#if HAVE_DUNE_UGGRID
42#include<dune/grid/uggrid.hh>
43#endif
44#if HAVE_DUNE_ALUGRID
45#include<dune/alugrid/grid.hh>
46#include <dune/alugrid/dgf.hh>
47#endif
49#include <dune/grid/io/file/gmshreader.hh>
50
51#include <dune/istl/bvector.hh>
53#include <dune/istl/solvers.hh>
54#include <dune/istl/solvercategory.hh>
56#include <dune/istl/io.hh>
57
59#include <dune/pdelab/common/function.hh>
60#include <dune/pdelab/common/functionutilities.hh>
61#include <dune/pdelab/common/topologyutility.hh>
62#include <dune/pdelab/common/vtkexport.hh>
63#include <dune/pdelab/backend/istl.hh>
64#include <dune/pdelab/constraints/conforming.hh>
65#include <dune/pdelab/constraints/hangingnode.hh>
66#include <dune/pdelab/constraints/p0.hh>
67#include <dune/pdelab/constraints/p0ghost.hh>
68#include <dune/pdelab/constraints/common/constraints.hh>
69#include <dune/pdelab/gridfunctionspace/gridfunctionspace.hh>
70#include <dune/pdelab/gridfunctionspace/gridfunctionspaceutilities.hh>
71#include <dune/pdelab/gridfunctionspace/interpolate.hh>
72#include <dune/pdelab/gridoperator/gridoperator.hh>
73#include <dune/pdelab/gridoperator/onestep.hh>
74#include <dune/pdelab/stationary/linearproblem.hh>
75#include <dune/pdelab/finiteelementmap/pkfem.hh>
76#include <dune/pdelab/finiteelementmap/p0fem.hh>
77#include <dune/pdelab/finiteelementmap/opbfem.hh>
78#include <dune/pdelab/finiteelementmap/qkfem.hh>
79#include <dune/pdelab/finiteelementmap/qkdg.hh>
80#include <dune/pdelab/adaptivity/adaptivity.hh>
81#include <dune/pdelab/instationary/onestep.hh>
82#include <dune/pdelab/common/instationaryfilenamehelper.hh>
83#include <dune/pdelab/newton/newton.hh>
84
85namespace Dune {
86 namespace PDELab {
87
88 // make grids
89 template<typename T>
90 class StructuredGrid
91 {
92 public:
93 // export types
94 typedef T Grid;
95 typedef typename T::ctype ctype;
96 static const int dim = T::dimension;
97 static const int dimworld = T::dimensionworld;
98
99 // constructors
100 StructuredGrid (Dune::GeometryType::BasicType meshtype, unsigned int cells)
101 {
102 FieldVector<ctype,dimworld> lowerLeft(0.0);
103 FieldVector<ctype,dimworld> upperRight(1.0);
104 std::array<unsigned int,dim> elements; elements.fill(cells);
105
106 StructuredGridFactory<T> factory;
107
108 if (meshtype==Dune::GeometryType::cube)
109 gridp = factory.createCubeGrid(lowerLeft,upperRight,elements);
110 else if (meshtype==Dune::GeometryType::simplex)
111 gridp = factory.createSimplexGrid(lowerLeft,upperRight,elements);
112 else
113 {
114 DUNE_THROW(GridError, className<StructuredGrid>()
115 << "::StructuredGrid(): grid type must be simplex or cube ");
116 }
117 }
118
119
120 StructuredGrid (Dune::GeometryType::BasicType meshtype,
121 std::array<double,dimworld> lower_left, std::array<double,dimworld> upper_right,
122 std::array<unsigned int,dim> cells)
123 {
124 FieldVector<ctype,dimworld> lowerLeft;
125 FieldVector<ctype,dimworld> upperRight;
126 std::array<unsigned int,dim> elements;
127
128 // copy data to correct types for StructuredGridFactory
129 for (size_t i=0; i<dimworld; i++)
130 {
131 lowerLeft[i] = lower_left[i];
132 upperRight[i] = upper_right[i];
133 }
134 for (size_t i=0; i<dim; i++)
135 {
136 elements[i] = cells[i];
137 }
138
139 StructuredGridFactory<T> factory;
140
141 if (meshtype==Dune::GeometryType::cube)
142 gridp = factory.createCubeGrid(lowerLeft,upperRight,elements);
143 else if (meshtype==Dune::GeometryType::simplex)
144 gridp = factory.createSimplexGrid(lowerLeft,upperRight,elements);
145 else
146 {
147 DUNE_THROW(GridError, className<StructuredGrid>()
148 << "::StructuredGrid(): grid type must be simplex or cube ");
149 }
150 }
151
152 // return shared pointer
153 std::shared_ptr<T> getSharedPtr ()
154 {
155 return gridp;
156 }
157
158 // return grid reference
159 T& getGrid ()
160 {
161 return *gridp;
162 }
163
164 // return grid reference const version
165 const T& getGrid () const
166 {
167 return *gridp;
168 }
169
170 T& operator*()
171 {
172 return *gridp;
173 }
174
175 T* operator->()
176 {
177 return gridp.operator->();
178 }
179
180 const T& operator*() const
181 {
182 return *gridp;
183 }
184
185 const T* operator->() const
186 {
187 return gridp.operator->();
188 }
189
190
191 private:
192 std::shared_ptr<T> gridp; // hold a shared pointer to a grid
193 };
194
195 // specialization for yaspgrid; treats paralle case right
196 template<int dim>
197 class StructuredGrid<YaspGrid<dim> >
198 {
199 public:
200
201 // export types
202 typedef YaspGrid<dim> Grid;
203 typedef typename Grid::ctype ctype;
204 static const int dimworld = Grid::dimensionworld;
205
206 // simple constructor for the unit cube
207 StructuredGrid (Dune::GeometryType::BasicType meshtype, unsigned int cells, int overlap=1)
208 {
209 // check element type
210 if (meshtype!=Dune::GeometryType::cube)
211 std::cout << "StructuredGrid(): element type " << meshtype << " is ignored" << std::endl;
212
213 // copy data to correct types for YaspGrid
215 std::array<int,dimworld> N(Dune::filledArray<dimworld, int>(cells));
216 std::bitset<dimworld> B(false);
217
218 // instantiate the grid
219 gridp = std::shared_ptr<Grid>(new Grid(L,N,B,overlap,Dune::MPIHelper::getCommunication()));
220 }
221
222 // constructor with sizes given
223 StructuredGrid (Dune::GeometryType::BasicType meshtype,
224 std::array<double,dimworld> lower_left, std::array<double,dimworld> upper_right,
225 std::array<unsigned int,dim> cells, int overlap=1)
226 {
227 // check that lower right corner is the origin
228 for(int d = 0; d < dimworld; ++d)
229 if(std::abs(lower_left[d]) > std::abs(upper_right[d])*1e-10)
230 DUNE_THROW(GridError, className<StructuredGrid>()
231 << "::createCubeGrid(): The lower coordinates "
232 "must be at the origin for YaspGrid.");
233
234 // check element type
235 if (meshtype!=Dune::GeometryType::cube)
236 std::cout << "StructuredGrid(): element type " << meshtype << " is ignored" << std::endl;
237
238 // copy data to correct types for YaspGrid
240 std::array<int,dimworld> N;
241 std::bitset<dimworld> B(false);
242 for (size_t i=0; i<dimworld; i++)
243 {
244 L[i] = upper_right[i];
245 N[i] = cells[i];
246 }
247
248 // instantiate the grid
249 gridp = std::shared_ptr<Grid>(new Grid(L,N,B,overlap,Dune::MPIHelper::getCommunication()));
250 }
251
252 // constructor with periodicity argument
253 StructuredGrid (Dune::GeometryType::BasicType meshtype,
254 std::array<double,dimworld> lower_left, std::array<double,dimworld> upper_right,
255 std::array<unsigned int,dim> cells, std::array<bool,dim> periodic, int overlap=1)
256 {
257 // check that lower right corner is the origin
258 for(int d = 0; d < dimworld; ++d)
259 if(std::abs(lower_left[d]) > std::abs(upper_right[d])*1e-10)
260 DUNE_THROW(GridError, className<StructuredGrid>()
261 << "::createCubeGrid(): The lower coordinates "
262 "must be at the origin for YaspGrid.");
263
264 // check element type
265 if (meshtype!=Dune::GeometryType::cube)
266 std::cout << "StructuredGrid(): element type " << meshtype << " is ignored" << std::endl;
267
268 // copy data to correct types for YaspGrid
270 std::array<int,dimworld> N;
271 std::bitset<dimworld> B(false);
272 for (size_t i=0; i<dimworld; i++)
273 {
274 L[i] = upper_right[i];
275 N[i] = cells[i];
276 B[i] = periodic[i];
277 }
278
279 // instantiate the grid
280 gridp = std::shared_ptr<Grid>(new Grid(L,N,B,overlap,Dune::MPIHelper::getCommunication()));
281 }
282
283 // return shared pointer
284 std::shared_ptr<Grid> getSharedPtr ()
285 {
286 return gridp;
287 }
288
289 // return grid reference
290 Grid& getGrid ()
291 {
292 return *gridp;
293 }
294
295 // return grid reference const version
296 const Grid& getGrid () const
297 {
298 return *gridp;
299 }
300
301 Grid& operator*()
302 {
303 return *gridp;
304 }
305
306 Grid* operator->()
307 {
308 return gridp.operator->();
309 }
310
311 const Grid& operator*() const
312 {
313 return *gridp;
314 }
315
316 const Grid* operator->() const
317 {
318 return gridp.operator->();
319 }
320
321 private:
322 std::shared_ptr<Grid> gridp; // hold a shared pointer to a grid
323 };
324
325 // unstructured grid read from gmsh file
326 template<typename T>
327 class UnstructuredGrid
328 {
329 public:
330 // export types
331 typedef T Grid;
332 typedef typename T::ctype ctype;
333 static const int dim = T::dimension;
334 static const int dimworld = T::dimensionworld;
335
336 // constructors
337 UnstructuredGrid (std::string filename, bool verbose = true, bool insert_boundary_segments=true)
338 {
339 Dune::GridFactory<T> factory;
340 Dune::GmshReader<T>::read(factory,filename,verbose,insert_boundary_segments);
341 gridp = std::shared_ptr<T>(factory.createGrid());
342 }
343
344 // return shared pointer
345 std::shared_ptr<T> getSharedPtr ()
346 {
347 return gridp;
348 }
349
350 // return grid reference
351 T& getGrid ()
352 {
353 return *gridp;
354 }
355
356 // return grid reference const version
357 const T& getGrid () const
358 {
359 return *gridp;
360 }
361
362 T& operator*()
363 {
364 return *gridp;
365 }
366
367 T* operator->()
368 {
369 return gridp.operator->();
370 }
371
372 const T& operator*() const
373 {
374 return *gridp;
375 }
376
377 const T* operator->() const
378 {
379 return gridp.operator->();
380 }
381
382 private:
383 std::shared_ptr<T> gridp; // hold a shared pointer to a grid
384 };
385
386
387 //============================================================================
388 // Continuous Lagrange Finite Element Space
389 //============================================================================
390
391 // finite element map base template
392 template<typename GV, typename C, typename R, unsigned int degree, unsigned int dim, Dune::GeometryType::BasicType gt>
393 class CGFEMBase
394 {};
395
396 template<typename GV, typename C, typename R, unsigned int degree, unsigned int dim>
397 class CGFEMBase<GV,C,R,degree,dim,Dune::GeometryType::simplex>
398 {
399 public:
400 typedef PkLocalFiniteElementMap<GV,C,R,degree> FEM;
401
402 CGFEMBase (const GV& gridview)
403 {
404 femp = std::shared_ptr<FEM>(new FEM(gridview));
405 }
406
407 FEM& getFEM() {return *femp;}
408 const FEM& getFEM() const {return *femp;}
409
410 private:
411 std::shared_ptr<FEM> femp;
412 };
413
414 template<typename GV, typename C, typename R, unsigned int degree, unsigned int dim>
415 class CGFEMBase<GV,C,R,degree,dim,Dune::GeometryType::cube>
416 {
417 public:
418 typedef QkLocalFiniteElementMap<GV,C,R,degree> FEM;
419
420 CGFEMBase (const GV& gridview)
421 {
422 femp = std::shared_ptr<FEM>(new FEM(gridview));
423 }
424
425 FEM& getFEM() {return *femp;}
426 const FEM& getFEM() const {return *femp;}
427
428 private:
429 std::shared_ptr<FEM> femp;
430 };
431
432 //============================================================================
433
434 // define enumeration type that differentiate conforming and nonconforming meshes
435 enum MeshType {
438 };
439
440 // constraints base template
441 template<typename Grid, unsigned int degree, Dune::GeometryType::BasicType gt, MeshType mt, SolverCategory::Category st, typename BCType, typename GV = typename Grid::LeafGridView>
442 class CGCONBase
443 {};
444
445 template<typename Grid, typename BCType, typename GV>
446 class CGCONBase<Grid,1,Dune::GeometryType::simplex,MeshType::nonconforming,SolverCategory::sequential,BCType,GV>
447 {
448 public:
449 typedef HangingNodesDirichletConstraints<Grid,HangingNodesConstraintsAssemblers::SimplexGridP1Assembler,BCType> CON;
450
451 CGCONBase (Grid& grid, const BCType& bctype, const GV& gv)
452 {
453 conp = std::shared_ptr<CON>(new CON(grid,true,bctype));
454 }
455
456 CGCONBase (Grid& grid, const BCType& bctype)
457 {
458 conp = std::shared_ptr<CON>(new CON(grid,true,bctype));
459 }
460
461 template<typename GFS>
462 void postGFSHook (const GFS& gfs) {}
463 CON& getCON() {return *conp;}
464 const CON& getCON() const {return *conp;}
465 template<typename GFS, typename DOF>
466 void make_consistent (const GFS& gfs, DOF& x) const {}
467 private:
468 std::shared_ptr<CON> conp;
469 };
470
471 template<typename Grid, typename BCType, typename GV>
472 class CGCONBase<Grid,1,Dune::GeometryType::cube,MeshType::nonconforming,SolverCategory::sequential,BCType,GV>
473 {
474 public:
475 typedef HangingNodesDirichletConstraints<Grid,HangingNodesConstraintsAssemblers::CubeGridQ1Assembler,BCType> CON;
476
477 CGCONBase (Grid& grid, const BCType& bctype, const GV& gv)
478 {
479 conp = std::shared_ptr<CON>(new CON(grid,true,bctype));
480 }
481
482 CGCONBase (Grid& grid, const BCType& bctype)
483 {
484 conp = std::shared_ptr<CON>(new CON(grid,true,bctype));
485 }
486
487 template<typename GFS>
488 void postGFSHook (const GFS& gfs) {}
489 CON& getCON() {return *conp;}
490 const CON& getCON() const {return *conp;}
491 template<typename GFS, typename DOF>
492 void make_consistent (const GFS& gfs, DOF& x) const {}
493 private:
494 std::shared_ptr<CON> conp;
495 };
496
497 template<typename Grid, unsigned int degree, Dune::GeometryType::BasicType gt,typename BCType, typename GV>
498 class CGCONBase<Grid,degree,gt,MeshType::conforming,SolverCategory::sequential,BCType,GV>
499 {
500 public:
501 typedef ConformingDirichletConstraints CON;
502
503 CGCONBase (Grid& grid, const BCType& bctype, const GV& gv)
504 {
505 conp = std::shared_ptr<CON>(new CON());
506 }
507
508 CGCONBase (Grid& grid, const BCType& bctype)
509 {
510 conp = std::shared_ptr<CON>(new CON());
511 }
512
513 template<typename GFS>
514 void postGFSHook (const GFS& gfs) {}
515 CON& getCON() {return *conp;}
516 const CON& getCON() const {return *conp;}
517 template<typename GFS, typename DOF>
518 void make_consistent (const GFS& gfs, DOF& x) const {}
519 private:
520 std::shared_ptr<CON> conp;
521 };
522
523 template<typename Grid, unsigned int degree, Dune::GeometryType::BasicType gt,typename BCType, typename GV>
524 class CGCONBase<Grid,degree,gt,MeshType::conforming,SolverCategory::overlapping,BCType,GV>
525 {
526 public:
527 typedef OverlappingConformingDirichletConstraints CON;
528
529 CGCONBase (Grid& grid, const BCType& bctype, const GV& gv)
530 {
531 conp = std::shared_ptr<CON>(new CON());
532 }
533
534 CGCONBase (Grid& grid, const BCType& bctype)
535 {
536 conp = std::shared_ptr<CON>(new CON());
537 }
538
539 template<typename GFS>
540 void postGFSHook (const GFS& gfs) {}
541 CON& getCON() {return *conp;}
542 const CON& getCON() const {return *conp;}
543 template<typename GFS, typename DOF>
544 void make_consistent (const GFS& gfs, DOF& x) const
545 {
546 // make vector consistent; this is needed for all overlapping solvers
547 ISTL::ParallelHelper<GFS> helper(gfs);
548 helper.maskForeignDOFs(Backend::native(x));
549 Dune::PDELab::AddDataHandle<GFS,DOF> adddh(gfs,x);
550 if (gfs.gridView().comm().size()>1)
551 gfs.gridView().communicate(adddh,Dune::InteriorBorder_All_Interface,Dune::ForwardCommunication);
552 }
553 private:
554 std::shared_ptr<CON> conp;
555 };
556
557 template<typename Grid, unsigned int degree, Dune::GeometryType::BasicType gt,typename BCType, typename GV>
558 class CGCONBase<Grid,degree,gt,MeshType::conforming,SolverCategory::nonoverlapping,BCType,GV>
559 {
560 public:
562 CGCONBase (Grid& grid, const BCType& bctype)
563 {
564 conp = std::shared_ptr<CON>(new CON);
565 }
566
567 template<typename GFS>
568 CON& getCON() {return *conp;}
569 const CON& getCON() const {return *conp;}
570 template<typename GFS, typename DOF>
571 void make_consistent (const GFS& gfs, DOF& x) const {}
572 private:
573 std::shared_ptr<CON> conp;
574 };
575
576
577 // continuous Lagrange finite elements
578 template<typename T, typename N, unsigned int degree, typename BCType,
580 typename VBET=ISTL::VectorBackend<> >
581 class CGSpace {
582 public:
583
584 // export types
585 typedef T Grid;
586 typedef typename T::LeafGridView GV;
587 typedef typename T::ctype ctype;
588 static const int dim = T::dimension;
589 static const int dimworld = T::dimensionworld;
590
591 typedef CGFEMBase<GV,ctype,N,degree,dim,gt> FEMB;
592 typedef CGCONBase<Grid,degree,gt,mt,st,BCType> CONB;
593
594 typedef typename FEMB::FEM FEM;
595 typedef typename CONB::CON CON;
596
597 typedef VBET VBE;
598 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
599
600 typedef N NT;
601 using DOF = Backend::Vector<GFS,N>;
603 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
604 typedef VTKGridFunctionAdapter<DGF> VTKF;
605
606 // constructor making the grid function space an all that is needed
607 CGSpace (Grid& grid, const BCType& bctype)
608 : gv(grid.leafGridView()), femb(gv), conb(grid,bctype)
609 {
610 gfsp = std::shared_ptr<GFS>(new GFS(gv,femb.getFEM(),conb.getCON()));
611 gfsp->name("cgspace");
612 // initialize ordering
613 gfsp->update();
614 conb.postGFSHook(*gfsp);
615 ccp = std::shared_ptr<CC>(new CC());
616 }
617
618 FEM& getFEM()
619 {
620 return femb.getFEM();
621 }
622
623 const FEM& getFEM() const
624 {
625 return femb.getFEM();
626 }
627
628 // return gfs reference
629 GFS& getGFS ()
630 {
631 return *gfsp;
632 }
633
634 // return gfs reference const version
635 const GFS& getGFS () const
636 {
637 return *gfsp;
638 }
639
640 // return gfs reference
641 CC& getCC ()
642 {
643 return *ccp;
644 }
645
646 // return gfs reference const version
647 const CC& getCC () const
648 {
649 return *ccp;
650 }
651
652 void assembleConstraints (const BCType& bctype)
653 {
654 ccp->clear();
655 constraints(bctype,*gfsp,*ccp);
656 }
657
658 void clearConstraints ()
659 {
660 ccp->clear();
661 }
662
663 void setConstrainedDOFS (DOF& x, NT nt) const
664 {
665 set_constrained_dofs(*ccp,nt,x);
666 conb.make_consistent(*gfsp,x);
667 }
668
669 void setNonConstrainedDOFS (DOF& x, NT nt) const
670 {
671 set_nonconstrained_dofs(*ccp,nt,x);
672 conb.make_consistent(*gfsp,x);
673 }
674
675 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
676 {
677 copy_constrained_dofs(*ccp,xin,xout);
678 conb.make_consistent(*gfsp,xout);
679 }
680
681 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
682 {
683 copy_nonconstrained_dofs(*ccp,xin,xout);
684 conb.make_consistent(*gfsp,xout);
685 }
686
687 private:
688 GV gv; // need this object here because FEM and GFS store a const reference !!
689 FEMB femb;
690 CONB conb;
691 std::shared_ptr<GFS> gfsp;
692 std::shared_ptr<CC> ccp;
693 };
694
695 // template specialization for nonoverlapping case
696 template<typename T, typename N, unsigned int degree, typename BCType,
698 typename VBET>
699 class CGSpace<T, N, degree, BCType, gt, mt, SolverCategory::nonoverlapping, VBET> {
700 public:
701
702 // export types
703 typedef T Grid;
704 typedef typename T::LeafGridView GV;
705 typedef typename T::ctype ctype;
707 static const int dim = T::dimension;
708 static const int dimworld = T::dimensionworld;
709
710 typedef CGFEMBase<ES,ctype,N,degree,dim,gt> FEMB;
711 typedef CGCONBase<Grid,degree,gt,mt,SolverCategory::nonoverlapping,BCType> CONB;
712
713 typedef typename FEMB::FEM FEM;
714 typedef typename CONB::CON CON;
715
716 typedef VBET VBE;
717 typedef GridFunctionSpace<ES,FEM,CON,VBE> GFS;
718
719 typedef N NT;
720 using DOF = Backend::Vector<GFS,N>;
722 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
723 typedef VTKGridFunctionAdapter<DGF> VTKF;
724
725 // constructor making the grid function space an all that is needed
726 CGSpace (Grid& grid, const BCType& bctype)
727 : gv(grid.leafGridView()), es(gv), femb(es), conb(grid,bctype)
728 {
729 gfsp = std::shared_ptr<GFS>(new GFS(es,femb.getFEM(),conb.getCON()));
730 gfsp->name("cgspace");
731 // initialize ordering
732 gfsp->update();
733 // conb.postGFSHook(*gfsp);
734 ccp = std::shared_ptr<CC>(new CC());
735 }
736
737 FEM& getFEM()
738 {
739 return femb.getFEM();
740 }
741
742 const FEM& getFEM() const
743 {
744 return femb.getFEM();
745 }
746
747 // return gfs reference
748 GFS& getGFS ()
749 {
750 return *gfsp;
751 }
752
753 // return gfs reference const version
754 const GFS& getGFS () const
755 {
756 return *gfsp;
757 }
758
759 // return gfs reference
760 CC& getCC ()
761 {
762 return *ccp;
763 }
764
765 // return gfs reference const version
766 const CC& getCC () const
767 {
768 return *ccp;
769 }
770
771 void assembleConstraints (const BCType& bctype)
772 {
773 ccp->clear();
774 constraints(bctype,*gfsp,*ccp);
775 }
776
777 void clearConstraints ()
778 {
779 ccp->clear();
780 }
781
782 void setConstrainedDOFS (DOF& x, NT nt) const
783 {
784 set_constrained_dofs(*ccp,nt,x);
785 conb.make_consistent(*gfsp,x);
786 }
787
788 void setNonConstrainedDOFS (DOF& x, NT nt) const
789 {
790 set_nonconstrained_dofs(*ccp,nt,x);
791 conb.make_consistent(*gfsp,x);
792 }
793
794 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
795 {
796 copy_constrained_dofs(*ccp,xin,xout);
797 conb.make_consistent(*gfsp,xout);
798 }
799
800 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
801 {
802 copy_nonconstrained_dofs(*ccp,xin,xout);
803 conb.make_consistent(*gfsp,xout);
804 }
805
806 private:
807 GV gv; // need this object here because FEM and GFS store a const reference !!
808 ES es;
809 FEMB femb;
810 CONB conb;
811 std::shared_ptr<GFS> gfsp;
812 std::shared_ptr<CC> ccp;
813 };
814
815
816
817 //============================================================================
818 // Discontinuous Finite Element Space
819 //============================================================================
820
821 // constraints base template
822 template<SolverCategory::Category st>
823 class DGCONBase
824 {};
825
826 template<>
827 class DGCONBase<SolverCategory::sequential>
828 {
829 public:
830 typedef NoConstraints CON;
831 DGCONBase ()
832 {
833 conp = std::shared_ptr<CON>(new CON());
834 }
835 CON& getCON() {return *conp;}
836 const CON& getCON() const {return *conp;}
837 template<typename GFS, typename DOF>
838 void make_consistent (const GFS& gfs, DOF& x) const {}
839 private:
840 std::shared_ptr<CON> conp;
841 };
842
843 template<>
844 class DGCONBase<SolverCategory::nonoverlapping>
845 {
846 public:
847 typedef P0ParallelGhostConstraints CON;
848 DGCONBase ()
849 {
850 conp = std::shared_ptr<CON>(new CON());
851 }
852 CON& getCON() {return *conp;}
853 const CON& getCON() const {return *conp;}
854 template<typename GFS, typename DOF>
855 void make_consistent (const GFS& gfs, DOF& x) const {}
856 private:
857 std::shared_ptr<CON> conp;
858 };
859
860 template<>
861 class DGCONBase<SolverCategory::overlapping>
862 {
863 public:
864 typedef P0ParallelConstraints CON;
865 DGCONBase ()
866 {
867 conp = std::shared_ptr<CON>(new CON());
868 }
869 CON& getCON() {return *conp;}
870 const CON& getCON() const {return *conp;}
871 template<typename GFS, typename DOF>
872 void make_consistent (const GFS& gfs, DOF& x) const
873 {
874 // make vector consistent; this is needed for all overlapping solvers
875 ISTL::ParallelHelper<GFS> helper(gfs);
876 helper.maskForeignDOFs(Backend::native(x));
877 Dune::PDELab::AddDataHandle<GFS,DOF> adddh(gfs,x);
878 if (gfs.gridView().comm().size()>1)
879 gfs.gridView().communicate(adddh,Dune::InteriorBorder_All_Interface,Dune::ForwardCommunication);
880 }
881 private:
882 std::shared_ptr<CON> conp;
883 };
884
885 // Discontinuous space
886 // default implementation, use only specializations below
887 template<typename T, typename N, unsigned int degree,
889 typename VBET=ISTL::VectorBackend<ISTL::Blocking::fixed,Dune::PB::PkSize<degree,T::dimension>::value> >
890 class DGPkSpace
891 {
892 public:
893
894 // export types
895 typedef T Grid;
896 typedef typename T::LeafGridView GV;
897 typedef typename T::ctype ctype;
898 static const int dim = T::dimension;
899 static const int dimworld = T::dimensionworld;
900 typedef N NT;
901#if HAVE_GMP
902 typedef OPBLocalFiniteElementMap<ctype,NT,degree,dim,gt,Dune::GMPField<512>,Dune::PB::BasisType::Pk> FEM;
903#else
904 typedef OPBLocalFiniteElementMap<ctype,NT,degree,dim,gt> FEM;
905#endif
906 typedef DGCONBase<st> CONB;
907 typedef typename CONB::CON CON;
908 typedef VBET VBE;
909 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
910 using DOF = Backend::Vector<GFS,N>;
912 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
913 typedef VTKGridFunctionAdapter<DGF> VTKF;
914
915 // constructor making the grid function space an all that is needed
916 DGPkSpace (const GV& gridview) : gv(gridview), conb()
917 {
918 femp = std::shared_ptr<FEM>(new FEM());
919 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
920 // initialize ordering
921 gfsp->update();
922 ccp = std::shared_ptr<CC>(new CC());
923 }
924
925 FEM& getFEM() { return *femp; }
926 const FEM& getFEM() const { return *femp; }
927
928 // return gfs reference
929 GFS& getGFS () { return *gfsp; }
930
931 // return gfs reference const version
932 const GFS& getGFS () const {return *gfsp;}
933
934 // return gfs reference
935 CC& getCC () { return *ccp;}
936
937 // return gfs reference const version
938 const CC& getCC () const { return *ccp;}
939
940 template<class BCTYPE>
941 void assembleConstraints (const BCTYPE& bctype)
942 {
943 ccp->clear();
944 constraints(bctype,*gfsp,*ccp);
945 }
946
947 void clearConstraints ()
948 {
949 ccp->clear();
950 }
951
952 void setConstrainedDOFS (DOF& x, NT nt) const
953 {
954 set_constrained_dofs(*ccp,nt,x);
955 conb.make_consistent(*gfsp,x);
956 }
957
958 void setNonConstrainedDOFS (DOF& x, NT nt) const
959 {
960 set_nonconstrained_dofs(*ccp,nt,x);
961 conb.make_consistent(*gfsp,x);
962 }
963
964 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
965 {
966 copy_constrained_dofs(*ccp,xin,xout);
967 conb.make_consistent(*gfsp,xout);
968 }
969
970 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
971 {
972 copy_nonconstrained_dofs(*ccp,xin,xout);
973 conb.make_consistent(*gfsp,xout);
974 }
975
976 private:
977 GV gv; // need this object here because FEM and GFS store a const reference !!
978 CONB conb;
979 std::shared_ptr<FEM> femp;
980 std::shared_ptr<GFS> gfsp;
981 std::shared_ptr<CC> ccp;
982 };
983
984 // Discontinuous space
985 // default implementation, use only specializations below
986 template<typename T, typename N, unsigned int degree,
988 //typename VBET=ISTL::VectorBackend<ISTL::Blocking::fixed,Dune::PB::PkSize<degree,T::dimension>::value> >
989 typename VBET=ISTL::VectorBackend<> >
990 class DGQkOPBSpace
991 {
992 public:
993
994 // export types
995 typedef T Grid;
996 typedef typename T::LeafGridView GV;
997 typedef typename T::ctype ctype;
998 static const int dim = T::dimension;
999 static const int dimworld = T::dimensionworld;
1000 typedef N NT;
1001#if HAVE_GMP
1002 typedef OPBLocalFiniteElementMap<ctype,NT,degree,dim,gt,Dune::GMPField<512>,Dune::PB::BasisType::Qk> FEM;
1003#else
1004 typedef OPBLocalFiniteElementMap<ctype,NT,degree,dim,gt,N,Dune::PB::BasisType::Qk> FEM;
1005#endif
1006 typedef DGCONBase<st> CONB;
1007 typedef typename CONB::CON CON;
1008 typedef VBET VBE;
1009 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
1010 using DOF = Backend::Vector<GFS,N>;
1012 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
1013 typedef VTKGridFunctionAdapter<DGF> VTKF;
1014
1015 // constructor making the grid function space an all that is needed
1016 DGQkOPBSpace (const GV& gridview) : gv(gridview), conb()
1017 {
1018 femp = std::shared_ptr<FEM>(new FEM());
1019 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
1020 // initialize ordering
1021 gfsp->update();
1022 ccp = std::shared_ptr<CC>(new CC());
1023 }
1024
1025 FEM& getFEM() { return *femp; }
1026 const FEM& getFEM() const { return *femp; }
1027
1028 // return gfs reference
1029 GFS& getGFS () { return *gfsp; }
1030
1031 // return gfs reference const version
1032 const GFS& getGFS () const {return *gfsp;}
1033
1034 // return gfs reference
1035 CC& getCC () { return *ccp;}
1036
1037 // return gfs reference const version
1038 const CC& getCC () const { return *ccp;}
1039
1040 template<class BCTYPE>
1041 void assembleConstraints (const BCTYPE& bctype)
1042 {
1043 ccp->clear();
1044 constraints(bctype,*gfsp,*ccp);
1045 }
1046
1047 void clearConstraints ()
1048 {
1049 ccp->clear();
1050 }
1051
1052 void setConstrainedDOFS (DOF& x, NT nt) const
1053 {
1054 set_constrained_dofs(*ccp,nt,x);
1055 conb.make_consistent(*gfsp,x);
1056 }
1057
1058 void setNonConstrainedDOFS (DOF& x, NT nt) const
1059 {
1060 set_nonconstrained_dofs(*ccp,nt,x);
1061 conb.make_consistent(*gfsp,x);
1062 }
1063
1064 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
1065 {
1066 copy_constrained_dofs(*ccp,xin,xout);
1067 conb.make_consistent(*gfsp,xout);
1068 }
1069
1070 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
1071 {
1072 copy_nonconstrained_dofs(*ccp,xin,xout);
1073 conb.make_consistent(*gfsp,xout);
1074 }
1075
1076 private:
1077 GV gv; // need this object here because FEM and GFS store a const reference !!
1078 CONB conb;
1079 std::shared_ptr<FEM> femp;
1080 std::shared_ptr<GFS> gfsp;
1081 std::shared_ptr<CC> ccp;
1082 };
1083
1084 // Discontinuous space
1085 // default implementation, use only specializations below
1086 template<typename T, typename N, unsigned int degree,
1088 typename VBET=ISTL::VectorBackend<ISTL::Blocking::fixed,Dune::QkStuff::QkSize<degree,T::dimension>::value> >
1089 class DGQkSpace
1090 {
1091 public:
1092
1093 // export types
1094 typedef T Grid;
1095 typedef typename T::LeafGridView GV;
1096 typedef typename T::ctype ctype;
1097 static const int dim = T::dimension;
1098 static const int dimworld = T::dimensionworld;
1099 typedef N NT;
1100 typedef QkDGLocalFiniteElementMap<ctype,NT,degree,dim> FEM;
1101 typedef DGCONBase<st> CONB;
1102 typedef typename CONB::CON CON;
1103 typedef VBET VBE;
1104 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
1105 using DOF = Backend::Vector<GFS,N>;
1107 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
1108 typedef VTKGridFunctionAdapter<DGF> VTKF;
1109
1110 // constructor making the grid function space an all that is needed
1111 DGQkSpace (const GV& gridview) : gv(gridview), conb()
1112 {
1113 femp = std::shared_ptr<FEM>(new FEM());
1114 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
1115 // initialize ordering
1116 gfsp->update();
1117 ccp = std::shared_ptr<CC>(new CC());
1118 }
1119
1120 FEM& getFEM() { return *femp; }
1121 const FEM& getFEM() const { return *femp; }
1122
1123 // return gfs reference
1124 GFS& getGFS () { return *gfsp; }
1125
1126 // return gfs reference const version
1127 const GFS& getGFS () const {return *gfsp;}
1128
1129 // return gfs reference
1130 CC& getCC () { return *ccp;}
1131
1132 // return gfs reference const version
1133 const CC& getCC () const { return *ccp;}
1134
1135 template<class BCTYPE>
1136 void assembleConstraints (const BCTYPE& bctype)
1137 {
1138 ccp->clear();
1139 constraints(bctype,*gfsp,*ccp);
1140 }
1141
1142 void clearConstraints ()
1143 {
1144 ccp->clear();
1145 }
1146
1147 void setConstrainedDOFS (DOF& x, NT nt) const
1148 {
1149 set_constrained_dofs(*ccp,nt,x);
1150 conb.make_consistent(*gfsp,x);
1151 }
1152
1153 void setNonConstrainedDOFS (DOF& x, NT nt) const
1154 {
1155 set_nonconstrained_dofs(*ccp,nt,x);
1156 conb.make_consistent(*gfsp,x);
1157 }
1158
1159 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
1160 {
1161 copy_constrained_dofs(*ccp,xin,xout);
1162 conb.make_consistent(*gfsp,xout);
1163 }
1164
1165 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
1166 {
1167 copy_nonconstrained_dofs(*ccp,xin,xout);
1168 conb.make_consistent(*gfsp,xout);
1169 }
1170
1171 private:
1172 GV gv; // need this object here because FEM and GFS store a const reference !!
1173 CONB conb;
1174 std::shared_ptr<FEM> femp;
1175 std::shared_ptr<GFS> gfsp;
1176 std::shared_ptr<CC> ccp;
1177 };
1178
1179
1180 // Discontinuous space using QK with Gauss Lobatto points (use only for cube elements)
1181 template<typename T, typename N, unsigned int degree,
1183 //typename VBET=ISTL::VectorBackend<ISTL::Blocking::fixed,Dune::QkStuff::QkSize<degree,T::dimension>::value> >
1184 typename VBET=ISTL::VectorBackend<> >
1185 class DGQkGLSpace
1186 {
1187 public:
1188
1189 // export types
1190 typedef T Grid;
1191 typedef typename T::LeafGridView GV;
1192 typedef typename T::ctype ctype;
1193 static const int dim = T::dimension;
1194 static const int dimworld = T::dimensionworld;
1195 typedef N NT;
1196 typedef QkDGLocalFiniteElementMap<ctype,NT,degree,dim,QkDGBasisPolynomial::lobatto> FEM;
1197 typedef DGCONBase<st> CONB;
1198 typedef typename CONB::CON CON;
1199 typedef VBET VBE;
1200 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
1201 using DOF = Backend::Vector<GFS,N>;
1203 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
1204 typedef VTKGridFunctionAdapter<DGF> VTKF;
1205
1206 // constructor making the grid function space an all that is needed
1207 DGQkGLSpace (const GV& gridview) : gv(gridview), conb()
1208 {
1209 femp = std::shared_ptr<FEM>(new FEM());
1210 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
1211 // initialize ordering
1212 gfsp->update();
1213 ccp = std::shared_ptr<CC>(new CC());
1214 }
1215
1216 FEM& getFEM() { return *femp; }
1217 const FEM& getFEM() const { return *femp; }
1218
1219 // return gfs reference
1220 GFS& getGFS () { return *gfsp; }
1221
1222 // return gfs reference const version
1223 const GFS& getGFS () const {return *gfsp;}
1224
1225 // return gfs reference
1226 CC& getCC () { return *ccp;}
1227
1228 // return gfs reference const version
1229 const CC& getCC () const { return *ccp;}
1230
1231 template<class BCTYPE>
1232 void assembleConstraints (const BCTYPE& bctype)
1233 {
1234 ccp->clear();
1235 constraints(bctype,*gfsp,*ccp);
1236 }
1237
1238 void clearConstraints ()
1239 {
1240 ccp->clear();
1241 }
1242
1243 void setConstrainedDOFS (DOF& x, NT nt) const
1244 {
1245 set_constrained_dofs(*ccp,nt,x);
1246 conb.make_consistent(*gfsp,x);
1247 }
1248
1249 void setNonConstrainedDOFS (DOF& x, NT nt) const
1250 {
1251 set_nonconstrained_dofs(*ccp,nt,x);
1252 conb.make_consistent(*gfsp,x);
1253 }
1254
1255 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
1256 {
1257 copy_constrained_dofs(*ccp,xin,xout);
1258 conb.make_consistent(*gfsp,xout);
1259 }
1260
1261 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
1262 {
1263 copy_nonconstrained_dofs(*ccp,xin,xout);
1264 conb.make_consistent(*gfsp,xout);
1265 }
1266
1267 private:
1268 GV gv; // need this object here because FEM and GFS store a const reference !!
1269 CONB conb;
1270 std::shared_ptr<FEM> femp;
1271 std::shared_ptr<GFS> gfsp;
1272 std::shared_ptr<CC> ccp;
1273 };
1274
1275
1276 // Discontinuous space using Legendre polynomials (use only for cube elements)
1277 template<typename T, typename N, unsigned int degree,
1279 //typename VBET=ISTL::VectorBackend<ISTL::Blocking::fixed,Dune::QkStuff::QkSize<degree,T::dimension>::value> >
1280 typename VBET=ISTL::VectorBackend<> >
1281 class DGLegendreSpace
1282 {
1283 public:
1284
1285 // export types
1286 typedef T Grid;
1287 typedef typename T::LeafGridView GV;
1288 typedef typename T::ctype ctype;
1289 static const int dim = T::dimension;
1290 static const int dimworld = T::dimensionworld;
1291 typedef N NT;
1292 typedef QkDGLocalFiniteElementMap<ctype,NT,degree,dim,QkDGBasisPolynomial::legendre> FEM;
1293 typedef DGCONBase<st> CONB;
1294 typedef typename CONB::CON CON;
1295 typedef VBET VBE;
1296 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
1297 using DOF = Backend::Vector<GFS,N>;
1299 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
1300 typedef VTKGridFunctionAdapter<DGF> VTKF;
1301
1302 // constructor making the grid function space an all that is needed
1303 DGLegendreSpace (const GV& gridview) : gv(gridview), conb()
1304 {
1305 femp = std::shared_ptr<FEM>(new FEM());
1306 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
1307 // initialize ordering
1308 gfsp->update();
1309 ccp = std::shared_ptr<CC>(new CC());
1310 }
1311
1312 FEM& getFEM() { return *femp; }
1313 const FEM& getFEM() const { return *femp; }
1314
1315 // return gfs reference
1316 GFS& getGFS () { return *gfsp; }
1317
1318 // return gfs reference const version
1319 const GFS& getGFS () const {return *gfsp;}
1320
1321 // return gfs reference
1322 CC& getCC () { return *ccp;}
1323
1324 // return gfs reference const version
1325 const CC& getCC () const { return *ccp;}
1326
1327 template<class BCTYPE>
1328 void assembleConstraints (const BCTYPE& bctype)
1329 {
1330 ccp->clear();
1331 constraints(bctype,*gfsp,*ccp);
1332 }
1333
1334 void clearConstraints ()
1335 {
1336 ccp->clear();
1337 }
1338
1339 void setConstrainedDOFS (DOF& x, NT nt) const
1340 {
1341 set_constrained_dofs(*ccp,nt,x);
1342 conb.make_consistent(*gfsp,x);
1343 }
1344
1345 void setNonConstrainedDOFS (DOF& x, NT nt) const
1346 {
1347 set_nonconstrained_dofs(*ccp,nt,x);
1348 conb.make_consistent(*gfsp,x);
1349 }
1350
1351 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
1352 {
1353 copy_constrained_dofs(*ccp,xin,xout);
1354 conb.make_consistent(*gfsp,xout);
1355 }
1356
1357 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
1358 {
1359 copy_nonconstrained_dofs(*ccp,xin,xout);
1360 conb.make_consistent(*gfsp,xout);
1361 }
1362
1363 private:
1364 GV gv; // need this object here because FEM and GFS store a const reference !!
1365 CONB conb;
1366 std::shared_ptr<FEM> femp;
1367 std::shared_ptr<GFS> gfsp;
1368 std::shared_ptr<CC> ccp;
1369 };
1370
1371
1372 // Discontinuous P0 space
1373 template<typename T, typename N,
1375 typename VBET=ISTL::VectorBackend<> >
1376 class P0Space
1377 {
1378 public:
1379
1380 // export types
1381 typedef T Grid;
1382 typedef typename T::LeafGridView GV;
1383 typedef typename T::ctype ctype;
1384 static const int dim = T::dimension;
1385 static const int dimworld = T::dimensionworld;
1386 typedef N NT;
1388 typedef DGCONBase<st> CONB;
1389 typedef typename CONB::CON CON;
1390 typedef VBET VBE;
1391 typedef GridFunctionSpace<GV,FEM,CON,VBE> GFS;
1392 using DOF = Backend::Vector<GFS,N>;
1394 typedef typename GFS::template ConstraintsContainer<N>::Type CC;
1395 typedef VTKGridFunctionAdapter<DGF> VTKF;
1396
1397 // constructor making the grid function space an all that is needed
1398 P0Space (const GV& gridview) : gv(gridview), conb()
1399 {
1400 auto geometryType = geometryTypeFromBasicType(gt, dim);
1401 femp = std::shared_ptr<FEM>(new FEM(geometryType));
1402 gfsp = std::shared_ptr<GFS>(new GFS(gv,*femp));
1403 // initialize ordering
1404 gfsp->update();
1405 ccp = std::shared_ptr<CC>(new CC());
1406 }
1407
1408 FEM& getFEM() { return *femp; }
1409 const FEM& getFEM() const { return *femp; }
1410
1411 // return gfs reference
1412 GFS& getGFS () { return *gfsp; }
1413
1414 // return gfs reference const version
1415 const GFS& getGFS () const {return *gfsp;}
1416
1417 // return gfs reference
1418 CC& getCC () { return *ccp;}
1419
1420 // return gfs reference const version
1421 const CC& getCC () const { return *ccp;}
1422
1423 template<class BCTYPE>
1424 void assembleConstraints (const BCTYPE& bctype)
1425 {
1426 ccp->clear();
1427 constraints(bctype,*gfsp,*ccp);
1428 }
1429
1430 void clearConstraints ()
1431 {
1432 ccp->clear();
1433 }
1434
1435 void setConstrainedDOFS (DOF& x, NT nt) const
1436 {
1437 set_constrained_dofs(*ccp,nt,x);
1438 conb.make_consistent(*gfsp,x);
1439 }
1440
1441 void setNonConstrainedDOFS (DOF& x, NT nt) const
1442 {
1443 set_nonconstrained_dofs(*ccp,nt,x);
1444 conb.make_consistent(*gfsp,x);
1445 }
1446
1447 void copyConstrainedDOFS (const DOF& xin, DOF& xout) const
1448 {
1449 copy_constrained_dofs(*ccp,xin,xout);
1450 conb.make_consistent(*gfsp,xout);
1451 }
1452
1453 void copyNonConstrainedDOFS (const DOF& xin, DOF& xout) const
1454 {
1455 copy_nonconstrained_dofs(*ccp,xin,xout);
1456 conb.make_consistent(*gfsp,xout);
1457 }
1458
1459 private:
1460 GV gv; // need this object here because FEM and GFS store a const reference !!
1461 CONB conb;
1462 std::shared_ptr<FEM> femp;
1463 std::shared_ptr<GFS> gfsp;
1464 std::shared_ptr<CC> ccp;
1465 };
1466
1467
1468 // how can we most easily specify a grid function
1469 // pass a function space as parameter
1470 template<typename FS, typename Functor>
1471 class UserFunction
1472 : public GridFunctionBase<GridFunctionTraits<typename FS::GV, typename FS::NT,
1473 1,FieldVector<typename FS::NT,1> >
1474 ,UserFunction<FS,Functor> >
1475 {
1476 public:
1477 typedef GridFunctionTraits<typename FS::GV, typename FS::NT,
1478 1,FieldVector<typename FS::NT,1> > Traits;
1479
1481 UserFunction (const FS& fs_, const Functor& f_)
1482 : fs(fs_), f(f_)
1483 {}
1484
1486 inline void evaluate (const typename Traits::ElementType& e,
1487 const typename Traits::DomainType& x,
1488 typename Traits::RangeType& y) const
1489 {
1490 typename Traits::DomainType x_ = e.geometry().global(x);
1491 std::vector<double> x__(x.size());
1492 for (size_t i=0; i<x.size(); ++i) x__[i]=x_[i];
1493 y = f(x__);
1494 }
1495
1496 inline const typename FS::GV& getGridView () const
1497 {
1498 return fs.getGFS().gridView();
1499 }
1500
1501 private:
1502 const FS fs; // store a copy of the function space
1503 const Functor f;
1504 };
1505
1506
1507 template<typename FS, typename LOP, SolverCategory::Category st = SolverCategory::sequential>
1508 class GalerkinGlobalAssembler
1509 {
1510 public:
1511 // export types
1512 typedef ISTL::BCRSMatrixBackend<> MBE;
1513 typedef Dune::PDELab::GridOperator<typename FS::GFS,typename FS::GFS,LOP,MBE,
1514 typename FS::NT,typename FS::NT,typename FS::NT,
1515 typename FS::CC,typename FS::CC> GO;
1516 typedef typename GO::Jacobian MAT;
1517
1518 GalerkinGlobalAssembler (const FS& fs, LOP& lop, const std::size_t nonzeros)
1519 {
1520 gop = std::shared_ptr<GO>(new GO(fs.getGFS(),fs.getCC(),fs.getGFS(),fs.getCC(),lop,MBE(nonzeros)));
1521 }
1522
1523 // return grid reference
1524 GO& getGO ()
1525 {
1526 return *gop;
1527 }
1528
1529 // return grid reference const version
1530 const GO& getGO () const
1531 {
1532 return *gop;
1533 }
1534
1535 GO& operator*()
1536 {
1537 return *gop;
1538 }
1539
1540 GO* operator->()
1541 {
1542 return gop.operator->();
1543 }
1544
1545 const GO& operator*() const
1546 {
1547 return *gop;
1548 }
1549
1550 const GO* operator->() const
1551 {
1552 return gop.operator->();
1553 }
1554
1555 private:
1556 std::shared_ptr<GO> gop;
1557 };
1558
1559
1560 template<typename FS, typename LOP, SolverCategory::Category st = SolverCategory::sequential>
1561 class GalerkinGlobalAssemblerNewBackend
1562 {
1563 public:
1564 // export types
1566 typedef Dune::PDELab::GridOperator<typename FS::GFS,typename FS::GFS,LOP,MBE,
1567 typename FS::NT,typename FS::NT,typename FS::NT,
1568 typename FS::CC,typename FS::CC> GO;
1569 typedef typename GO::Jacobian MAT;
1570
1571 GalerkinGlobalAssemblerNewBackend (const FS& fs, LOP& lop, const MBE& mbe)
1572 {
1573 gop = std::shared_ptr<GO>(new GO(fs.getGFS(),fs.getCC(),fs.getGFS(),fs.getCC(),lop,mbe));
1574 }
1575
1576 // return grid reference
1577 GO& getGO ()
1578 {
1579 return *gop;
1580 }
1581
1582 // return grid reference const version
1583 const GO& getGO () const
1584 {
1585 return *gop;
1586 }
1587
1588 GO& operator*()
1589 {
1590 return *gop;
1591 }
1592
1593 GO* operator->()
1594 {
1595 return gop.operator->();
1596 }
1597
1598 const GO& operator*() const
1599 {
1600 return *gop;
1601 }
1602
1603 const GO* operator->() const
1604 {
1605 return gop.operator->();
1606 }
1607
1608 private:
1609 std::shared_ptr<GO> gop;
1610 };
1611
1612
1613 // variant with two different function spaces
1614 template<typename FSU, typename FSV, typename LOP, SolverCategory::Category st>
1615 class GlobalAssembler
1616 {
1617 public:
1618 // export types
1619 typedef ISTL::BCRSMatrixBackend<> MBE;
1620 typedef Dune::PDELab::GridOperator<typename FSU::GFS,typename FSV::GFS,LOP,MBE,
1621 typename FSU::NT,typename FSU::NT,typename FSU::NT,
1622 typename FSU::CC,typename FSV::CC> GO;
1623 typedef typename GO::Jacobian MAT;
1624
1625 GlobalAssembler (const FSU& fsu, const FSV& fsv, LOP& lop, const std::size_t nonzeros)
1626 {
1627 gop = std::shared_ptr<GO>(new GO(fsu.getGFS(),fsu.getCC(),fsv.getGFS(),fsv.getCC(),lop,MBE(nonzeros)));
1628 }
1629
1630 // return grid reference
1631 GO& getGO ()
1632 {
1633 return *gop;
1634 }
1635
1636 // return grid reference const version
1637 const GO& getGO () const
1638 {
1639 return *gop;
1640 }
1641
1642 GO& operator*()
1643 {
1644 return *gop;
1645 }
1646
1647 GO* operator->()
1648 {
1649 return gop.operator->();
1650 }
1651
1652 const GO& operator*() const
1653 {
1654 return *gop;
1655 }
1656
1657 const GO* operator->() const
1658 {
1659 return gop.operator->();
1660 }
1661
1662 private:
1663 std::shared_ptr<GO> gop;
1664 };
1665
1666
1667 template<typename GO1, typename GO2, bool implicit = true>
1668 class OneStepGlobalAssembler
1669 {
1670 public:
1671 // export types
1672 typedef ISTL::BCRSMatrixBackend<> MBE;
1673 typedef Dune::PDELab::OneStepGridOperator<typename GO1::GO,typename GO2::GO,implicit> GO;
1674 typedef typename GO::Jacobian MAT;
1675
1676 OneStepGlobalAssembler (GO1& go1, GO2& go2)
1677 {
1678 gop = std::shared_ptr<GO>(new GO(*go1,*go2));
1679 }
1680
1681 // return grid reference
1682 GO& getGO ()
1683 {
1684 return *gop;
1685 }
1686
1687 // return grid reference const version
1688 const GO& getGO () const
1689 {
1690 return *gop;
1691 }
1692
1693 GO& operator*()
1694 {
1695 return *gop;
1696 }
1697
1698 GO* operator->()
1699 {
1700 return gop.operator->();
1701 }
1702
1703 const GO& operator*() const
1704 {
1705 return *gop;
1706 }
1707
1708 const GO* operator->() const
1709 {
1710 return gop.operator->();
1711 }
1712
1713 private:
1714 std::shared_ptr<GO> gop;
1715 };
1716
1717
1718 // packaging of the CG_AMG_SSOR solver: default version is sequential
1719 template<typename FS, typename ASS, SolverCategory::Category st = SolverCategory::sequential>
1720 class ISTLSolverBackend_CG_AMG_SSOR
1721 {
1722 public:
1723 // types exported
1725
1726 ISTLSolverBackend_CG_AMG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1727 int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
1728 {
1729 lsp = std::shared_ptr<LS>(new LS(maxiter_,verbose_,reuse_,usesuperlu_));
1730 }
1731
1732 LS& getLS () {return *lsp;}
1733 const LS& getLS () const { return *lsp;}
1734 LS& operator*(){return *lsp;}
1735 LS* operator->() { return lsp.operator->(); }
1736 const LS& operator*() const{return *lsp;}
1737 const LS* operator->() const{ return lsp.operator->();}
1738
1739 private:
1740 std::shared_ptr<LS> lsp;
1741 };
1742
1743 // packaging of the CG_AMG_SSOR solver: nonoverlapping version
1744 template<typename FS, typename ASS>
1745 class ISTLSolverBackend_CG_AMG_SSOR<FS,ASS, SolverCategory::nonoverlapping>
1746 {
1747 public:
1748 // types exported
1750
1751 ISTLSolverBackend_CG_AMG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1752 int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
1753 {
1754 lsp = std::shared_ptr<LS>(new LS(fs.getGFS(),maxiter_,verbose_,reuse_,usesuperlu_));
1755 }
1756
1757 LS& getLS () {return *lsp;}
1758 const LS& getLS () const { return *lsp;}
1759 LS& operator*(){return *lsp;}
1760 LS* operator->() { return lsp.operator->(); }
1761 const LS& operator*() const{return *lsp;}
1762 const LS* operator->() const{ return lsp.operator->();}
1763
1764 private:
1765 std::shared_ptr<LS> lsp;
1766 };
1767
1768 // packaging of the CG_AMG_SSOR solver: overlapping version
1769 template<typename FS, typename ASS>
1770 class ISTLSolverBackend_CG_AMG_SSOR<FS,ASS, SolverCategory::overlapping>
1771 {
1772 public:
1773 // types exported
1775
1776 ISTLSolverBackend_CG_AMG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1777 int verbose_=1, bool reuse_=false, bool usesuperlu_=true)
1778 {
1779 lsp = std::shared_ptr<LS>(new LS(fs.getGFS(),maxiter_,verbose_,reuse_,usesuperlu_));
1780 }
1781
1782 LS& getLS () {return *lsp;}
1783 const LS& getLS () const { return *lsp;}
1784 LS& operator*(){return *lsp;}
1785 LS* operator->() { return lsp.operator->(); }
1786 const LS& operator*() const{return *lsp;}
1787 const LS* operator->() const{ return lsp.operator->();}
1788
1789 private:
1790 std::shared_ptr<LS> lsp;
1791 };
1792
1793 // packaging of the CG_SSOR solver: default version is sequential
1794 template<typename FS, typename ASS, SolverCategory::Category st = SolverCategory::sequential>
1795 class ISTLSolverBackend_CG_SSOR
1796 {
1797 public:
1798 // types exported
1799 typedef ISTLBackend_SEQ_CG_SSOR LS;
1800
1801 ISTLSolverBackend_CG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1802 int steps_=5, int verbose_=1)
1803 {
1804 lsp = std::shared_ptr<LS>(new LS(maxiter_,verbose_));
1805 }
1806
1807 LS& getLS () {return *lsp;}
1808 const LS& getLS () const { return *lsp;}
1809 LS& operator*(){return *lsp;}
1810 LS* operator->() { return lsp.operator->(); }
1811 const LS& operator*() const{return *lsp;}
1812 const LS* operator->() const{ return lsp.operator->();}
1813
1814 private:
1815 std::shared_ptr<LS> lsp;
1816 };
1817
1818 // packaging of the CG_SSOR solver: nonoverlapping version
1819 template<typename FS, typename ASS>
1820 class ISTLSolverBackend_CG_SSOR<FS,ASS,SolverCategory::nonoverlapping>
1821 {
1822 public:
1823 // types exported
1824 typedef ISTLBackend_NOVLP_CG_SSORk<typename ASS::GO> LS;
1825
1826 ISTLSolverBackend_CG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1827 int steps_=5, int verbose_=1)
1828 {
1829 lsp = std::shared_ptr<LS>(new LS(fs.getGFS(),maxiter_,steps_,verbose_));
1830 }
1831
1832 LS& getLS () {return *lsp;}
1833 const LS& getLS () const { return *lsp;}
1834 LS& operator*(){return *lsp;}
1835 LS* operator->() { return lsp.operator->(); }
1836 const LS& operator*() const{return *lsp;}
1837 const LS* operator->() const{ return lsp.operator->();}
1838
1839 private:
1840 std::shared_ptr<LS> lsp;
1841 };
1842
1843 // packaging of the CG_SSOR solver: overlapping version
1844 template<typename FS, typename ASS>
1845 class ISTLSolverBackend_CG_SSOR<FS,ASS,SolverCategory::overlapping>
1846 {
1847 public:
1848 // types exported
1849 typedef ISTLBackend_OVLP_CG_SSORk<typename FS::GFS, typename FS::CC> LS;
1850
1851 ISTLSolverBackend_CG_SSOR (const FS& fs, const ASS& ass, unsigned maxiter_=5000,
1852 int steps_=5, int verbose_=1)
1853 {
1854 lsp = std::shared_ptr<LS>(new LS(fs.getGFS(),fs.getCC(),maxiter_,steps_,verbose_));
1855 }
1856
1857 LS& getLS () {return *lsp;}
1858 const LS& getLS () const { return *lsp;}
1859 LS& operator*(){return *lsp;}
1860 LS* operator->() { return lsp.operator->(); }
1861 const LS& operator*() const{return *lsp;}
1862 const LS* operator->() const{ return lsp.operator->();}
1863
1864 private:
1865 std::shared_ptr<LS> lsp;
1866 };
1867
1868
1869 // packaging of a default solver that should always work
1870 // in the sequential case : BCGS SSOR
1871 template<typename FS, typename ASS, SolverCategory::Category st = SolverCategory::sequential>
1872 class ISTLSolverBackend_IterativeDefault
1873 {
1874 public:
1875 // types exported
1876 typedef ISTLBackend_SEQ_BCGS_SSOR LS;
1877
1878 ISTLSolverBackend_IterativeDefault (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
1879 {
1880 lsp = std::shared_ptr<LS>(new LS(maxiter_,verbose_));
1881 }
1882
1883 LS& getLS () {return *lsp;}
1884 const LS& getLS () const { return *lsp;}
1885 LS& operator*(){return *lsp;}
1886 LS* operator->() { return lsp.operator->(); }
1887 const LS& operator*() const{return *lsp;}
1888 const LS* operator->() const{ return lsp.operator->();}
1889
1890 private:
1891 std::shared_ptr<LS> lsp;
1892 };
1893
1894 // in the nonoverlapping case : BCGS SSORk
1895 template<typename FS, typename ASS>
1896 class ISTLSolverBackend_IterativeDefault<FS,ASS,SolverCategory::nonoverlapping>
1897 {
1898 public:
1899 // types exported
1901
1902 ISTLSolverBackend_IterativeDefault (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
1903 {
1904 lsp = std::shared_ptr<LS>(new LS(ass.getGO(),maxiter_,3,verbose_));
1905 }
1906
1907 LS& getLS () {return *lsp;}
1908 const LS& getLS () const { return *lsp;}
1909 LS& operator*(){return *lsp;}
1910 LS* operator->() { return lsp.operator->(); }
1911 const LS& operator*() const{return *lsp;}
1912 const LS* operator->() const{ return lsp.operator->();}
1913
1914 private:
1915 std::shared_ptr<LS> lsp;
1916 };
1917
1918 // in the overlapping case : BCGS SSORk
1919 template<typename FS, typename ASS>
1920 class ISTLSolverBackend_IterativeDefault<FS,ASS,SolverCategory::overlapping>
1921 {
1922 public:
1923 // types exported
1924 typedef ISTLBackend_OVLP_BCGS_SSORk<typename FS::GFS, typename FS::CC> LS;
1925
1926 ISTLSolverBackend_IterativeDefault (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
1927 {
1928 lsp = std::shared_ptr<LS>(new LS(fs.getGFS(),fs.getCC(),maxiter_,3,verbose_));
1929 }
1930
1931 LS& getLS () {return *lsp;}
1932 const LS& getLS () const { return *lsp;}
1933 LS& operator*(){return *lsp;}
1934 LS* operator->() { return lsp.operator->(); }
1935 const LS& operator*() const{return *lsp;}
1936 const LS* operator->() const{ return lsp.operator->();}
1937
1938 private:
1939 std::shared_ptr<LS> lsp;
1940 };
1941
1942 // packaging of a default solver that should always work
1943 // in the sequential case : BCGS SSOR
1944 template<typename FS, typename ASS, SolverCategory::Category st = SolverCategory::sequential>
1945 class ISTLSolverBackend_ExplicitDiagonal
1946 {
1947 public:
1948 // types exported
1950
1951 ISTLSolverBackend_ExplicitDiagonal (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
1952 {
1953 lsp = std::shared_ptr<LS>(new LS());
1954 }
1955
1956 LS& getLS () {return *lsp;}
1957 const LS& getLS () const { return *lsp;}
1958 LS& operator*(){return *lsp;}
1959 LS* operator->() { return lsp.operator->(); }
1960 const LS& operator*() const{return *lsp;}
1961 const LS* operator->() const{ return lsp.operator->();}
1962
1963 private:
1964 std::shared_ptr<LS> lsp;
1965 };
1966
1967 // packaging of a default solver that should always work
1968 // in the sequential case : BCGS SSOR
1969 template<typename FS, typename ASS>
1970 class ISTLSolverBackend_ExplicitDiagonal<FS,ASS,SolverCategory::overlapping>
1971 {
1972 public:
1973 // types exported
1975
1976 ISTLSolverBackend_ExplicitDiagonal (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
1977 {
1978 lsp = std::shared_ptr<LS>(new LS(fs.getGFS()));
1979 }
1980
1981 LS& getLS () {return *lsp;}
1982 const LS& getLS () const { return *lsp;}
1983 LS& operator*(){return *lsp;}
1984 LS* operator->() { return lsp.operator->(); }
1985 const LS& operator*() const{return *lsp;}
1986 const LS* operator->() const{ return lsp.operator->();}
1987
1988 private:
1989 std::shared_ptr<LS> lsp;
1990 };
1991
1992 // packaging of a default solver that should always work
1993 // in the sequential case : BCGS SSOR
1994 template<typename FS, typename ASS>
1995 class ISTLSolverBackend_ExplicitDiagonal<FS,ASS,SolverCategory::nonoverlapping>
1996 {
1997 public:
1998 // types exported
2000
2001 ISTLSolverBackend_ExplicitDiagonal (const FS& fs, const ASS& ass, unsigned maxiter_=5000, int verbose_=1)
2002 {
2003 lsp = std::shared_ptr<LS>(new LS(fs.getGFS()));
2004 }
2005
2006 LS& getLS () {return *lsp;}
2007 const LS& getLS () const { return *lsp;}
2008 LS& operator*(){return *lsp;}
2009 LS* operator->() { return lsp.operator->(); }
2010 const LS& operator*() const{return *lsp;}
2011 const LS* operator->() const{ return lsp.operator->();}
2012
2013 private:
2014 std::shared_ptr<LS> lsp;
2015 };
2016
2017
2018} // end namespace PDELab
2019 } // end namespace Dune
2020
2021#endif // DUNE_PDELAB_BOILERPLATE_PDELAB_HH
The AMG preconditioner.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
BasicType
Each entity can be tagged by one of these basic types plus its space dimension.
Definition: type.hh:120
@ cube
Cube element in any nonnegative dimension.
Definition: type.hh:122
@ simplex
Simplicial element in any nonnegative dimension.
Definition: type.hh:121
static std::unique_ptr< Grid > read(const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:935
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:275
virtual std::unique_ptr< GridType > createGrid()
Finalize grid creation and hand over the grid.
Definition: gridfactory.hh:333
static constexpr int dimensionworld
The dimension of the world the grid lives in.
Definition: grid.hh:390
ct ctype
Define type used for coordinates in grid module.
Definition: grid.hh:518
Dirichlet Constraints construction.
Definition: conforming.hh:38
convert a grid function space and a coefficient vector into a grid function
Definition: gridfunctionspaceutilities.hh:76
Standard grid operator implementation.
Definition: gridoperator.hh:36
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperator.hh:47
Overlapping parallel conjugate gradient solver preconditioned with AMG smoothed by SSOR.
Definition: ovlpistlsolverbackend.hh:1259
Nonoverlapping parallel BiCGSTAB solver preconditioned by block SSOR.
Definition: novlpistlsolverbackend.hh:838
Nonoverlapping parallel CG solver preconditioned with AMG smoothed by SSOR.
Definition: novlpistlsolverbackend.hh:1076
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: novlpistlsolverbackend.hh:637
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: ovlpistlsolverbackend.hh:1011
Sequential conjugate gradient solver preconditioned with AMG smoothed by SSOR.
Definition: seqistlsolverbackend.hh:856
Solver to be used for explicit time-steppers with (block-)diagonal mass matrix.
Definition: seqistlsolverbackend.hh:661
Partition view (or entity set) of a grid view.
Definition: partitionviewentityset.hh:120
A free function to provide the demangled class name of a given object or type as a string.
A few common exception classes.
@ conforming
Output conforming data.
Definition: common.hh:73
@ nonconforming
Output non-conforming data.
Definition: common.hh:81
GeometryType geometryType(const Dune::GeometryType &t)
mapping from GeometryType to VTKGeometryType
Definition: common.hh:151
Convenience header which includes all available VTK writers.
Some generic functions for pretty printing vectors and matrices.
Implementations of the inverse operator interface.
Utility to generate an array with a certain value.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
auto periodic(RawPreBasisIndicator &&rawPreBasisIndicator, PIS &&periodicIndexSet)
Create a pre-basis factory that can create a periodic pre-basis.
Definition: periodicbasis.hh:199
Grid< dim, dimworld, ct, GridFamily >::LeafGridView leafGridView(const Grid< dim, dimworld, ct, GridFamily > &grid)
leaf grid view for the given grid
Definition: grid.hh:805
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:171
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:88
void constraints(const GFS &gfs, CG &cg, const bool verbose=false)
construct constraints
Definition: constraints.hh:749
void set_nonconstrained_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
Definition: constraints.hh:960
void set_constrained_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
construct constraints from given boundary condition function
Definition: constraints.hh:796
void copy_constrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:936
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:987
std::size_t degree(const Node &node)
Returns the degree of node as run time information.
Definition: nodeinterface.hh:79
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:13
The OneDGrid class.
Define general, extensible interface for operators. The available implementation wraps a matrix.
Various parser methods to get data into a ParameterTree object.
Define general preconditioner interface.
Dune::FieldVector< GV::Grid::ctype, GV::dimension > DomainType
domain type in dim-size coordinates
Definition: function.hh:50
Backend using (possibly nested) ISTL BCRSMatrices.
Definition: bcrsmatrixbackend.hh:188
GV::Traits::template Codim< 0 >::Entity ElementType
codim 0 entity
Definition: function.hh:119
Category
Definition: solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:25
@ nonoverlapping
Category for non-overlapping solvers.
Definition: solvercategory.hh:27
@ overlapping
Category for overlapping solvers.
Definition: solvercategory.hh:29
A unique label for each type of element that can occur in a grid.
The UGGrid class.
A class to construct structured cube and simplex grids using the grid factory.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)