ExaDG
Loading...
Searching...
No Matches
calculate_characteristic_element_length.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 INCLUDE_FUNCTIONALITIES_CALCULATE_CHARACTERISTIC_ELEMENT_LENGTH_H_
23#define INCLUDE_FUNCTIONALITIES_CALCULATE_CHARACTERISTIC_ELEMENT_LENGTH_H_
24
25// deal.II
26#include <deal.II/base/mpi.h>
27#include <deal.II/fe/mapping.h>
28#include <deal.II/grid/tria.h>
29
30// ExaDG
31#include <exadg/grid/grid_data.h>
32
33namespace ExaDG
34{
40template<int dim>
41inline double
42calculate_minimum_vertex_distance(dealii::Triangulation<dim> const & triangulation,
43 dealii::Mapping<dim> const & mapping,
44 MPI_Comm const & mpi_comm)
45{
46 double min_cell_diameter = std::numeric_limits<double>::max();
47
48 for(auto const & cell : triangulation.active_cell_iterators())
49 {
50 if(cell->is_locally_owned())
51 {
52 auto const vertices = mapping.get_vertices(cell);
53 for(unsigned int i = 0; i < vertices.size(); ++i)
54 for(unsigned int j = i + 1; j < vertices.size(); ++j)
55 min_cell_diameter = std::min(min_cell_diameter, vertices[i].distance(vertices[j]));
56 }
57 }
58
59 return dealii::Utilities::MPI::min(min_cell_diameter, mpi_comm);
60}
61
66template<typename Number>
67inline Number
68calculate_characteristic_element_length(Number const element_volume,
69 unsigned int const dim,
70 ElementType const & element_type)
71{
72 if(element_type == ElementType::Hypercube)
73 {
74 // calculate the length h of a dim-dimensional cube that has the volume "element_volume".
75 return std::exp(std::log(element_volume) / (double)dim);
76 }
77 else if(element_type == ElementType::Simplex)
78 {
79 // calculate the length h of an equilateral triangle/tetrahedron that has the volume
80 // "element_volume" (V = h^dim / factor).
81 double factor = 1.0;
82
83 if(dim == 2)
84 factor = 4.0 / std::sqrt(3.0);
85 else if(dim == 3)
86 factor = 6.0 * std::sqrt(2.0);
87 else
88 AssertThrow(false, dealii::ExcMessage("Only dim = 2, 3 implemented."));
89
90 return std::exp(std::log(element_volume * factor) / (double)dim);
91 }
92 else
93 {
94 AssertThrow(false,
95 dealii::ExcMessage("This function is not implemented for the given ElementType."));
96 return std::exp(std::log(element_volume) / (double)dim);
97 }
98}
99
106template<typename Number>
107inline Number
108calculate_high_order_element_length(Number const element_length,
109 unsigned int const fe_degree,
110 bool const is_dg)
111{
112 unsigned int n_nodes_1d = fe_degree;
113
114 if(is_dg)
115 n_nodes_1d += 1;
116
117 return element_length / ((double)n_nodes_1d);
118}
119
120} // namespace ExaDG
121
122#endif /* INCLUDE_FUNCTIONALITIES_CALCULATE_CHARACTERISTIC_ELEMENT_LENGTH_H_ */
Definition driver.cpp:33
Number calculate_high_order_element_length(Number const element_length, unsigned int const fe_degree, bool const is_dg)
Definition calculate_characteristic_element_length.h:108
double calculate_minimum_vertex_distance(dealii::Triangulation< dim > const &triangulation, dealii::Mapping< dim > const &mapping, MPI_Comm const &mpi_comm)
Definition calculate_characteristic_element_length.h:42
Number calculate_characteristic_element_length(Number const element_volume, unsigned int const dim, ElementType const &element_type)
Definition calculate_characteristic_element_length.h:68