41 template <
int dim,
int spacedim>
52 if (
const auto *p_tria =
dynamic_cast<
64 const std::vector<Point<spacedim>> &vertices = tria.
get_vertices();
65 std::vector<bool> boundary_vertices(vertices.size(),
false);
71 for (; cell != endc; ++cell)
72 for (
const unsigned int face : cell->face_indices())
73 if (cell->face(face)->at_boundary())
74 for (
unsigned int i = 0; i < cell->face(face)->n_vertices(); ++i)
75 boundary_vertices[cell->face(face)->vertex_index(i)] =
true;
79 double max_distance_sqr = 0;
80 std::vector<bool>::const_iterator pi = boundary_vertices.begin();
81 const unsigned int N = boundary_vertices.size();
82 for (
unsigned int i = 0; i < N; ++i, ++pi)
84 std::vector<bool>::const_iterator pj = pi + 1;
85 for (
unsigned int j = i + 1; j < N; ++j, ++pj)
86 if ((*pi ==
true) && (*pj ==
true) &&
87 ((vertices[i] - vertices[j]).norm_square() > max_distance_sqr))
88 max_distance_sqr = (vertices[i] - vertices[j]).norm_square();
96 template <
int dim,
int spacedim>
105 reference_cell.template get_default_linear_mapping<dim, spacedim>());
110 template <
int dim,
int spacedim>
116 unsigned int mapping_degree = 1;
118 mapping_degree = p->get_degree();
119 else if (
const auto *p =
121 mapping_degree = p->get_degree();
128 reference_cell.template get_gauss_type_quadrature<dim>(mapping_degree +
130 const unsigned int n_q_points = quadrature_formula.
size();
141 double local_volume = 0;
144 for (
const auto &cell :
triangulation.active_cell_iterators())
145 if (cell->is_locally_owned())
148 for (
unsigned int q = 0; q < n_q_points; ++q)
149 local_volume += fe_values.
JxW(q);
152 const double global_volume =
155 return global_volume;
160 template <
int dim,
int spacedim>
161 std::pair<unsigned int, double>
165 double max_ratio = 1;
166 unsigned int index = 0;
168 for (
unsigned int i = 0; i < dim; ++i)
169 for (
unsigned int j = i + 1; j < dim; ++j)
171 unsigned int ax = i % dim;
172 unsigned int next_ax = j % dim;
175 cell->extent_in_direction(ax) / cell->extent_in_direction(next_ax);
177 if (ratio > max_ratio)
182 else if (1.0 / ratio > max_ratio)
184 max_ratio = 1.0 / ratio;
188 return std::make_pair(index, max_ratio);
210 struct TransformR2UAffine
225 [1] = {{-1.000000}, {1.000000}};
229 {1.000000, 0.000000};
240 [2] = {{-0.500000, -0.500000},
241 {0.500000, -0.500000},
242 {-0.500000, 0.500000},
243 {0.500000, 0.500000}};
253 {0.750000, 0.250000, 0.250000, -0.250000};
259 {-0.250000, -0.250000, -0.250000},
260 {0.250000, -0.250000, -0.250000},
261 {-0.250000, 0.250000, -0.250000},
262 {0.250000, 0.250000, -0.250000},
263 {-0.250000, -0.250000, 0.250000},
264 {0.250000, -0.250000, 0.250000},
265 {-0.250000, 0.250000, 0.250000},
266 {0.250000, 0.250000, 0.250000}
285 template <
int dim,
int spacedim>
294 for (
unsigned int d = 0; d < spacedim; ++d)
295 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
296 for (
unsigned int e = 0; e < dim; ++e)
297 A[d][e] += vertices[v][d] * TransformR2UAffine<dim>::KA[v][e];
301 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
302 b += vertices[v] * TransformR2UAffine<dim>::Kb[v];
304 return std::make_pair(A, b);
321 for (
const auto &cell :
triangulation.active_cell_iterators())
323 if (cell->is_locally_owned())
325 double aspect_ratio_cell = 0.0;
330 for (
unsigned int q = 0; q < quadrature.
size(); ++q)
341 aspect_ratio_cell = std::numeric_limits<double>::infinity();
346 for (
unsigned int i = 0; i < dim; ++i)
347 for (
unsigned int j = 0; j < dim; ++j)
348 J(i, j) = jacobian[i][j];
354 const double ar = max_sv / min_sv;
361 aspect_ratio_cell =
std::max(aspect_ratio_cell, ar);
366 aspect_ratio_vector(cell->active_cell_index()) = aspect_ratio_cell;
370 return aspect_ratio_vector;
391 template <
int dim,
int spacedim>
397 const auto predicate = [](
const iterator &) {
return true; };
400 tria, std::function<
bool(
const iterator &)>(predicate));
405 template <
int dim,
int spacedim>
410 double min_diameter = std::numeric_limits<double>::max();
411 for (
const auto &cell :
triangulation.active_cell_iterators())
412 if (!cell->is_artificial())
413 min_diameter =
std::min(min_diameter, cell->diameter(mapping));
415 const double global_min_diameter =
417 return global_min_diameter;
422 template <
int dim,
int spacedim>
427 double max_diameter = 0.;
428 for (
const auto &cell :
triangulation.active_cell_iterators())
429 if (!cell->is_artificial())
430 max_diameter =
std::max(max_diameter, cell->diameter(mapping));
432 const double global_max_diameter =
434 return global_max_diameter;
440#include "grid_tools_geometry.inst"
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int q_point) const
double JxW(const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
number singular_value(const size_type i) const
Abstract base class for mapping classes.
unsigned int size() const
const std::vector< Point< spacedim > > & get_vertices() const
cell_iterator end() const
unsigned int n_cells() const
active_cell_iterator begin_active(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
T sum(const T &t, const MPI_Comm mpi_communicator)
T max(const T &t, const MPI_Comm mpi_communicator)
T min(const T &t, const MPI_Comm mpi_communicator)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)