-
Notifications
You must be signed in to change notification settings - Fork 737
/
point_evaluation_13.cc
154 lines (116 loc) · 4.38 KB
/
point_evaluation_13.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
141
142
143
144
145
146
147
148
149
150
151
152
153
154
// ---------------------------------------------------------------------
//
// Copyright (C) 2021 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.
//
// ---------------------------------------------------------------------
// Check FEPointEvaluation for vector-valued FE_Q and MappingCartesian by
// comparing to the output of FEValues with the same settings.
// This is a modified version of point_evaluation_03.cc.
#include <deal.II/base/function_lib.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping_cartesian.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria.h>
#include <deal.II/lac/vector.h>
#include <deal.II/matrix_free/fe_point_evaluation.h>
#include <deal.II/numerics/vector_tools.h>
#include <iostream>
#include "../tests.h"
template <int dim>
class MyFunction : public Function<dim>
{
public:
MyFunction()
: Function<dim>(dim)
{}
double
value(const Point<dim> &p, const unsigned int component) const override
{
return component + p[component];
}
};
template <int dim>
void
test()
{
using namespace dealii;
Triangulation<dim> tria;
GridGenerator::subdivided_hyper_cube(tria, 2, 0, 1);
MappingCartesian<dim> mapping;
deallog << "Cartesian linear mapping" << std::endl;
std::vector<Point<dim>> unit_points;
for (unsigned int i = 0; i < 13; ++i)
{
Point<dim> p;
for (unsigned int d = 0; d < dim; ++d)
p[d] = static_cast<double>(i) / 17. + 0.015625 * d;
unit_points.push_back(p);
}
FESystem<dim> fe(FE_Q<dim>(2), dim);
FEValues<dim> fe_values(mapping,
fe,
Quadrature<dim>(unit_points),
update_values | update_gradients);
DoFHandler<dim> dof_handler(tria);
dof_handler.distribute_dofs(fe);
Vector<double> vector(dof_handler.n_dofs());
FEPointEvaluation<dim, dim> evaluator(mapping,
fe,
update_values | update_gradients);
VectorTools::interpolate(mapping, dof_handler, MyFunction<dim>(), vector);
std::vector<double> solution_values(fe.dofs_per_cell);
std::vector<Tensor<1, dim>> function_values(unit_points.size());
std::vector<Tensor<2, dim>> function_gradients(unit_points.size());
FEValuesExtractors::Vector extractor(0);
for (const auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
fe_values[extractor].get_function_values(vector, function_values);
fe_values[extractor].get_function_gradients(vector, function_gradients);
cell->get_dof_values(vector,
solution_values.begin(),
solution_values.end());
evaluator.reinit(cell, unit_points);
evaluator.evaluate(solution_values,
EvaluationFlags::values | EvaluationFlags::gradients);
deallog << "Cell with center " << cell->center(true) << std::endl;
for (unsigned int i = 0; i < function_values.size(); ++i)
deallog << mapping.transform_unit_to_real_cell(cell, unit_points[i])
<< ": " << evaluator.get_value(i) << " error value "
<< (function_values[i] - evaluator.get_value(i)).norm()
<< " error grad "
<< (evaluator.get_gradient(i) - function_gradients[i]).norm()
<< std::endl;
deallog << std::endl;
for (unsigned int i = 0; i < unit_points.size(); ++i)
{
evaluator.submit_value(evaluator.get_value(i), i);
evaluator.submit_gradient(evaluator.get_gradient(i), i);
}
evaluator.integrate(solution_values,
EvaluationFlags::values | EvaluationFlags::gradients);
for (const auto i : solution_values)
deallog << i << " ";
deallog << std::endl;
}
}
int
main()
{
initlog();
deallog << std::setprecision(10);
test<2>();
test<3>();
}