ExaDG
Loading...
Searching...
No Matches
container_interface_data.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_FUNCTIONS_AND_BOUNDARY_CONDITIONS_CONTAINER_INTERFACE_DATA_H_
23#define INCLUDE_EXADG_FUNCTIONS_AND_BOUNDARY_CONDITIONS_CONTAINER_INTERFACE_DATA_H_
24
25// deal.II
26#include <deal.II/base/tensor.h>
27
28#include <map>
29#include <tuple>
30#include <vector>
31
32// ExaDG
33#include <exadg/matrix_free/integrators.h>
34#include <exadg/utilities/tensor_utilities.h>
35
36namespace ExaDG
37{
43template<int rank, int dim, typename number_type>
45{
46public:
47 typedef dealii::Tensor<rank, dim, number_type> data_type;
48
49private:
50 static unsigned int const n_components = rank_to_n_components<rank, dim>();
51
52 using SetBoundaryIDs = std::set<dealii::types::boundary_id>;
53
54 using quad_index = unsigned int;
55
56 using Id = std::tuple<unsigned int /*face*/, unsigned int /*q*/, unsigned int /*v*/>;
57
58 using MapVectorIndex = std::map<Id, dealii::types::global_dof_index>;
59
60 using ArrayQuadraturePoints = std::vector<dealii::Point<dim>>;
61
62 using ArraySolutionValues = std::vector<data_type>;
63
64public:
66
67 template<typename Number>
68 void
69 setup(dealii::MatrixFree<dim, Number> const & matrix_free_,
70 unsigned int const dof_index_,
71 std::vector<quad_index> const & quad_indices_,
72 SetBoundaryIDs const & set_bids_)
73 {
74 quad_indices = quad_indices_;
75
76 for(auto q_index : quad_indices)
77 {
78 // initialize maps
79 map_vector_index.emplace(q_index, MapVectorIndex());
80 map_q_points.emplace(q_index, ArrayQuadraturePoints());
81 map_solution.emplace(q_index, ArraySolutionValues());
82
83 MapVectorIndex & map_index = map_vector_index.find(q_index)->second;
84 ArrayQuadraturePoints & array_q_points_dst = map_q_points.find(q_index)->second;
85 ArraySolutionValues & array_solution_dst = map_solution.find(q_index)->second;
86
87 // create map "ID = {face, q, v} <-> vector_index" and fill array of quadrature points
88 for(unsigned int face = matrix_free_.n_inner_face_batches();
89 face < matrix_free_.n_inner_face_batches() + matrix_free_.n_boundary_face_batches();
90 ++face)
91 {
92 // only consider relevant boundary IDs
93 if(set_bids_.find(matrix_free_.get_boundary_id(face)) != set_bids_.end())
94 {
95 FaceIntegrator<dim, n_components, Number> integrator(matrix_free_,
96 true,
97 dof_index_,
98 q_index);
99 integrator.reinit(face);
100
101 for(unsigned int q = 0; q < integrator.n_q_points; ++q)
102 {
103 dealii::Point<dim, dealii::VectorizedArray<Number>> q_points =
104 integrator.quadrature_point(q);
105
106 for(unsigned int v = 0; v < dealii::VectorizedArray<Number>::size(); ++v)
107 {
108 dealii::Point<dim> q_point;
109 for(unsigned int d = 0; d < dim; ++d)
110 q_point[d] = q_points[d][v];
111
112 Id id = std::make_tuple(face, q, v);
113 dealii::types::global_dof_index index = array_q_points_dst.size();
114 map_index.emplace(id, index);
115 array_q_points_dst.push_back(q_point);
116 }
117 }
118 }
119 }
120
121 array_solution_dst.resize(array_q_points_dst.size(), data_type());
122 }
123 }
124
125 std::vector<quad_index> const &
126 get_quad_indices();
127
128 ArrayQuadraturePoints &
129 get_array_q_points(quad_index const & q_index);
130
131 ArraySolutionValues &
132 get_array_solution(quad_index const & q_index);
133
134 data_type
135 get_data(unsigned int const q_index,
136 unsigned int const face,
137 unsigned int const q,
138 unsigned int const v) const;
139
140private:
141 std::vector<quad_index> quad_indices;
142
143 mutable std::map<quad_index, MapVectorIndex> map_vector_index;
144 mutable std::map<quad_index, ArrayQuadraturePoints> map_q_points;
145 mutable std::map<quad_index, ArraySolutionValues> map_solution;
146};
147} // namespace ExaDG
148
149
150
151#endif /* INCLUDE_EXADG_FUNCTIONS_AND_BOUNDARY_CONDITIONS_CONTAINER_INTERFACE_DATA_H_ */
Definition container_interface_data.h:45
Definition driver.cpp:33