4#include <deal.II/base/point.h>
5#include <deal.II/dofs/dof_handler.h>
6#include <deal.II/dofs/dof_renumbering.h>
7#include <deal.II/dofs/dof_tools.h>
8#include <deal.II/fe/fe_dgq.h>
9#include <deal.II/fe/fe_system.h>
10#include <deal.II/fe/mapping_q1.h>
11#include <deal.II/lac/affine_constraints.h>
21 using namespace dealii;
29 template <
typename Components_,
typename NumberType_,
typename Mesh_>
class Discretization
41 static_assert(Components::count_fe_subsystems() > 1,
42 "LDG must have a defined submodel of the Model with index 1.");
43 for (
uint i = 0; i < Components::count_fe_subsystems(); ++i) {
44 fe.emplace_back(std::make_shared<FESystem<dim>>(FE_DGQ<dim>(
json.
get_uint(
"/discretization/fe_order")),
45 Components::count_fe_functions(i)));
46 dof_handler.emplace_back(std::make_shared<DoFHandler<dim>>(
mesh.get_triangulation()));
62 std::vector<const DoFHandler<dim> *> ret;
63 for (
uint i = 0; i < Components::count_fe_subsystems(); ++i)
80 double min_dist = std::numeric_limits<double>::max();
83 if (dist < min_dist) {
93 std::vector<uint> block_structure{
dof_handler[0]->n_dofs()};
94 if (Components::count_variables() > 0) block_structure.push_back(Components::count_variables());
95 return block_structure;
101 spdlog::get(
"log")->info(
"FEM: Number of active cells: {}",
mesh.get_triangulation().n_active_cells());
103 for (
uint i = 0; i < Components::count_fe_subsystems(); ++i) {
110 if (i == 0) spdlog::get(
"log")->info(
"FEM: Number of degrees of freedom: {}",
dof_handler[0]->n_dofs());
120 std::vector<std::shared_ptr<FESystem<dim>>>
fe;
A wrapper around the boost json value class.
Definition json.hh:19
uint get_uint(const std::string &key) const
Get the value of a key in the json object.
Class to manage the system on which we solve, i.e. fe spaces, grids, etc. This class is a System for ...
Definition ldg.hh:30
uint get_closest_dof(const Point< dim > &p) const
Definition ldg.hh:77
const auto & get_constraints(const uint i=0) const
Definition ldg.hh:52
const auto & get_support_points() const
Definition ldg.hh:72
std::vector< uint > get_block_structure() const
Definition ldg.hh:91
NumberType_ NumberType
Definition ldg.hh:33
std::vector< std::shared_ptr< FESystem< dim > > > fe
Definition ldg.hh:120
Components_ Components
Definition ldg.hh:32
const auto get_dof_handler_list() const
Definition ldg.hh:60
Mesh_ Mesh
Definition ldg.hh:36
const auto & get_mapping() const
Definition ldg.hh:68
std::vector< Point< dim > > support_points
Definition ldg.hh:124
void setup_dofs()
Definition ldg.hh:99
void reinit()
Definition ldg.hh:75
static constexpr uint dim
Definition ldg.hh:37
const auto & get_json() const
Definition ldg.hh:73
const auto & get_dof_handler(const uint i=0) const
Definition ldg.hh:58
auto & get_dof_handler(const uint i=0)
Definition ldg.hh:59
const auto & get_triangulation() const
Definition ldg.hh:69
auto & get_triangulation()
Definition ldg.hh:70
std::vector< AffineConstraints< NumberType > > constraints
Definition ldg.hh:122
auto & get_constraints(const uint i=0)
Definition ldg.hh:53
const auto & get_fe(uint i=0) const
Definition ldg.hh:67
Vector< NumberType > VectorType
Definition ldg.hh:34
MappingQ1< dim > mapping
Definition ldg.hh:123
JSONValue json
Definition ldg.hh:118
Discretization(Mesh &mesh, const JSONValue &json)
Definition ldg.hh:39
BlockSparseMatrix< NumberType > SparseMatrixType
Definition ldg.hh:35
std::vector< std::shared_ptr< DoFHandler< dim > > > dof_handler
Definition ldg.hh:121
Mesh & mesh
Definition ldg.hh:117
const Point< dim > & get_support_point(const uint &dof) const
Definition ldg.hh:71
Definition complex_math.hh:14
unsigned int uint
Definition utils.hh:22