/
fe_interface_values_02.cc
140 lines (111 loc) · 4.15 KB
/
fe_interface_values_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
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
// ---------------------------------------------------------------------
//
// Copyright (C) 2018 - 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.
//
// ---------------------------------------------------------------------
// evaluate jump(), average(), shape_value() of FEInterfaceValues
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_interface_values.h>
#include <deal.II/fe/mapping_q.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_refinement.h>
#include <deal.II/grid/tria.h>
#include <fstream>
#include <iostream>
#include "../tests.h"
template <int dim>
void
make_2_cells(Triangulation<dim> &tria);
template <>
void make_2_cells<2>(Triangulation<2> &tria)
{
const unsigned int dim = 2;
std::vector<unsigned int> repetitions = {2, 1};
Point<dim> p1;
Point<dim> p2(2.0, 1.0);
GridGenerator::subdivided_hyper_rectangle(tria, repetitions, p1, p2);
}
template <>
void make_2_cells<3>(Triangulation<3> &tria)
{
const unsigned int dim = 3;
std::vector<unsigned int> repetitions = {2, 1, 1};
Point<dim> p1;
Point<dim> p2(2.0, 1.0, 1.0);
GridGenerator::subdivided_hyper_rectangle(tria, repetitions, p1, p2);
}
template <int dim>
void
test(unsigned int fe_degree)
{
Triangulation<dim> tria;
make_2_cells(tria);
DoFHandler<dim> dofh(tria);
FE_DGQ<dim> fe(fe_degree);
deallog << fe.get_name() << std::endl;
dofh.distribute_dofs(fe);
MappingQ<dim> mapping(1);
UpdateFlags update_flags = update_values | update_gradients |
update_quadrature_points | update_JxW_values;
FEInterfaceValues<dim> fiv(mapping,
fe,
QGauss<dim - 1>(fe.degree + 1),
update_flags);
auto cell = dofh.begin();
for (const unsigned int f : GeometryInfo<dim>::face_indices())
if (!cell->at_boundary(f))
{
fiv.reinit(cell,
f,
numbers::invalid_unsigned_int,
cell->neighbor(f),
cell->neighbor_of_neighbor(f),
numbers::invalid_unsigned_int);
const unsigned int n_dofs = fiv.n_current_interface_dofs();
Vector<double> cell_vector(n_dofs);
const auto &q_points = fiv.get_quadrature_points();
cell_vector = 0.0;
for (unsigned int qpoint = 0; qpoint < q_points.size(); ++qpoint)
for (unsigned int i = 0; i < n_dofs; ++i)
cell_vector(i) +=
fiv.shape_value(true, i, qpoint) * fiv.get_JxW_values()[qpoint];
deallog << "shape_value(true): " << cell_vector << std::endl;
cell_vector = 0.0;
for (unsigned int qpoint = 0; qpoint < q_points.size(); ++qpoint)
for (unsigned int i = 0; i < n_dofs; ++i)
cell_vector(i) +=
fiv.shape_value(false, i, qpoint) * fiv.get_JxW_values()[qpoint];
deallog << "shape_value(false): " << cell_vector << std::endl;
cell_vector = 0.0;
for (unsigned int qpoint = 0; qpoint < q_points.size(); ++qpoint)
for (unsigned int i = 0; i < n_dofs; ++i)
cell_vector(i) +=
fiv.jump(i, qpoint) * fiv.get_JxW_values()[qpoint];
deallog << "jump(): " << cell_vector << std::endl;
cell_vector = 0.0;
for (unsigned int qpoint = 0; qpoint < q_points.size(); ++qpoint)
for (unsigned int i = 0; i < n_dofs; ++i)
cell_vector(i) +=
fiv.average(i, qpoint) * fiv.get_JxW_values()[qpoint];
deallog << "average(): " << cell_vector << std::endl;
}
}
int
main()
{
initlog();
test<2>(0);
test<2>(1);
test<3>(0);
test<3>(1);
}