ExaDG
Loading...
Searching...
No Matches
deformed_cube_manifold.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_GRID_DEFORMED_CUBE_MANIFOLD_H_
23#define EXADG_GRID_DEFORMED_CUBE_MANIFOLD_H_
24
25// deal.II
26#include <deal.II/grid/manifold_lib.h>
27
28namespace ExaDG
29{
35template<int dim>
36class DeformedCubeManifold : public dealii::ChartManifold<dim, dim, dim>
37{
38public:
39 DeformedCubeManifold(double const left,
40 double const right,
41 double const deformation,
42 unsigned int const frequency = 1)
43 : left(left), right(right), deformation(deformation), frequency(frequency)
44 {
45 }
46
47 dealii::Point<dim>
48 push_forward(dealii::Point<dim> const & chart_point) const override
49 {
50 double sinval = deformation;
51 for(unsigned int d = 0; d < dim; ++d)
52 sinval *=
53 std::sin(frequency * dealii::numbers::PI * (chart_point(d) - left) / (right - left));
54 dealii::Point<dim> space_point;
55 for(unsigned int d = 0; d < dim; ++d)
56 space_point(d) = chart_point(d) + sinval;
57
58 return space_point;
59 }
60
61 dealii::Point<dim>
62 pull_back(dealii::Point<dim> const & space_point) const override
63 {
64 dealii::Point<dim> x = space_point;
65 dealii::Point<dim> one;
66 for(unsigned int d = 0; d < dim; ++d)
67 one(d) = 1.;
68
69 // Newton iteration to solve the nonlinear equation given by the point
70 dealii::Tensor<1, dim> sinvals, residual;
71
72 for(unsigned int d = 0; d < dim; ++d)
73 sinvals[d] = std::sin(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
74 double sinval = deformation;
75 for(unsigned int d = 0; d < dim; ++d)
76 sinval *= sinvals[d];
77 residual = space_point - (x + sinval * one);
78
79 unsigned int its = 0;
80 while(residual.norm() > 1e-12 and its < 100)
81 {
82 dealii::Tensor<2, dim> jacobian;
83 for(unsigned int d = 0; d < dim; ++d)
84 jacobian[d][d] = 1.;
85 for(unsigned int d = 0; d < dim; ++d)
86 {
87 double sinval_der =
88 deformation * frequency / (right - left) * dealii::numbers::PI *
89 std::cos(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
90 for(unsigned int e = 0; e < dim; ++e)
91 if(e != d)
92 sinval_der *= sinvals[e];
93 for(unsigned int e = 0; e < dim; ++e)
94 jacobian[e][d] += sinval_der;
95 }
96
97 x += invert(jacobian) * residual;
98
99 for(unsigned int d = 0; d < dim; ++d)
100 sinvals[d] = std::sin(frequency * dealii::numbers::PI * (x(d) - left) / (right - left));
101 sinval = deformation;
102 for(unsigned int d = 0; d < dim; ++d)
103 sinval *= sinvals[d];
104 residual = space_point - (x + sinval * one);
105
106 ++its;
107 }
108
109 AssertThrow(residual.norm() < 1e-12, dealii::ExcMessage("Newton for point did not converge."));
110
111 return x;
112 }
113
114 std::unique_ptr<dealii::Manifold<dim>>
115 clone() const override
116 {
117 return std::make_unique<DeformedCubeManifold<dim>>(left, right, deformation, frequency);
118 }
119
120private:
121 double const left;
122 double const right;
123 double const deformation;
124 unsigned int const frequency;
125};
126
127template<int dim>
128void
129apply_deformed_cube_manifold(dealii::Triangulation<dim> & triangulation,
130 double const left,
131 double const right,
132 double const deformation,
133 unsigned int const frequency)
134{
135 DeformedCubeManifold<dim> manifold(left, right, deformation, frequency);
136 triangulation.set_all_manifold_ids(1);
137 triangulation.set_manifold(1, manifold);
138
139 // the vertices need to be placed correctly according to the manifold description such that mesh
140 // refinements and high-order mappings are done correctly (which invoke pull-back and push-forward
141 // operations)
142 std::vector<bool> vertex_touched(triangulation.n_vertices(), false);
143
144 for(auto const & cell : triangulation.cell_iterators())
145 {
146 for(unsigned int const v : cell->vertex_indices())
147 {
148 if(vertex_touched[cell->vertex_index(v)] == false)
149 {
150 dealii::Point<dim> & vertex = cell->vertex(v);
151 dealii::Point<dim> new_point = manifold.push_forward(vertex);
152 vertex = new_point;
153 vertex_touched[cell->vertex_index(v)] = true;
154 }
155 }
156 }
157}
158
159} // namespace ExaDG
160
161#endif /* EXADG_GRID_DEFORMED_CUBE_MANIFOLD_H_ */
Definition deformed_cube_manifold.h:37
Definition driver.cpp:33