Dune Core Modules (2.9.0)

nedelecsimplexinterpolation.hh
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 // SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5 #ifndef DUNE_LOCALFUNCTIONS_NEDELEC_NEDELECSIMPLEX_NEDELECSIMPLEXINTERPOLATION_HH
6 #define DUNE_LOCALFUNCTIONS_NEDELEC_NEDELECSIMPLEX_NEDELECSIMPLEXINTERPOLATION_HH
7 
8 #include <fstream>
9 #include <utility>
10 #include <numeric>
11 
13 
15 #include <dune/geometry/referenceelements.hh>
16 #include <dune/geometry/type.hh>
17 
18 #include <dune/localfunctions/common/localkey.hh>
19 #include <dune/localfunctions/utility/interpolationhelper.hh>
20 #include <dune/localfunctions/utility/polynomialbasis.hh>
21 #include <dune/localfunctions/orthonormal/orthonormalbasis.hh>
22 
23 namespace Dune
24 {
25 
26  // Internal Forward Declarations
27  // -----------------------------
28 
29  template < unsigned int dim, class Field >
30  struct NedelecL2InterpolationFactory;
31 
32 
33 
34  // LocalCoefficientsContainer
35  // --------------------------
36 
37  class LocalCoefficientsContainer
38  {
39  typedef LocalCoefficientsContainer This;
40 
41  public:
42  template <class Setter>
43  LocalCoefficientsContainer ( const Setter &setter )
44  {
45  setter.setLocalKeys(localKey_);
46  }
47 
48  const LocalKey &localKey ( const unsigned int i ) const
49  {
50  assert( i < size() );
51  return localKey_[ i ];
52  }
53 
54  std::size_t size () const
55  {
56  return localKey_.size();
57  }
58 
59  private:
60  std::vector< LocalKey > localKey_;
61  };
62 
63 
64 
65  // NedelecCoefficientsFactory
66  // --------------------------------
67 
68  template < unsigned int dim >
69  struct NedelecCoefficientsFactory
70  {
71  typedef std::size_t Key;
72  typedef const LocalCoefficientsContainer Object;
73 
74  template< GeometryType::Id geometryId >
75  static Object *create( const Key &key )
76  {
77  typedef NedelecL2InterpolationFactory< dim, double > InterpolationFactory;
78  if( !supports< geometryId >( key ) )
79  return nullptr;
80  typename InterpolationFactory::Object *interpolation = InterpolationFactory::template create< geometryId >( key );
81  Object *localKeys = new Object( *interpolation );
82  InterpolationFactory::release( interpolation );
83  return localKeys;
84  }
85 
86  template< GeometryType::Id geometryId >
87  static bool supports ( const Key &key )
88  {
89  GeometryType gt = geometryId;
90  return gt.isTriangle() || gt.isTetrahedron() ;
91  }
92  static void release( Object *object ) { delete object; }
93  };
94 
95 
96 
97  // NedelecL2InterpolationBuilder
98  // ------------------------
99 
100  // L2 Interpolation requires:
101  // - for element
102  // - test basis
103  // - for each face (dynamic)
104  // - test basis
105  // - tangents
106  // - for each edge (dynamic)
107  // - test basis
108  // - tangent
109  template< unsigned int dim, class Field >
110  struct NedelecL2InterpolationBuilder
111  {
112  static const unsigned int dimension = dim;
113 
114  // for the dofs associated to the element
115  typedef OrthonormalBasisFactory< dimension, Field > TestBasisFactory;
116  typedef typename TestBasisFactory::Object TestBasis;
117 
118  // for the dofs associated to the faces
119  typedef OrthonormalBasisFactory< dimension-1, Field > TestFaceBasisFactory;
120  typedef typename TestFaceBasisFactory::Object TestFaceBasis;
121 
122  // for the dofs associated to the edges
123  typedef OrthonormalBasisFactory< 1, Field > TestEdgeBasisFactory;
124  typedef typename TestEdgeBasisFactory::Object TestEdgeBasis;
125 
126  // the tangent of the edges
127  typedef FieldVector< Field, dimension > Tangent;
128 
129  // the normal and the tangents of the faces
130  typedef FieldVector< Field, dimension > Normal;
131  typedef std::array<FieldVector< Field, dimension >,dim-1> FaceTangents;
132 
133  NedelecL2InterpolationBuilder () = default;
134 
135  NedelecL2InterpolationBuilder ( const NedelecL2InterpolationBuilder & ) = delete;
136  NedelecL2InterpolationBuilder ( NedelecL2InterpolationBuilder && ) = delete;
137 
138  ~NedelecL2InterpolationBuilder ()
139  {
140  TestBasisFactory::release( testBasis_ );
141  for( FaceStructure &f : faceStructure_ )
142  TestFaceBasisFactory::release( f.basis_ );
143  for( EdgeStructure& e : edgeStructure_ )
144  TestEdgeBasisFactory::release( e.basis_ );
145  }
146 
147  unsigned int topologyId () const
148  {
149  return geometry_.id();
150  }
151 
152  GeometryType type () const
153  {
154  return geometry_;
155  }
156 
157  std::size_t order () const
158  {
159  return order_;
160  }
161 
162  // number of faces
163  unsigned int faceSize () const
164  {
165  return numberOfFaces_;
166  }
167 
168  // number of edges
169  unsigned int edgeSize () const
170  {
171  return numberOfEdges_;
172  }
173 
174  // basis associated to the element
175  TestBasis *testBasis () const
176  {
177  return testBasis_;
178  }
179 
180  // basis associated to face f
181  TestFaceBasis *testFaceBasis ( unsigned int f ) const
182  {
183  assert( f < faceSize() );
184  return faceStructure_[ f ].basis_;
185  }
186 
187  // basis associated to edge e
188  TestEdgeBasis *testEdgeBasis ( unsigned int e ) const
189  {
190  assert( e < edgeSize() );
191  return edgeStructure_[ e ].basis_;
192  }
193 
194  const Tangent& edgeTangent ( unsigned int e ) const
195  {
196  assert( e < edgeSize() );
197  return edgeStructure_[ e ].tangent_;
198  }
199 
200  const FaceTangents& faceTangents ( unsigned int f ) const
201  {
202  assert( f < faceSize() );
203  return faceStructure_[ f ].faceTangents_;
204  }
205 
206  const Normal &normal ( unsigned int f ) const
207  {
208  assert( f < faceSize() );
209  return faceStructure_[ f ].normal_;
210  }
211 
212  template< GeometryType::Id geometryId >
213  void build ( std::size_t order )
214  {
215  constexpr GeometryType geometry = geometryId;
216  order_ = order;
217  geometry_ = geometry;
218 
219  /*
220  * The Nedelec parameter begins at 1.
221  * This is the numbering used by J.C. Nedelec himself.
222  * See "Mixed Finite Elements in \R^3" published in 1980.
223  *
224  * This construction is based on the construction of Raviart-Thomas elements.
225  * There the numbering starts at 0.
226  * Because of this we reduce the order internally by 1.
227  */
228  order--;
229 
230  // if dimension == 2: order-1 on element
231  // if dimension == 3: order-2 on element
232  int requiredOrder = static_cast<int>(dimension==3);
233  testBasis_ = (order > requiredOrder ? TestBasisFactory::template create< geometry >( order-1-requiredOrder ) : nullptr);
234 
235  const auto &refElement = ReferenceElements< Field, dimension >::general( type() );
236 
237  numberOfFaces_ = refElement.size( 1 );
238  faceStructure_.reserve( numberOfFaces_ );
239 
240  // compute the basis, tangents and normals of each face
241  for (std::size_t i=0; i<numberOfFaces_; i++)
242  {
243  FieldVector<Field,dimension> zero(0);
244  FaceTangents faceTangents;
245  faceTangents.fill(zero);
246 
247  // use the first dim-1 vertices of a face to compute the tangents
248  auto vertices = refElement.subEntities(i,1,dim).begin();
249  auto vertex1 = *vertices;
250  for(int j=1; j<dim;j++)
251  {
252  auto vertex2 = vertices[j];
253 
254  faceTangents[j-1] = refElement.position(vertex2,dim) - refElement.position(vertex1,dim);
255 
256  // By default, edges point from the vertex with the smaller index
257  // to the vertex with the larger index.
258  if (vertex1>vertex2)
259  faceTangents[j-1] *=-1;
260 
261  vertex1 = vertex2;
262  }
263 
264  /* For simplices or cubes of arbitrary dimension you could just use
265  *
266  * ```
267  * GeometryType faceGeometry = Impl::getBase(geometry_);
268  * TestFaceBasis *faceBasis = ( dim == 3 && order > 0 ? TestFaceBasisFactory::template create< faceGeometry >( order-1 ) : nullptr);
269  * ```
270  *
271  * For i.e. Prisms and Pyramids in 3d this does not work because they contain squares and triangles as faces.
272  * And depending on the dynamic face index a different face geometry is needed.
273  *
274  */
275  TestFaceBasis *faceBasis = ( dim == 3 && order > 0 ? Impl::IfGeometryType< CreateFaceBasis, dimension-1 >::apply( refElement.type( i, 1 ), order-1 ) : nullptr);
276  faceStructure_.emplace_back( faceBasis, refElement.integrationOuterNormal(i), faceTangents );
277  }
278  assert( faceStructure_.size() == numberOfFaces_ );
279 
280  numberOfEdges_ = refElement.size( dim-1 );
281  edgeStructure_.reserve( numberOfEdges_ );
282 
283  // compute the basis and tangent of each edge
284  for (std::size_t i=0; i<numberOfEdges_; i++)
285  {
286  auto vertexIterator = refElement.subEntities(i,dim-1,dim).begin();
287  auto v0 = *vertexIterator;
288  auto v1 = *(++vertexIterator);
289 
290  // By default, edges point from the vertex with the smaller index
291  // to the vertex with the larger index.
292  if (v0>v1)
293  std::swap(v0,v1);
294  auto tangent = std::move(refElement.position(v1,dim) - refElement.position(v0,dim));
295 
296  TestEdgeBasis *edgeBasis = Impl::IfGeometryType< CreateEdgeBasis, 1 >::apply( refElement.type( i, dim-1 ), order );
297  edgeStructure_.emplace_back( edgeBasis, tangent );
298  }
299  assert( edgeStructure_.size() == numberOfEdges_ );
300  }
301 
302  private:
303 
304  // helper struct for edges
305  struct EdgeStructure
306  {
307  EdgeStructure( TestEdgeBasis *teb, const Tangent &t )
308  : basis_( teb ), tangent_( t )
309  {}
310 
311  TestEdgeBasis *basis_;
313  };
314 
315  template< GeometryType::Id edgeGeometryId >
316  struct CreateEdgeBasis
317  {
318  static TestEdgeBasis *apply ( std::size_t order ) { return TestEdgeBasisFactory::template create< edgeGeometryId >( order ); }
319  };
320 
321  // helper struct for faces
322  struct FaceStructure
323  {
324  FaceStructure( TestFaceBasis *tfb, const Normal& normal, const FaceTangents& faceTangents )
325  : basis_( tfb ), normal_(normal), faceTangents_( faceTangents )
326  {}
327 
328  TestFaceBasis *basis_;
330  const FaceTangents faceTangents_;
331  };
332 
333  template< GeometryType::Id faceGeometryId >
334  struct CreateFaceBasis
335  {
336  static TestFaceBasis *apply ( std::size_t order ) { return TestFaceBasisFactory::template create< faceGeometryId >( order ); }
337  };
338 
339  TestBasis *testBasis_ = nullptr;
340  std::vector< FaceStructure > faceStructure_;
341  unsigned int numberOfFaces_;
342  std::vector< EdgeStructure > edgeStructure_;
343  unsigned int numberOfEdges_;
344  GeometryType geometry_;
345  std::size_t order_;
346  };
347 
348 
349 
350  // NedelecL2Interpolation
351  // ----------------------------
352 
358  template< unsigned int dimension, class F>
360  : public InterpolationHelper< F ,dimension >
361  {
363  typedef InterpolationHelper<F,dimension> Base;
364 
365  public:
366  typedef F Field;
367  typedef NedelecL2InterpolationBuilder<dimension,Field> Builder;
368  typedef typename Builder::FaceTangents FaceTangents;
369 
371  : order_(0),
372  size_(0)
373  {}
374 
375  template< class Function, class Vector >
376  auto interpolate ( const Function &function, Vector &coefficients ) const
377  -> std::enable_if_t< std::is_same< decltype(std::declval<Vector>().resize(1) ),void >::value,void>
378  {
379  coefficients.resize(size());
380  typename Base::template Helper<Function,Vector,true> func( function,coefficients );
381  interpolate(func);
382  }
383 
384  template< class Basis, class Matrix >
385  auto interpolate ( const Basis &basis, Matrix &matrix ) const
386  -> std::enable_if_t< std::is_same<
387  decltype(std::declval<Matrix>().rowPtr(0)), typename Matrix::Field* >::value,void>
388  {
389  matrix.resize( size(), basis.size() );
390  typename Base::template Helper<Basis,Matrix,false> func( basis,matrix );
391  interpolate(func);
392  }
393 
394  std::size_t order() const
395  {
396  return order_;
397  }
398  std::size_t size() const
399  {
400  return size_;
401  }
402 
403  template <GeometryType::Id geometryId>
404  void build( std::size_t order )
405  {
406  size_ = 0;
407  order_ = order;
408  builder_.template build<geometryId>(order_);
409  if (builder_.testBasis())
410  size_ += dimension*builder_.testBasis()->size();
411 
412  for ( unsigned int f=0; f<builder_.faceSize(); ++f )
413  if (builder_.testFaceBasis(f))
414  size_ += (dimension-1)*builder_.testFaceBasis(f)->size();
415 
416  for ( unsigned int e=0; e<builder_.edgeSize(); ++e )
417  if (builder_.testEdgeBasis(e))
418  size_ += builder_.testEdgeBasis(e)->size();
419  }
420 
421  void setLocalKeys(std::vector< LocalKey > &keys) const
422  {
423  keys.resize(size());
424  unsigned int row = 0;
425  for (unsigned int e=0; e<builder_.edgeSize(); ++e)
426  {
427  if (builder_.edgeSize())
428  for (unsigned int i=0; i<builder_.testEdgeBasis(e)->size(); ++i,++row)
429  keys[row] = LocalKey(e,dimension-1,i);
430  }
431  for (unsigned int f=0; f<builder_.faceSize(); ++f)
432  {
433  if (builder_.faceSize())
434  for (unsigned int i=0; i<builder_.testFaceBasis(f)->size()*(dimension-1); ++i,++row)
435  keys[row] = LocalKey(f,1,i);
436  }
437 
438  if (builder_.testBasis())
439  for (unsigned int i=0; i<builder_.testBasis()->size()*dimension; ++i,++row)
440  keys[row] = LocalKey(0,0,i);
441  assert( row == size() );
442  }
443 
444  protected:
445  template< class Func, class Container, bool type >
446  void interpolate ( typename Base::template Helper<Func,Container,type> &func ) const
447  {
448  const Dune::GeometryType geoType( builder_.topologyId(), dimension );
449 
450  std::vector<Field> testBasisVal;
451 
452  for (unsigned int i=0; i<size(); ++i)
453  for (unsigned int j=0; j<func.size(); ++j)
454  func.set(i,j,0);
455 
456  unsigned int row = 0;
457 
458  // edge dofs:
459  typedef Dune::QuadratureRule<Field, 1> EdgeQuadrature;
460  typedef Dune::QuadratureRules<Field, 1> EdgeQuadratureRules;
461 
462  const auto &refElement = Dune::ReferenceElements< Field, dimension >::general( geoType );
463 
464  for (unsigned int e=0; e<builder_.edgeSize(); ++e)
465  {
466  if (!builder_.testEdgeBasis(e))
467  continue;
468  testBasisVal.resize(builder_.testEdgeBasis(e)->size());
469 
470  const auto &geometry = refElement.template geometry< dimension-1 >( e );
471  const Dune::GeometryType subGeoType( geometry.type().id(), 1 );
472  const EdgeQuadrature &edgeQuad = EdgeQuadratureRules::rule( subGeoType, 2*order_+2 );
473 
474  const unsigned int quadratureSize = edgeQuad.size();
475  for( unsigned int qi = 0; qi < quadratureSize; ++qi )
476  {
477  if (dimension>1)
478  builder_.testEdgeBasis(e)->template evaluate<0>(edgeQuad[qi].position(),testBasisVal);
479  else
480  testBasisVal[0] = 1.;
481  computeEdgeDofs(row,
482  testBasisVal,
483  func.evaluate( geometry.global( edgeQuad[qi].position() ) ),
484  builder_.edgeTangent(e),
485  edgeQuad[qi].weight(),
486  func);
487  }
488 
489  row += builder_.testEdgeBasis(e)->size();
490  }
491 
492  // face dofs:
493  typedef Dune::QuadratureRule<Field, dimension-1> FaceQuadrature;
494  typedef Dune::QuadratureRules<Field, dimension-1> FaceQuadratureRules;
495 
496  for (unsigned int f=0; f<builder_.faceSize(); ++f)
497  {
498  if (builder_.testFaceBasis(f))
499  {
500  testBasisVal.resize(builder_.testFaceBasis(f)->size());
501 
502  const auto &geometry = refElement.template geometry< 1 >( f );
503  const Dune::GeometryType subGeoType( geometry.type().id(), dimension-1 );
504  const FaceQuadrature &faceQuad = FaceQuadratureRules::rule( subGeoType, 2*order_+2 );
505 
506  const unsigned int quadratureSize = faceQuad.size();
507  for( unsigned int qi = 0; qi < quadratureSize; ++qi )
508  {
509  if (dimension>1)
510  builder_.testFaceBasis(f)->template evaluate<0>(faceQuad[qi].position(),testBasisVal);
511  else
512  testBasisVal[0] = 1.;
513 
514  computeFaceDofs( row,
515  testBasisVal,
516  func.evaluate( geometry.global( faceQuad[qi].position() ) ),
517  builder_.faceTangents(f),
518  builder_.normal(f),
519  faceQuad[qi].weight(),
520  func);
521  }
522 
523  row += builder_.testFaceBasis(f)->size()*(dimension-1);
524  }
525  }
526 
527  // element dofs
528  if (builder_.testBasis())
529  {
530  testBasisVal.resize(builder_.testBasis()->size());
531 
532  typedef Dune::QuadratureRule<Field, dimension> Quadrature;
534  const Quadrature &elemQuad = QuadratureRules::rule( geoType, 2*order_+1 );
535 
536  const unsigned int quadratureSize = elemQuad.size();
537  for( unsigned int qi = 0; qi < quadratureSize; ++qi )
538  {
539  builder_.testBasis()->template evaluate<0>(elemQuad[qi].position(),testBasisVal);
540  computeInteriorDofs(row,
541  testBasisVal,
542  func.evaluate(elemQuad[qi].position()),
543  elemQuad[qi].weight(),
544  func );
545  }
546 
547  row += builder_.testBasis()->size()*dimension;
548  }
549  assert(row==size());
550  }
551 
552  private:
562  template <class MVal, class NedVal,class Matrix>
563  void computeEdgeDofs (unsigned int startRow,
564  const MVal &mVal,
565  const NedVal &nedVal,
566  const FieldVector<Field,dimension> &tangent,
567  const Field &weight,
568  Matrix &matrix) const
569  {
570  const unsigned int endRow = startRow+mVal.size();
571  typename NedVal::const_iterator nedIter = nedVal.begin();
572  for ( unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
573  {
574  Field cFactor = (*nedIter)*tangent;
575  typename MVal::const_iterator mIter = mVal.begin();
576  for (unsigned int row = startRow; row!=endRow; ++mIter, ++row )
577  matrix.add(row,col, (weight*cFactor)*(*mIter) );
578 
579  assert( mIter == mVal.end() );
580  }
581  }
582 
593  template <class MVal, class NedVal,class Matrix>
594  void computeFaceDofs (unsigned int startRow,
595  const MVal &mVal,
596  const NedVal &nedVal,
597  const FaceTangents& faceTangents,
598  const FieldVector<Field,dimension> &normal,
599  const Field &weight,
600  Matrix &matrix) const
601  {
602  const unsigned int endRow = startRow+mVal.size()*(dimension-1);
603  typename NedVal::const_iterator nedIter = nedVal.begin();
604  for ( unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
605  {
606  auto const& u=*nedIter;
607  auto const& n=normal;
608  FieldVector<Field,dimension> nedTimesNormal = { u[1]*n[2]-u[2]*n[1],
609  u[2]*n[0]-u[0]*n[2],
610  u[0]*n[1]-u[1]*n[0]};
611  typename MVal::const_iterator mIter = mVal.begin();
612  for (unsigned int row = startRow; row!=endRow; ++mIter)
613  {
614  for(int i=0; i<dimension-1;i++)
615  {
616  auto test = *mIter*faceTangents[i];
617  matrix.add(row+i,col, weight*(nedTimesNormal*test) );
618  }
619  row += dimension-1;
620  }
621 
622  assert( mIter == mVal.end() );
623  }
624  }
625 
634  template <class MVal, class NedVal,class Matrix>
635  void computeInteriorDofs (unsigned int startRow,
636  const MVal &mVal,
637  const NedVal &nedVal,
638  Field weight,
639  Matrix &matrix) const
640  {
641  const unsigned int endRow = startRow+mVal.size()*dimension;
642  typename NedVal::const_iterator nedIter = nedVal.begin();
643  for ( unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
644  {
645  typename MVal::const_iterator mIter = mVal.begin();
646  for (unsigned int row = startRow; row!=endRow; ++mIter,row+=dimension )
647  for (unsigned int i=0; i<dimension; ++i)
648  matrix.add(row+i,col, (weight*(*mIter))*(*nedIter)[i] );
649 
650  assert( mIter == mVal.end() );
651  }
652  }
653 
654  public:
655  Builder builder_;
656  std::size_t order_;
657  std::size_t size_;
658  };
659 
660  template < unsigned int dim, class Field >
661  struct NedelecL2InterpolationFactory
662  {
663  typedef NedelecL2InterpolationBuilder<dim,Field> Builder;
664  typedef const NedelecL2Interpolation<dim,Field> Object;
665  typedef std::size_t Key;
666  typedef typename std::remove_const<Object>::type NonConstObject;
667 
668  template <GeometryType::Id geometryId>
669  static Object *create( const Key &key )
670  {
671  if ( !supports<geometryId>(key) )
672  return 0;
673  NonConstObject *interpol = new NonConstObject();
674  interpol->template build<geometryId>(key);
675  return interpol;
676  }
677 
678  template <GeometryType::Id geometryId>
679  static bool supports( const Key &key )
680  {
681  GeometryType gt = geometryId;
682  return gt.isTriangle() || gt.isTetrahedron() ;
683  }
684  static void release( Object *object ) { delete object; }
685  };
686 
687 } // namespace Dune
688 
689 #endif // #ifndef DUNE_LOCALFUNCTIONS_NEDELEC_NEDELECSIMPLEX_NEDELECSIMPLEXINTERPOLATION_HH
Iterator begin()
begin iterator
Definition: densevector.hh:347
Base class template for function classes.
Definition: function.hh:41
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Describe position of one degree of freedom.
Definition: localkey.hh:23
A generic dynamic dense matrix.
Definition: matrix.hh:561
An L2-based interpolation for Nedelec.
Definition: nedelecsimplexinterpolation.hh:361
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:154
A container for all quadrature rules of dimension dim
Definition: quadraturerules.hh:200
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:266
A few common exception classes.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
Dune namespace.
Definition: alignedallocator.hh:13
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:198
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 4, 22:30, 2024)