ExaDG
Loading...
Searching...
No Matches
periodic_box.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 APPLICATIONS_GRID_TOOLS_PERIODIC_BOX_H_
23#define APPLICATIONS_GRID_TOOLS_PERIODIC_BOX_H_
24
25// deal.II
26#include <deal.II/distributed/tria.h>
27#include <deal.II/grid/grid_generator.h>
28#include <deal.II/grid/grid_tools.h>
29
30// ExaDG
31#include <exadg/grid/deformed_cube_manifold.h>
32
33namespace ExaDG
34{
35template<int dim>
36void
37create_periodic_box(dealii::Triangulation<dim> & triangulation,
38 unsigned int const n_refine_space,
39 std::vector<dealii::GridTools::PeriodicFacePair<
40 typename dealii::Triangulation<dim>::cell_iterator>> & periodic_faces,
41 unsigned int const n_subdivisions,
42 double const left,
43 double const right,
44 bool const curvilinear_mesh = false,
45 double const deformation = 0.1)
46{
47 dealii::GridGenerator::subdivided_hyper_cube(triangulation, n_subdivisions, left, right);
48
49 if(curvilinear_mesh)
50 {
51 unsigned int const frequency = 2;
52 apply_deformed_cube_manifold(triangulation, left, right, deformation, frequency);
53 }
54
55 for(auto const & cell : triangulation.cell_iterators())
56 {
57 for(unsigned int const face_number : cell->face_indices())
58 {
59 // x-direction
60 if((std::fabs(cell->face(face_number)->center()(0) - left) < 1e-12))
61 cell->face(face_number)->set_all_boundary_ids(0);
62 else if((std::fabs(cell->face(face_number)->center()(0) - right) < 1e-12))
63 cell->face(face_number)->set_all_boundary_ids(1);
64 // y-direction
65 else if((std::fabs(cell->face(face_number)->center()(1) - left) < 1e-12))
66 cell->face(face_number)->set_all_boundary_ids(2);
67 else if((std::fabs(cell->face(face_number)->center()(1) - right) < 1e-12))
68 cell->face(face_number)->set_all_boundary_ids(3);
69 // z-direction
70 else if(dim == 3 and (std::fabs(cell->face(face_number)->center()(2) - left) < 1e-12))
71 cell->face(face_number)->set_all_boundary_ids(4);
72 else if(dim == 3 and (std::fabs(cell->face(face_number)->center()(2) - right) < 1e-12))
73 cell->face(face_number)->set_all_boundary_ids(5);
74 }
75 }
76
77 dealii::GridTools::collect_periodic_faces(triangulation, 0, 1, 0 /*x-direction*/, periodic_faces);
78 dealii::GridTools::collect_periodic_faces(triangulation, 2, 3, 1 /*y-direction*/, periodic_faces);
79 if(dim == 3)
80 dealii::GridTools::collect_periodic_faces(
81 triangulation, 4, 5, 2 /*z-direction*/, periodic_faces);
82
83 triangulation.add_periodicity(periodic_faces);
84
85 // perform global refinements
86 triangulation.refine_global(n_refine_space);
87}
88
89} // namespace ExaDG
90
91#endif /* APPLICATIONS_GRID_TOOLS_PERIODIC_BOX_H_ */
Definition driver.cpp:33