/
hp_unify_dof_indices_02.cc
108 lines (85 loc) · 3.09 KB
/
hp_unify_dof_indices_02.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
// ---------------------------------------------------------------------
//
// Copyright (C) 2017 - 2020 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------
// have a 2x2 coarse mesh (or 2x2x1) and verify DoF indices in the hp-
// case with an FECollection that contains multiple copies of the same
// FE_Q(2) element. the hp-code will unify DoF indices on boundaries
// between all subdomains.
//
// in this testcase, three FE objects are distributed on four cells,
// from which two, which are diagonally opposite, get the same
// active_fe_index assigned.
#include <deal.II/base/tensor.h>
#include <deal.II/base/utilities.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/intergrid_map.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/hp/fe_collection.h>
#include <numeric>
#include "../tests.h"
#include "hp_unify_dof_indices.h"
template <int dim>
void
test()
{
parallel::distributed::Triangulation<dim> triangulation(
MPI_COMM_WORLD, Triangulation<dim>::limit_level_difference_at_vertices);
std::vector<unsigned int> reps(dim, 1U);
reps[0] = 2;
reps[1] = 2;
Point<dim> top_right;
for (unsigned int d = 0; d < dim; ++d)
top_right[d] = (d == 0 ? 2 : 1);
GridGenerator::subdivided_hyper_rectangle(triangulation,
reps,
Point<dim>(),
top_right);
Assert(triangulation.n_global_active_cells() == 4, ExcInternalError());
Assert(triangulation.n_active_cells() == 4, ExcInternalError());
hp::FECollection<dim> fe;
fe.push_back(FE_Q<dim>(2));
fe.push_back(FE_Q<dim>(2));
fe.push_back(FE_Q<dim>(2));
DoFHandler<dim> dof_handler(triangulation);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
{
if (cell->id().to_string() == "0_0:")
cell->set_active_fe_index(0);
if (cell->id().to_string() == "1_0:")
cell->set_active_fe_index(1);
if (cell->id().to_string() == "2_0:")
cell->set_active_fe_index(2);
if (cell->id().to_string() == "3_0:")
cell->set_active_fe_index(0);
}
dof_handler.distribute_dofs(fe);
log_dof_diagnostics(dof_handler);
}
int
main(int argc, char *argv[])
{
Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
MPILogInitAll log;
deallog.push("2d");
test<2>();
deallog.pop();
deallog.push("3d");
test<3>();
deallog.pop();
}