16#ifndef dealii_cuda_hanging_nodes_internal_h
17#define dealii_cuda_hanging_nodes_internal_h
21#ifdef DEAL_II_COMPILER_CUDA_AWARE
60 template <
typename CellIterator>
63 std::vector<types::global_dof_index> & dof_indices,
64 const CellIterator & cell,
65 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
66 unsigned int & mask)
const;
80 unsigned int n_dofs_1d,
81 std::vector<types::global_dof_index> &dofs)
const;
86 unsigned int n_dofs_1d)
const;
98 std::vector<std::vector<std::pair<cell_iterator, unsigned int>>>
143 const unsigned int face_no =
151 interpolation_matrix,
154 std::vector<unsigned int> mapping =
155 FETools::lexicographic_to_hierarchic_numbering<1>(fe_degree);
161 mapped_matrix(i, j) = interpolation_matrix(mapping[i], mapping[j]);
163 cudaError_t error_code =
165 &mapped_matrix[0][0],
176 unsigned int fe_degree,
178 const std::vector<unsigned int> &lexicographic_mapping)
179 : n_raw_lines(dof_handler.get_triangulation().n_raw_lines())
180 , line_to_cells(dim == 3 ? n_raw_lines : 0)
181 , lexicographic_mapping(lexicographic_mapping)
182 , fe_degree(fe_degree)
183 , dof_handler(dof_handler)
188 internal::setup_constraint_weigths<dim>(
fe_degree);
210 const unsigned int line_to_children[12][2] = {{0, 2},
223 std::vector<std::vector<std::pair<cell_iterator, unsigned int>>>
224 line_to_inactive_cells(n_raw_lines);
227 for (
const auto &cell : dof_handler.cell_iterators())
229 for (
unsigned int line = 0; line < GeometryInfo<3>::lines_per_cell;
232 const unsigned int line_idx = cell->line(line)->index();
233 if (cell->is_active())
234 line_to_cells[line_idx].push_back(std::make_pair(cell, line));
236 line_to_inactive_cells[line_idx].push_back(
237 std::make_pair(cell, line));
244 for (
unsigned int line_idx = 0; line_idx < n_raw_lines; ++line_idx)
246 if ((line_to_cells[line_idx].size() > 0) &&
247 line_to_inactive_cells[line_idx].size() > 0)
252 line_to_inactive_cells[line_idx][0].first;
253 const unsigned int neighbor_line =
254 line_to_inactive_cells[line_idx][0].second;
256 for (
unsigned int c = 0; c < 2; ++c)
259 inactive_cell->child(line_to_children[neighbor_line][c]);
260 const unsigned int child_line_idx =
261 child->line(neighbor_line)->index();
264 for (
const auto cl : line_to_cells[line_idx])
265 line_to_cells[child_line_idx].push_back(cl);
274 template <
typename CellIterator>
277 std::vector<types::global_dof_index> & dof_indices,
278 const CellIterator & cell,
279 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
280 unsigned int & mask)
const
283 const unsigned int n_dofs_1d = fe_degree + 1;
284 const unsigned int dofs_per_face =
287 std::vector<types::global_dof_index> neighbor_dofs(dofs_per_face);
289 const auto lex_face_mapping =
294 if ((!cell->at_boundary(face)) &&
295 (cell->neighbor(face)->has_children() ==
false))
300 if (neighbor->level() < cell->level())
302 const unsigned int neighbor_face =
303 cell->neighbor_face_no(face);
306 unsigned int subface = 0;
307 for (; subface < GeometryInfo<dim>::max_children_per_face;
309 if (neighbor->neighbor_child_on_subface(neighbor_face,
314 neighbor->face(neighbor_face)->get_dof_indices(neighbor_dofs);
318 for (
auto &index : neighbor_dofs)
319 index = partitioner->global_to_local(index);
344 unsigned int offset = (face % 2 == 1) ? fe_degree : 0;
346 for (
unsigned int i = 0; i < n_dofs_1d; ++i)
348 unsigned int idx = 0;
351 idx = n_dofs_1d * offset + i;
354 idx = n_dofs_1d * i + offset;
356 dof_indices[idx] = neighbor_dofs[lex_face_mapping[i]];
361 const bool transpose = !(cell->face_orientation(face));
365 if (cell->face_rotation(face))
367 if (cell->face_flip(face))
370 rotate_face(
rotate, n_dofs_1d, neighbor_dofs);
371 rotate_subface_index(
rotate, subface);
375 transpose_face(neighbor_dofs);
376 transpose_subface_index(subface);
385 if (subface % 2 == 0)
387 if (subface / 2 == 0)
396 if (subface % 2 == 0)
398 if (subface / 2 == 0)
407 if (subface % 2 == 0)
409 if (subface / 2 == 0)
414 unsigned int offset = (face % 2 == 1) ? fe_degree : 0;
416 for (
unsigned int i = 0; i < n_dofs_1d; ++i)
418 for (
unsigned int j = 0; j < n_dofs_1d; ++j)
420 unsigned int idx = 0;
424 idx = n_dofs_1d * n_dofs_1d * i +
425 n_dofs_1d * j + offset;
429 idx = n_dofs_1d * n_dofs_1d * j +
430 n_dofs_1d * offset + i;
434 idx = n_dofs_1d * n_dofs_1d * offset +
438 neighbor_dofs[lex_face_mapping[n_dofs_1d * i +
456 const unsigned int line_to_edge[12][4] = {
506 for (
unsigned int local_line = 0;
507 local_line < GeometryInfo<dim>::lines_per_cell;
511 if (!(mask & line_to_edge[local_line][0]))
514 const unsigned int line = cell->line(local_line)->index();
515 for (
const auto edge_neighbor : line_to_cells[line])
519 if (neighbor_cell->level() < cell->level())
521 const unsigned int local_line_neighbor =
522 edge_neighbor.second;
523 mask |= line_to_edge[local_line][1] |
524 line_to_edge[local_line][2];
526 bool flipped =
false;
527 if (cell->line(local_line)->vertex_index(0) ==
528 neighbor_cell->line(local_line_neighbor)
533 mask |= line_to_edge[local_line][3];
535 else if (cell->line(local_line)->vertex_index(1) ==
536 neighbor_cell->line(local_line_neighbor)
541 else if (cell->line(local_line)->vertex_index(1) ==
542 neighbor_cell->line(local_line_neighbor)
548 else if (cell->line(local_line)->vertex_index(0) ==
549 neighbor_cell->line(local_line_neighbor)
553 mask |= line_to_edge[local_line][3];
560 neighbor_dofs.resize(n_dofs_1d * n_dofs_1d *
562 neighbor_cell->get_dof_indices(neighbor_dofs);
566 for (
auto &index : neighbor_dofs)
567 index = partitioner->global_to_local(index);
569 for (
unsigned int i = 0; i < n_dofs_1d; ++i)
572 const unsigned int idx =
573 line_dof_idx(local_line, i, n_dofs_1d);
574 dof_indices[idx] = neighbor_dofs
575 [lexicographic_mapping[line_dof_idx(
577 flipped ? fe_degree - i : i,
595 unsigned int &subface_index)
const
597 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
600 times = times < 0 ? times + 4 : times;
601 for (
int t = 0; t < times; ++t)
602 subface_index = rot_mapping[subface_index];
611 unsigned int n_dofs_1d,
612 std::vector<types::global_dof_index> &dofs)
const
614 const unsigned int rot_mapping[4] = {2, 0, 3, 1};
617 times = times < 0 ? times + 4 : times;
619 std::vector<types::global_dof_index>
copy(dofs.size());
620 for (
int t = 0; t < times; ++t)
625 for (
unsigned int i = 0; i < 4; ++i)
626 dofs[rot_mapping[i]] =
copy[i];
629 const unsigned int n_int = n_dofs_1d - 2;
630 unsigned int offset = 4;
631 for (
unsigned int i = 0; i < n_int; ++i)
634 dofs[offset + i] =
copy[offset + 2 * n_int + (n_int - 1 - i)];
636 dofs[offset + n_int + i] =
637 copy[offset + 3 * n_int + (n_int - 1 - i)];
639 dofs[offset + 2 * n_int + i] =
copy[offset + n_int + i];
641 dofs[offset + 3 * n_int + i] =
copy[offset + i];
647 for (
unsigned int i = 0; i < n_int; ++i)
648 for (
unsigned int j = 0; j < n_int; ++j)
649 dofs[offset + i * n_int + j] =
650 copy[offset + j * n_int + (n_int - 1 - i)];
660 unsigned int n_dofs_1d)
const
662 unsigned int x, y, z;
667 (local_line % 4 == 0) ? 0 : (local_line % 4 == 1) ? fe_degree : dof;
669 (local_line % 4 == 2) ? 0 : (local_line % 4 == 3) ? fe_degree : dof;
670 z = (local_line / 4) * fe_degree;
674 x = ((local_line - 8) % 2) * fe_degree;
675 y = ((local_line - 8) / 2) * fe_degree;
679 return n_dofs_1d * n_dofs_1d * z + n_dofs_1d * y + x;
687 std::vector<types::global_dof_index> &dofs)
const
689 const std::vector<types::global_dof_index>
copy(dofs);
696 const unsigned int n_int = fe_degree - 1;
697 unsigned int offset = 4;
698 for (
unsigned int i = 0; i < n_int; ++i)
701 dofs[offset + i] =
copy[offset + 2 * n_int + i];
703 dofs[offset + n_int + i] =
copy[offset + 3 * n_int + i];
705 dofs[offset + 2 * n_int + i] =
copy[offset + i];
707 dofs[offset + 3 * n_int + i] =
copy[offset + n_int + i];
712 for (
unsigned int i = 0; i < n_int; ++i)
713 for (
unsigned int j = 0; j < n_int; ++j)
714 dofs[offset + i * n_int + j] =
copy[offset + j * n_int + i];
725 else if (subface == 2)
736 template <
unsigned int size>
737 __device__
inline unsigned int
745 template <
unsigned int size>
746 __device__
inline unsigned int
747 index3(
unsigned int i,
unsigned int j,
unsigned int k)
749 return i + size * j + size * size * k;
754 template <
unsigned int fe_degree,
755 unsigned int direction,
758 __device__
inline void
762 const unsigned int x_idx = threadIdx.x % (fe_degree + 1);
763 const unsigned int y_idx = threadIdx.y;
765 const unsigned int this_type =
768 const unsigned int interp_idx = (direction == 0) ? x_idx : y_idx;
773 const bool constrained_face =
781 const bool constrained_dof =
784 (y_idx == fe_degree))) ||
787 (x_idx == fe_degree)));
789 if (constrained_face && constrained_dof)
791 const bool type = constraint_mask & this_type;
795 for (
unsigned int i = 0; i <= fe_degree; ++i)
797 const unsigned int real_idx =
798 (direction == 0) ? index2<fe_degree + 1>(i, y_idx) :
799 index2<fe_degree + 1>(x_idx, i);
813 for (
unsigned int i = 0; i <= fe_degree; ++i)
815 const unsigned int real_idx =
816 (direction == 0) ? index2<fe_degree + 1>(i, y_idx) :
817 index2<fe_degree + 1>(x_idx, i);
823 fe_degree - interp_idx] :
835 if (constrained_face && constrained_dof)
836 values[index2<fe_degree + 1>(x_idx, y_idx)] = t;
843 template <
unsigned int fe_degree,
844 unsigned int direction,
847 __device__
inline void
851 const unsigned int x_idx = threadIdx.x % (fe_degree + 1);
852 const unsigned int y_idx = threadIdx.y;
853 const unsigned int z_idx = threadIdx.z;
855 const unsigned int this_type =
859 const unsigned int face1_type =
863 const unsigned int face2_type =
882 const unsigned int constrained_face =
883 constraint_mask & (face1 | face2 | edge);
885 const unsigned int interp_idx =
886 (direction == 0) ? x_idx : (direction == 1) ? y_idx : z_idx;
887 const unsigned int face1_idx =
888 (direction == 0) ? y_idx : (direction == 1) ? z_idx : x_idx;
889 const unsigned int face2_idx =
890 (direction == 0) ? z_idx : (direction == 1) ? x_idx : y_idx;
893 const bool on_face1 = (constraint_mask & face1_type) ?
895 (face1_idx == fe_degree);
896 const bool on_face2 = (constraint_mask & face2_type) ?
898 (face2_idx == fe_degree);
899 const bool constrained_dof =
900 (((constraint_mask & face1) && on_face1) ||
901 ((constraint_mask & face2) && on_face2) ||
902 ((constraint_mask & edge) && on_face1 && on_face2));
904 if (constrained_face && constrained_dof)
906 const bool type = constraint_mask & this_type;
909 for (
unsigned int i = 0; i <= fe_degree; ++i)
911 const unsigned int real_idx =
913 index3<fe_degree + 1>(i, y_idx, z_idx) :
915 index3<fe_degree + 1>(x_idx, i, z_idx) :
916 index3<fe_degree + 1>(x_idx, y_idx, i);
930 for (
unsigned int i = 0; i <= fe_degree; ++i)
932 const unsigned int real_idx =
934 index3<fe_degree + 1>(i, y_idx, z_idx) :
936 index3<fe_degree + 1>(x_idx, i, z_idx) :
937 index3<fe_degree + 1>(x_idx, y_idx, i);
943 fe_degree - interp_idx] :
956 if (constrained_face && constrained_dof)
957 values[index3<fe_degree + 1>(x_idx, y_idx, z_idx)] = t;
973 template <
int dim,
int fe_degree,
bool transpose,
typename Number>
979 internal::interpolate_boundary_2d<fe_degree, 0, transpose>(
982 internal::interpolate_boundary_2d<fe_degree, 1, transpose>(
988 internal::interpolate_boundary_3d<fe_degree, 0, transpose>(
991 internal::interpolate_boundary_3d<fe_degree, 1, transpose>(
994 internal::interpolate_boundary_3d<fe_degree, 2, transpose>(
typename DoFHandler< dim >::active_cell_iterator active_cell_iterator
std::vector< std::vector< std::pair< cell_iterator, unsigned int > > > line_to_cells
const unsigned int fe_degree
void transpose_subface_index(unsigned int &subface) const
const std::vector< unsigned int > & lexicographic_mapping
void rotate_subface_index(int times, unsigned int &subface_index) const
void setup_constraints(std::vector< types::global_dof_index > &dof_indices, const CellIterator &cell, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner, unsigned int &mask) const
void rotate_face(int times, unsigned int n_dofs_1d, std::vector< types::global_dof_index > &dofs) const
const DoFHandler< dim > & dof_handler
unsigned int line_dof_idx(int local_line, unsigned int dof, unsigned int n_dofs_1d) const
void setup_line_to_cell()
void transpose_face(std::vector< types::global_dof_index > &dofs) const
const unsigned int n_raw_lines
HangingNodes(unsigned int fe_degree, const DoFHandler< dim > &dof_handler, const std::vector< unsigned int > &lexicographic_mapping)
typename DoFHandler< dim >::cell_iterator cell_iterator
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInternalError()
#define AssertCuda(error_code)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
constexpr unsigned int constr_edge_yz
void interpolate_boundary_2d(const unsigned int constraint_mask, Number *values)
void interpolate_boundary_3d(const unsigned int constraint_mask, Number *values)
void setup_constraint_weigths(unsigned int fe_degree)
__constant__ double constraint_weights[(mf_max_elem_degree+1) *(mf_max_elem_degree+1)]
constexpr unsigned int constr_edge_zx
constexpr unsigned int constr_type_y
unsigned int index2(unsigned int i, unsigned int j)
constexpr unsigned int constr_edge_xy
constexpr unsigned int constr_face_z
constexpr unsigned int constr_face_y
constexpr unsigned int constr_type_z
constexpr unsigned int constr_type_x
constexpr unsigned int constr_face_x
unsigned int index3(unsigned int i, unsigned int j, unsigned int k)
void resolve_hanging_nodes(const unsigned int constraint_mask, Number *values)
constexpr unsigned int mf_max_elem_degree
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void copy(const T *begin, const T *end, U *dest)