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