Doxygen API reference documentation for ideal.II
Loading...
Searching...
No Matches
spacetime_tria.hh
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2022 - 2023 by the ideal.II authors
4//
5// This file is part of the ideal.II library.
6//
7// The ideal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 3.0 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of ideal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef INCLUDE_IDEAL_II_GRID_SPACETIME_TRIA_HH_
17#define INCLUDE_IDEAL_II_GRID_SPACETIME_TRIA_HH_
18
19#include <ideal.II/grid/slab_tria.hh>
20
21#include <deal.II/grid/tria.h>
22
23#include <list>
24#include <memory>
25
26namespace idealii::spacetime
27{
35 template <int dim>
37 {
38 public:
44 dealii::types::global_cell_index max_N_intervals_per_slab = 0);
45
54 virtual void
55 generate(std::shared_ptr<dealii::Triangulation<dim>> space_tria,
56 unsigned int M,
57 double t0 = 0.,
58 double T = 1.) = 0;
59
63 unsigned int
64 M();
74 end();
75
81 virtual void
82 refine_global(const unsigned int times_space = 1,
83 const unsigned int times_time = 1) = 0;
84
85 protected:
86 dealii::types::global_cell_index max_N_intervals_per_slab;
87 std::list<slab::Triangulation<dim>> trias;
88 };
89} // namespace idealii::spacetime
90
91#endif /* INCLUDE_IDEAL_II_GRID_FIXED_TRIA_HH_ */
The base class for quadrature formulae in space and time.
Definition spacetime_quadrature.hh:35
The spacetime triangulation object.
Definition spacetime_tria.hh:37
slab::TriaIterator< dim > end()
An iterator pointing behind the first slab::Triangulation.
unsigned int M()
Return the number of slabs in the triangulation.
Triangulation(dealii::types::global_cell_index max_N_intervals_per_slab=0)
Constructor that initializes the underlying list object.
slab::TriaIterator< dim > begin()
An iterator pointing to the first slab::Triangulation.
virtual void refine_global(const unsigned int times_space=1, const unsigned int times_time=1)=0
Do uniform mesh refinement in time and space.
virtual void generate(std::shared_ptr< dealii::Triangulation< dim > > space_tria, unsigned int M, double t0=0., double T=1.)=0
Generate a list of M slab triangulations with matching temporal meshes and space_tria.
typename std::list< Triangulation< dim > >::iterator TriaIterator
A shortened type for Iterators over a list of shared pointers to Triangulation<dim> objects.
Definition slab_tria.hh:125
Namespace for general spacetime object and collections of slab objects.
Definition idealii.hh:89