Fluid structure interaction suite
compute_linear_transformation.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2022 by Luca Heltai
4 //
5 // This file is part of the FSI-suite platform, based on the deal.II library.
6 //
7 // The FSI-suite platform is free software; you can use it, redistribute it,
8 // and/or modify it under the terms of the GNU Lesser General Public License as
9 // published by the Free Software Foundation; either version 3.0 of the License,
10 // or (at your option) any later version. The full text of the license can be
11 // found in the file LICENSE at the top level of the FSI-suite platform
12 // distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef compute_linear_transformation_h
17 #define compute_linear_transformation_h
18 
20 
21 #include <deal.II/dofs/dof_tools.h>
22 
23 #include <deal.II/fe/fe_nothing.h>
24 #include <deal.II/fe/fe_values.h>
25 
27 
28 
41 template <int dim, int spacedim, int N>
42 dealii::Quadrature<spacedim>
44  const dealii::Quadrature<dim> &quadrature,
45  const std::array<dealii::Point<spacedim>, N> &vertices);
46 
47 
48 
49 // Template implementation
50 #ifndef DOXYGEN
51 template <int dim, int spacedim, int N>
52 dealii::Quadrature<spacedim>
54  const dealii::Quadrature<dim> &quadrature,
55  const std::array<dealii::Point<spacedim>, N> &vertices)
56 {
58  const auto CellType = dealii::ReferenceCell::n_vertices_to_type(
59  dim, N); // understand the kind of reference cell from vertices
60 
61  dealii::Triangulation<dim, spacedim> tria;
62  dealii::GridGenerator::reference_cell(
63  tria, CellType); // store reference cell stored in tria
64  dealii::FE_Nothing<dim, spacedim> dummy_fe(CellType);
65  dealii::DoFHandler<dim, spacedim> dh(tria);
66  dh.distribute_dofs(dummy_fe);
67  dealii::FEValues<dim, spacedim> fe_values(dummy_fe,
68  quadrature,
71  const auto &cell = dh.begin_active();
72  for (unsigned int i = 0; i < N; ++i)
73  cell->vertex(i) = vertices[i]; // the vertices of this real cell
74  fe_values.reinit(cell);
75 
76  return dealii::Quadrature<spacedim>(
77  fe_values.get_quadrature_points(),
78  fe_values.get_JxW_values()); // points and weigths in the real space
79 }
80 #endif
81 
82 #endif
Quadrature< spacedim > compute_linear_transformation(const Quadrature< dim > &quadrature, const std::array< Point< spacedim >, N > &vertices)
Given a dim-dimensional quadrature formula to integrate over vertices, returns a spacedim-dimensional...
Point< 3 > vertices[4]
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
update_JxW_values
update_quadrature_points