DiFfRG
Loading...
Searching...
No Matches
dg.hh
Go to the documentation of this file.
1#pragma once
2
3// external libraries
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>
12
13// DiFfRG
17
18namespace DiFfRG
19{
20 namespace DG
21 {
22 using namespace dealii;
23
30 template <typename Components_, typename NumberType_, typename Mesh_> class Discretization
31 {
32 public:
33 using Components = Components_;
34 using NumberType = NumberType_;
35 using VectorType = Vector<NumberType>;
36 using SparseMatrixType = SparseMatrix<NumberType>;
37 using Mesh = Mesh_;
38 static constexpr uint dim = Mesh::dim;
39
41 : mesh(mesh), json(json),
42 fe(std::make_shared<FESystem<dim>>(FE_DGQ<dim>(json.get_uint("/discretization/fe_order")),
43 Components::count_fe_functions(0))),
45 {
46 setup_dofs();
47 };
48
49 const auto &get_constraints(const uint i = 0) const
50 {
51 (void)i;
52 return constraints;
53 }
54 auto &get_constraints(const uint i = 0)
55 {
56 (void)i;
57 return constraints;
58 }
59 const auto &get_dof_handler(const uint i = 0) const
60 {
61 (void)i;
62 return dof_handler;
63 }
64 auto &get_dof_handler(const uint i = 0)
65 {
66 (void)i;
67 return dof_handler;
68 }
69 const auto &get_fe(uint i = 0) const
70 {
71 if (i != 0) throw std::runtime_error("Wrong FE index");
72 return *fe;
73 }
74 const auto &get_mapping() const { return mapping; }
75 const auto &get_triangulation() const { return mesh.get_triangulation(); }
76 auto &get_triangulation() { return mesh.get_triangulation(); }
77 const Point<dim> &get_support_point(const uint &dof) const { return support_points[dof]; }
78 const auto &get_support_points() const { return support_points; }
79 const auto &get_json() const { return json; }
80
81 void reinit() { setup_dofs(); }
82
83 uint get_closest_dof(const Point<dim> &p) const
84 {
85 uint dof = 0;
86 double min_dist = std::numeric_limits<double>::max();
87 for (uint i = 0; i < support_points.size(); ++i) {
88 const auto dist = p.distance(support_points[i]);
89 if (dist < min_dist) {
90 min_dist = dist;
91 dof = i;
92 }
93 }
94 return dof;
95 }
96
97 std::vector<uint> get_block_structure() const
98 {
99 std::vector<uint> block_structure{dof_handler.n_dofs()};
100 if (Components::count_variables() > 0) block_structure.push_back(Components::count_variables());
101 return block_structure;
102 }
103
104 protected:
106 {
107 dof_handler.distribute_dofs(*fe);
108
109 spdlog::get("log")->info("FEM: Number of active cells: {}", mesh.get_triangulation().n_active_cells());
110 spdlog::get("log")->info("FEM: Number of degrees of freedom: {}", dof_handler.n_dofs());
111
112 constraints.clear();
113 DoFTools::make_hanging_node_constraints(dof_handler, constraints);
114 constraints.close();
115
116 support_points.resize(dof_handler.n_dofs());
117 DoFTools::map_dofs_to_support_points(mapping, dof_handler, support_points);
118 }
119
122
123 std::shared_ptr<FESystem<dim>> fe;
124 DoFHandler<dim> dof_handler;
125 AffineConstraints<NumberType> constraints;
126 MappingQ1<dim> mapping;
127 std::vector<Point<dim>> support_points;
128 };
129 } // namespace DG
130} // namespace DiFfRG
Class to manage the system on which we solve, i.e. fe spaces, grids, etc. This class is a System for ...
Definition dg.hh:31
static constexpr uint dim
Definition dg.hh:38
AffineConstraints< NumberType > constraints
Definition dg.hh:125
DoFHandler< dim > dof_handler
Definition dg.hh:124
std::shared_ptr< FESystem< dim > > fe
Definition dg.hh:123
const auto & get_triangulation() const
Definition dg.hh:75
auto & get_constraints(const uint i=0)
Definition dg.hh:54
Mesh_ Mesh
Definition dg.hh:37
std::vector< Point< dim > > support_points
Definition dg.hh:127
void setup_dofs()
Definition dg.hh:105
SparseMatrix< NumberType > SparseMatrixType
Definition dg.hh:36
Vector< NumberType > VectorType
Definition dg.hh:35
const auto & get_fe(uint i=0) const
Definition dg.hh:69
Discretization(Mesh &mesh, const JSONValue &json)
Definition dg.hh:40
const auto & get_constraints(const uint i=0) const
Definition dg.hh:49
void reinit()
Definition dg.hh:81
const auto & get_dof_handler(const uint i=0) const
Definition dg.hh:59
const auto & get_mapping() const
Definition dg.hh:74
NumberType_ NumberType
Definition dg.hh:34
auto & get_triangulation()
Definition dg.hh:76
const auto & get_json() const
Definition dg.hh:79
const Point< dim > & get_support_point(const uint &dof) const
Definition dg.hh:77
uint get_closest_dof(const Point< dim > &p) const
Definition dg.hh:83
auto & get_dof_handler(const uint i=0)
Definition dg.hh:64
const auto & get_support_points() const
Definition dg.hh:78
std::vector< uint > get_block_structure() const
Definition dg.hh:97
MappingQ1< dim > mapping
Definition dg.hh:126
Components_ Components
Definition dg.hh:33
Mesh & mesh
Definition dg.hh:120
JSONValue json
Definition dg.hh:121
A wrapper around the boost json value class.
Definition json.hh:19
Definition complex_math.hh:14
unsigned int uint
Definition utils.hh:22
Definition tuples.hh:117