This module describes the experimental simplex support in deal.II.
More...
This module describes the experimental simplex support in deal.II.
Simplex and mixed meshes in deal.II are still experimental, i.e., work in progress. Large parts of the library have been ported to be able to operate on such kind of meshes. However, there are still many functions that need to be generalized. You can get a good overview of the ported functionalities by taking a look at the tests in the folder "tests/simplex". In the following, we provide two very basic examples to get started and provide some implementation details.
Example: simplex mesh
The following code shows how to work with simplex meshes:
#include <fstream>
#include <iostream>
class Step3
{
public:
Step3();
void run();
private:
void make_grid();
void setup_system();
void assemble_system();
void solve();
void output_results() const;
};
Step3::Step3()
, fe(2)
, quadrature_formula(3)
{}
void Step3::make_grid()
{
std::cout <<
"Number of active cells: " <<
triangulation.n_active_cells()
<< std::endl;
}
void Step3::setup_system()
{
dof_handler.distribute_dofs(fe);
std::cout << "Number of degrees of freedom: " << dof_handler.n_dofs()
<< std::endl;
sparsity_pattern.copy_from(dsp);
system_matrix.reinit(sparsity_pattern);
solution.reinit(dof_handler.n_dofs());
system_rhs.reinit(dof_handler.n_dofs());
}
void Step3::assemble_system()
{
fe,
quadrature_formula,
const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
for (const auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
cell_rhs = 0;
for (const unsigned int q_index : fe_values.quadrature_point_indices())
{
for (const unsigned int i : fe_values.dof_indices())
for (const unsigned int j : fe_values.dof_indices())
(fe_values.shape_grad(i, q_index) *
fe_values.shape_grad(j, q_index) *
fe_values.JxW(q_index));
for (const unsigned int i : fe_values.dof_indices())
cell_rhs(i) += (fe_values.shape_value(i, q_index) *
1. *
fe_values.JxW(q_index));
}
cell->get_dof_indices(local_dof_indices);
for (const unsigned int i : fe_values.dof_indices())
for (const unsigned int j : fe_values.dof_indices())
system_matrix.add(local_dof_indices[i],
local_dof_indices[j],
for (const unsigned int i : fe_values.dof_indices())
system_rhs(local_dof_indices[i]) += cell_rhs(i);
}
std::map<types::global_dof_index, double> boundary_values;
system_matrix,
solution,
system_rhs);
}
void Step3::solve()
{
}
void Step3::output_results() const
{
std::ofstream output("solution.vtk");
}
void Step3::run()
{
make_grid();
setup_system();
assemble_system();
solve();
output_results();
}
int main()
{
Step3 laplace_problem;
laplace_problem.run();
return 0;
}
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
void read(std::istream &in, Format format=Default)
unsigned int depth_console(const unsigned int n)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
bool write_higher_order_cells
void set_flags(const FlagType &flags)
void write_vtk(std::ostream &out) const
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Example: mixed mesh
The following code shows how to work with mixed meshes:
#include <fstream>
#include <iostream>
class Step3
{
public:
Step3();
void run();
private:
void make_grid();
void setup_system();
void assemble_system();
void solve();
void output_results() const;
};
Step3::Step3()
{}
void Step3::make_grid()
{
std::cout <<
"Number of active cells: " <<
triangulation.n_active_cells()
<< std::endl;
}
void Step3::setup_system()
{
for (const auto &cell : dof_handler.active_cell_iterators())
{
cell->set_active_fe_index(0);
cell->set_active_fe_index(1);
else
Assert(
false, ExcNotImplemented());
}
dof_handler.distribute_dofs(fe);
std::cout << "Number of degrees of freedom: " << dof_handler.n_dofs()
<< std::endl;
sparsity_pattern.copy_from(dsp);
system_matrix.reinit(sparsity_pattern);
solution.reinit(dof_handler.n_dofs());
system_rhs.reinit(dof_handler.n_dofs());
}
void Step3::assemble_system()
{
fe,
quadrature_formula,
std::vector<types::global_dof_index> local_dof_indices;
for (const auto &cell : dof_handler.active_cell_iterators())
{
const auto &fe_values = hp_fe_values.get_present_fe_values();
const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
cell_rhs.
reinit(dofs_per_cell);
local_dof_indices.resize(dofs_per_cell);
cell_rhs = 0;
for (const unsigned int q_index : fe_values.quadrature_point_indices())
{
for (const unsigned int i : fe_values.dof_indices())
for (const unsigned int j : fe_values.dof_indices())
(fe_values.shape_grad(i, q_index) *
fe_values.shape_grad(j, q_index) *
fe_values.JxW(q_index));
for (const unsigned int i : fe_values.dof_indices())
cell_rhs(i) += (fe_values.shape_value(i, q_index) *
1. *
fe_values.JxW(q_index));
}
cell->get_dof_indices(local_dof_indices);
for (const unsigned int i : fe_values.dof_indices())
for (const unsigned int j : fe_values.dof_indices())
system_matrix.
add(local_dof_indices[i],
local_dof_indices[j],
for (const unsigned int i : fe_values.dof_indices())
system_rhs(local_dof_indices[i]) += cell_rhs(i);
}
std::map<types::global_dof_index, double> boundary_values;
system_matrix,
solution,
system_rhs);
}
void Step3::solve()
{
}
void Step3::output_results() const
{
std::ofstream output("solution.vtk");
}
void Step3::run()
{
make_grid();
setup_system();
assemble_system();
solve();
output_results();
}
int main()
{
Step3 laplace_problem;
laplace_problem.run();
return 0;
}
#define Assert(cond, exc)
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Triangle
Reference cells
In 2D, we provide triangles and quadrilaterals with the following possible orientations in 3D:
2D: triangle and quadrilateral
Possible orientations of triangles and quadrilaterals in 3D
In 3D, tetrahedra, pyramids, wedges, and hexahedra are available:
Each surface of a 3D reference cell consists of 2D reference cells. The documentation of the enumeration of the numbering of their vertices and lines are given in the right columns.
◆ implicit_function()
Generate a Triangulation from the zero level set of an implicit function, using the CGAL library.
This function is only implemented for dim
equal to two or three, and requires that deal.II is configured using DEAL_II_WITH_CGAL
. When dim
is equal to three, the implicit_function
is supposed to be negative in the interior of the domain, positive outside, and to be entirely enclosed in a ball of radius outer_ball_radius
centered at the point interior_point
. The triangulation that is generated covers the volume bounded by the zero level set of the implicit function where the implicit_function
is negative.
When dim
is equal to two, the generated surface triangulation is the zero level set of the implicit_function
, oriented such that the surface triangulation has normals pointing towards the region where implicit_function
is positive.
The struct data
can be used to pass additional arguments to the CGAL::Mesh_criteria_3 class (see https://doc.cgal.org/latest/Mesh_3/index.html for more information.)
An example usage of this function is given by
Point<3>(.5, 0, 0), 1.0, cell_size = 0.2);
void implicit_function(Triangulation< dim, 3 > &tria, const Function< 3 > &implicit_function, const CGALWrappers::AdditionalData< dim > &data=CGALWrappers::AdditionalData< dim >{}, const Point< 3 > &interior_point=Point< 3 >(), const double &outer_ball_radius=1.0)
const ::Triangulation< dim, spacedim > & tria
The above snippet of code generates the following grid for dim
equal to two and three respectively
- Parameters
-
[out] | tria | The output triangulation |
[in] | implicit_function | The implicit function |
[in] | data | Additional parameters to pass to the CGAL::make_mesh_3 function and to the CGAL::make_surface_mesh functions |
[in] | interior_point | A point in the interior of the domain, for which implicit_function is negative |
[in] | outer_ball_radius | The radius of the ball that will contain the generated Triangulation object |