/
copy_serial_tria_06.cc
121 lines (97 loc) · 3.73 KB
/
copy_serial_tria_06.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
109
110
111
112
113
114
115
116
117
118
119
120
121
// ---------------------------------------------------------------------
//
// Copyright (C) 2019 - 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.
//
// ---------------------------------------------------------------------
// Create a serial triangulation, copy it, and create a partitioner
// by matrixfree.
#include <deal.II/base/mpi.h>
#include <deal.II/distributed/fully_distributed_tria.h>
#include <deal.II/distributed/shared_tria.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_out.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_description.h>
#include <deal.II/matrix_free/matrix_free.h>
#include "./tests.h"
using namespace dealii;
template <int dim>
void
test(int n_refinements, MPI_Comm comm)
{
// create serial triangulation
Triangulation<dim> basetria;
GridGenerator::hyper_cube(basetria);
basetria.refine_global(n_refinements);
GridTools::partition_triangulation_zorder(
Utilities::MPI::n_mpi_processes(comm), basetria);
// create instance of pft
parallel::fullydistributed::Triangulation<dim> tria_pft(comm);
// extract relevant information form serial triangulation
auto construction_data =
TriangulationDescription::Utilities::create_description_from_triangulation(
basetria, comm);
// actually create triangulation
tria_pft.create_triangulation(construction_data);
unsigned int degree = 2;
// test triangulation
FE_DGQ<dim> fe(degree);
DoFHandler<dim> dof_handler(tria_pft);
dof_handler.distribute_dofs(fe);
// test 1: print the indices of the dofs of each active cell
{
std::vector<types::global_dof_index> dofs(Utilities::pow(degree + 1, dim));
for (auto cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned() || cell->is_ghost())
{
cell->get_dof_indices(dofs);
deallog << std::setw(4) << cell->subdomain_id() << " : ";
for (const auto dof : dofs)
deallog << std::setw(4) << dof << " ";
deallog << std::endl;
}
}
// test 2: let matrixfree setup all partitioners
{
typename dealii::MatrixFree<dim, double>::AdditionalData additional_data;
additional_data.mapping_update_flags =
update_gradients | update_JxW_values | update_quadrature_points;
additional_data.mapping_update_flags_inner_faces =
update_gradients | update_JxW_values;
additional_data.mapping_update_flags_boundary_faces =
update_gradients | update_JxW_values | update_quadrature_points;
additional_data.hold_all_faces_to_owned_cells = true;
additional_data.mapping_update_flags_faces_by_cells =
update_gradients | update_JxW_values | update_quadrature_points;
MappingQ<dim> mapping(1);
QGauss<1> quad(degree + 1);
AffineConstraints<double> constraint;
MatrixFree<dim, double> matrix_free;
matrix_free.reinit(mapping, dof_handler, constraint, quad, additional_data);
}
}
int
main(int argc, char *argv[])
{
Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
MPILogInitAll all;
const int n_refinements = 3;
const MPI_Comm comm = MPI_COMM_WORLD;
{
deallog.push("2d");
test<2>(n_refinements, comm);
deallog.pop();
}
}