/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FV/discretization.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FV/discretization.hh Source File
DiFfRG
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
17
18namespace DiFfRG
19{
20 namespace FV
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>(0), Components::count_fe_functions(0))),
44 {
45 setup_dofs();
46 }
47
48 const auto &get_constraints(const uint i = 0) const
49 {
50 (void)i;
51 return constraints;
52 }
53 auto &get_constraints(const uint i = 0)
54 {
55 (void)i;
56 return constraints;
57 }
58 const auto &get_dof_handler(const uint i = 0) const
59 {
60 (void)i;
61 return dof_handler;
62 }
63 auto &get_dof_handler(const uint i = 0)
64 {
65 (void)i;
66 return dof_handler;
67 }
68 const auto &get_fe(uint i = 0) const
69 {
70 if (i != 0) throw std::runtime_error("Wrong FE index");
71 return *fe;
72 }
73 const auto &get_mapping() const { return mapping; }
74 const auto &get_triangulation() const { return mesh.get_triangulation(); }
75 auto &get_triangulation() { return mesh.get_triangulation(); }
76 const Point<dim> &get_support_point(const uint &dof) const { return support_points[dof]; }
77 const auto &get_support_points() const { return support_points; }
78 const auto &get_json() const { return json; }
79
80 void reinit() { setup_dofs(); }
81
82 uint get_closest_dof(const Point<dim> &p) const
83 {
84 uint dof = 0;
85 double min_dist = std::numeric_limits<double>::max();
86 for (uint i = 0; i < support_points.size(); ++i) {
87 const auto dist = p.distance(support_points[i]);
88 if (dist < min_dist) {
89 min_dist = dist;
90 dof = i;
91 }
92 }
93 return dof;
94 }
95
96 std::vector<uint> get_block_structure() const
97 {
98 std::vector<uint> block_structure{dof_handler.n_dofs()};
99 if (Components::count_variables() > 0) block_structure.push_back(Components::count_variables());
100 return block_structure;
101 }
102
103 protected:
105 {
106 dof_handler.distribute_dofs(*fe);
107
108 spdlog::get("log")->info("FV: Number of active cells: {}", mesh.get_triangulation().n_active_cells());
109 spdlog::get("log")->info("FV: Number of degrees of freedom: {}", dof_handler.n_dofs());
110
111 constraints.clear();
112 DoFTools::make_hanging_node_constraints(dof_handler, constraints);
113 constraints.close();
114
115 support_points.resize(dof_handler.n_dofs());
116 DoFTools::map_dofs_to_support_points(mapping, dof_handler, support_points);
117 }
118
121
122 std::shared_ptr<FESystem<dim>> fe;
123 DoFHandler<dim> dof_handler;
124 AffineConstraints<NumberType> constraints;
125 MappingQ1<dim> mapping;
126 std::vector<Point<dim>> support_points;
127 };
128 } // namespace FV
129} // 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:31
auto & get_dof_handler(const uint i=0)
Definition discretization.hh:63
MappingQ1< dim > mapping
Definition discretization.hh:125
std::shared_ptr< FESystem< dim > > fe
Definition discretization.hh:122
Mesh_ Mesh
Definition discretization.hh:37
std::vector< uint > get_block_structure() const
Definition discretization.hh:96
Mesh & mesh
Definition discretization.hh:119
JSONValue json
Definition discretization.hh:120
const auto & get_triangulation() const
Definition discretization.hh:74
SparseMatrix< NumberType > SparseMatrixType
Definition discretization.hh:36
void reinit()
Definition discretization.hh:80
const Point< dim > & get_support_point(const uint &dof) const
Definition discretization.hh:76
void setup_dofs()
Definition discretization.hh:104
const auto & get_dof_handler(const uint i=0) const
Definition discretization.hh:58
const auto & get_support_points() const
Definition discretization.hh:77
const auto & get_mapping() const
Definition discretization.hh:73
AffineConstraints< NumberType > constraints
Definition discretization.hh:124
uint get_closest_dof(const Point< dim > &p) const
Definition discretization.hh:82
std::vector< Point< dim > > support_points
Definition discretization.hh:126
static constexpr uint dim
Definition discretization.hh:38
Discretization(Mesh &mesh, const JSONValue &json)
Definition discretization.hh:40
DoFHandler< dim > dof_handler
Definition discretization.hh:123
Vector< NumberType > VectorType
Definition discretization.hh:35
NumberType_ NumberType
Definition discretization.hh:34
auto & get_constraints(const uint i=0)
Definition discretization.hh:53
const auto & get_constraints(const uint i=0) const
Definition discretization.hh:48
const auto & get_fe(uint i=0) const
Definition discretization.hh:68
const auto & get_json() const
Definition discretization.hh:78
Components_ Components
Definition discretization.hh:33
auto & get_triangulation()
Definition discretization.hh:75
A wrapper around the boost json value class.
Definition json.hh:19
Definition complex_math.hh:10
unsigned int uint
Definition utils.hh:24