PRISMS-PF Manual
Loading...
Searching...
No Matches
dof_manager.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/mg_level_object.h>
7#include <deal.II/dofs/dof_handler.h>
8#include <deal.II/fe/fe_system.h>
9
14#include <prismspf/core/types.h>
15
16#include <prismspf/config.h>
17
18#include <memory>
19
21
25template <unsigned int dim, unsigned int degree>
27{
28public:
32 explicit DoFManager(const std::vector<FieldAttributes> &field_attributes);
33
37 DoFManager(const std::vector<FieldAttributes> &field_attributes,
38 const TriangulationManager<dim> &triangulation_handler);
39
43 void
44 init(const TriangulationManager<dim> &triangulation_handler);
45
49 void
50 reinit(const TriangulationManager<dim> &triangulation_handler);
51
55 [[nodiscard]] std::vector<const dealii::DoFHandler<dim> *>
56 get_field_dof_handlers(const std::set<unsigned int> &field_indices,
57 unsigned int relative_level = 0) const;
58
62 [[nodiscard]] const dealii::DoFHandler<dim> &
63 get_field_dof_handler(Types::Index field_index, unsigned int relative_level = 0) const
64 {
65 return *dof_handlers[field_index][relative_level];
66 }
67
71 const std::vector<std::array<dealii::DoFHandler<dim>, 2>> &
73 {
74 return level_dof_handlers;
75 }
76
80 [[nodiscard]] const dealii::DoFHandler<dim> &
81 get_dof_handler(const unsigned int &rank, unsigned int relative_level = 0) const
82 {
83 return level_dof_handlers[relative_level][rank];
84 }
85
89 [[nodiscard]] dealii::types::global_dof_index
91 {
92 dealii::types::global_dof_index n_dofs = 0;
93 for (const auto &dof_handler_set : dof_handlers)
94 {
95 n_dofs += dof_handler_set[0]->n_dofs();
96 }
97 return n_dofs;
98 }
99
100private:
109 std::vector<std::vector<const dealii::DoFHandler<dim> *>> dof_handlers;
110
114 std::vector<std::array<dealii::DoFHandler<dim>, 2>> level_dof_handlers;
115};
116
117PRISMS_PF_END_NAMESPACE
Class that manages the deal.II DoFHandlers.
Definition dof_manager.h:27
void init(const TriangulationManager< dim > &triangulation_handler)
Initialize the DoFHandlers.
Definition dof_manager.cc:43
dealii::types::global_dof_index get_total_dofs() const
Get the total DoFs excluding multigrid DoFs.
Definition dof_manager.h:90
std::vector< const dealii::DoFHandler< dim > * > get_field_dof_handlers(const std::set< unsigned int > &field_indices, unsigned int relative_level=0) const
Getter function for the DoFHandlers (constant reference).
Definition dof_manager.cc:66
std::vector< std::vector< const dealii::DoFHandler< dim > * > > dof_handlers
Collection of the triangulation DoFs. The number of DoFHandlers should be equal to or less than the n...
Definition dof_manager.h:109
const dealii::DoFHandler< dim > & get_dof_handler(const unsigned int &rank, unsigned int relative_level=0) const
Getter function for a specific scalar or vector DoFHandler.
Definition dof_manager.h:81
void reinit(const TriangulationManager< dim > &triangulation_handler)
Reinitialize the DoFHandlers.
Definition dof_manager.cc:59
const dealii::DoFHandler< dim > & get_field_dof_handler(Types::Index field_index, unsigned int relative_level=0) const
Getter function for the DoFHandler (reference).
Definition dof_manager.h:63
const std::vector< std::array< dealii::DoFHandler< dim >, 2 > > & get_dof_handlers() const
Getter function for the scalar and vector DoFHandlers.
Definition dof_manager.h:72
std::vector< std::array< dealii::DoFHandler< dim >, 2 > > level_dof_handlers
A scalar and a vector dof handler for each level.
Definition dof_manager.h:114
This class handlers the generation and manipulation of triangulations.
Definition triangulation_manager.h:26
Definition conditional_ostreams.cc:20
unsigned int Index
Type for field indices.
Definition types.h:19