17 #ifndef dealii__matrix_free_mapping_data_on_the_fly_h 18 #define dealii__matrix_free_mapping_data_on_the_fly_h 21 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/subscriptor.h> 24 #include <deal.II/base/vectorization.h> 25 #include <deal.II/base/aligned_vector.h> 26 #include <deal.II/matrix_free/shape_info.h> 27 #include <deal.II/matrix_free/mapping_info.h> 28 #include <deal.II/fe/fe_values.h> 29 #include <deal.II/fe/fe_nothing.h> 32 DEAL_II_NAMESPACE_OPEN
37 namespace MatrixFreeFunctions
57 template <
int dim,
typename Number=
double>
82 void reinit(typename ::Triangulation<dim>::cell_iterator cell);
93 typename ::Triangulation<dim>::cell_iterator
get_cell ()
const;
188 template <
int dim,
typename Number>
194 fe_values(mapping, fe_dummy,
Quadrature<dim>(quadrature),
195 internal::MatrixFreeFunctions::
MappingInfo<dim,Number>::compute_update_flags(update_flags)),
196 quadrature_1d(quadrature),
197 inverse_jacobians(fe_values.get_quadrature().size()),
198 jxw_values(fe_values.get_quadrature().size()),
199 quadrature_points(fe_values.get_quadrature().size()),
200 normal_vectors(fe_values.get_quadrature().size())
208 template <
int dim,
typename Number>
213 fe_values(fe_dummy,
Quadrature<dim>(quadrature),
214 internal::MatrixFreeFunctions::
MappingInfo<dim,Number>::compute_update_flags(update_flags)),
215 quadrature_1d(quadrature),
216 inverse_jacobians(fe_values.get_quadrature().size()),
217 jxw_values(fe_values.get_quadrature().size()),
218 quadrature_points(fe_values.get_quadrature().size()),
219 normal_vectors(fe_values.get_quadrature().size())
227 template <
int dim,
typename Number>
232 if (present_cell == cell)
235 fe_values.reinit(present_cell);
236 for (
unsigned int q=0; q<fe_values.get_quadrature().size(); ++q)
239 for (
unsigned int d=0; d<dim; ++d)
240 for (
unsigned int e=0; e<dim; ++e)
241 inverse_jacobians[q][d][e] = fe_values.inverse_jacobian(q)[e][d];
243 for (
unsigned int d=0; d<dim; ++d)
244 quadrature_points[q][d] = fe_values.quadrature_point(q)[d];
246 for (
unsigned int d=0; d<dim; ++d)
247 normal_vectors[q][d] = fe_values.normal_vector(q)[d];
249 jxw_values[q] = fe_values.JxW(q);
255 template <
int dim,
typename Number>
260 return present_cell != typename ::Triangulation<dim>::cell_iterator();
265 template <
int dim,
typename Number>
267 typename ::Triangulation<dim>::cell_iterator
270 return fe_values.get_cell();
275 template <
int dim,
typename Number>
277 const ::FEValues<dim> &
285 template <
int dim,
typename Number>
290 return inverse_jacobians;
295 template <
int dim,
typename Number>
300 return normal_vectors;
305 template <
int dim,
typename Number>
310 return quadrature_points;
315 template <
int dim,
typename Number>
325 template <
int dim,
typename Number>
330 return quadrature_1d;
338 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
FE_Nothing< dim > fe_dummy
typename ::Triangulation< dim >::cell_iterator present_cell
const AlignedVector< VectorizedArray< Number > > & get_JxW_values() const
bool is_initialized() const
const AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > & get_inverse_jacobians() const
UpdateFlags get_update_flags() const
Transformed quadrature points.
::FEValues< dim > fe_values
typename ::Triangulation< dim >::cell_iterator get_cell() const
AlignedVector< Point< dim, VectorizedArray< Number > > > quadrature_points
AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > inverse_jacobians
#define Assert(cond, exc)
const AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > & get_normal_vectors() const
AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > normal_vectors
const Quadrature< 1 > quadrature_1d
const Quadrature< 1 > & get_quadrature() const
AlignedVector< VectorizedArray< Number > > jxw_values
const ::FEValues< dim > & get_fe_values() const
Gradient of volume element.
const AlignedVector< Point< dim, VectorizedArray< Number > > > & get_quadrature_points() const
static ::ExceptionBase & ExcNotImplemented()
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)