17#ifndef dealii_matrix_free_mapping_data_on_the_fly_h
18#define dealii_matrix_free_mapping_data_on_the_fly_h
41 namespace MatrixFreeFunctions
59 template <
int dim,
typename Number>
90 reinit(typename ::Triangulation<dim>::cell_iterator cell);
102 typename ::Triangulation<dim>::cell_iterator
111 const ::FEValues<dim> &
172 template <
int dim,
typename Number>
183 , quadrature_1d(quadrature)
209 template <
int dim,
typename Number>
220 template <
int dim,
typename Number>
223 typename ::Triangulation<dim>::cell_iterator cell)
225 if (present_cell == cell)
228 fe_values->reinit(present_cell);
229 for (
unsigned int q = 0; q < fe_values->get_quadrature().size(); ++q)
232 mapping_info_storage.JxW_values[q] = fe_values->JxW(q);
237 for (
unsigned int d = 0; d < dim; ++d)
238 for (
unsigned int e = 0; e < dim; ++e)
239 mapping_info_storage.jacobians[0][q][d][e] = jac[d][e];
242 for (
unsigned int d = 0; d < dim; ++d)
243 mapping_info_storage.quadrature_points[q][d] =
244 fe_values->quadrature_point(q)[d];
247 for (
unsigned int d = 0; d < dim; ++d)
248 mapping_info_storage.normal_vectors[q][d] =
249 fe_values->normal_vector(q)[d];
250 mapping_info_storage.normals_times_jacobians[0][q] =
251 mapping_info_storage.normal_vectors[q] *
252 mapping_info_storage.jacobians[0][q];
259 template <
int dim,
typename Number>
263 return present_cell !=
264 typename ::Triangulation<dim>::cell_iterator();
269 template <
int dim,
typename Number>
270 inline typename ::Triangulation<dim>::cell_iterator
273 return fe_values->get_cell();
278 template <
int dim,
typename Number>
279 inline const ::FEValues<dim> &
287 template <
int dim,
typename Number>
291 return mapping_info_storage;
296 template <
int dim,
typename Number>
300 return mapping_info_storage;
305 template <
int dim,
typename Number>
309 return quadrature_1d;
void resize(const size_type new_size)
Abstract base class for mapping classes.
const MappingInfoStorage< dim, dim, Number > & get_data_storage() const
const Quadrature< 1 > & get_quadrature() const
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
bool is_initialized() const
const ::FEValues< dim > & get_fe_values() const
MappingDataOnTheFly(const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
FE_Nothing< dim > fe_dummy
MappingInfoStorage< dim, dim, Number > & get_data_storage()
MappingDataOnTheFly()=default
typename::Triangulation< dim >::cell_iterator present_cell
const Quadrature< 1 > quadrature_1d
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
std::unique_ptr<::FEValues< dim > > fe_values
MappingInfoStorage< dim, dim, Number > mapping_info_storage
typename::Triangulation< dim >::cell_iterator get_cell() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
std::vector< QuadratureDescriptor > descriptor
AlignedVector< Tensor< 1, spacedim, Number > > normal_vectors
std::array< AlignedVector< Tensor< 1, spacedim, Number > >, 2 > normals_times_jacobians
AlignedVector< Number > JxW_values
AlignedVector< unsigned int > quadrature_point_offsets
AlignedVector< Point< spacedim, Number > > quadrature_points
AlignedVector< unsigned int > data_index_offsets
std::array< AlignedVector< Tensor< 2, spacedim, Number > >, 2 > jacobians
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)