 |
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
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
61 template <
int dim,
typename Number,
typename VectorizedArrayType>
66 "Type of Number and of VectorizedArrayType do not match.");
91 reinit(typename ::Triangulation<dim>::cell_iterator cell);
103 typename ::Triangulation<dim>::cell_iterator
112 const ::FEValues<dim> &
164 template <
int dim,
typename Number,
typename VectorizedArrayType>
175 , quadrature_1d(quadrature)
201 template <
int dim,
typename Number,
typename VectorizedArrayType>
212 template <
int dim,
typename Number,
typename VectorizedArrayType>
215 typename ::Triangulation<dim>::cell_iterator cell)
217 if (present_cell == cell)
220 fe_values.reinit(present_cell);
221 for (
unsigned int q = 0; q < fe_values.get_quadrature().size(); ++q)
224 mapping_info_storage.JxW_values[q] = fe_values.JxW(q);
229 for (
unsigned int d = 0;
d < dim; ++
d)
230 for (
unsigned int e = 0;
e < dim; ++
e)
231 mapping_info_storage.jacobians[0][q][
d][
e] = jac[
d][
e];
234 for (
unsigned int d = 0;
d < dim; ++
d)
235 mapping_info_storage.quadrature_points[q][
d] =
236 fe_values.quadrature_point(q)[
d];
239 for (
unsigned int d = 0;
d < dim; ++
d)
240 mapping_info_storage.normal_vectors[q][
d] =
241 fe_values.normal_vector(q)[
d];
242 mapping_info_storage.normals_times_jacobians[0][q] =
243 mapping_info_storage.normal_vectors[q] *
244 mapping_info_storage.jacobians[0][q];
251 template <
int dim,
typename Number,
typename VectorizedArrayType>
256 return present_cell !=
257 typename ::Triangulation<dim>::cell_iterator();
262 template <
int dim,
typename Number,
typename VectorizedArrayType>
263 inline typename ::Triangulation<dim>::cell_iterator
266 return fe_values.get_cell();
271 template <
int dim,
typename Number,
typename VectorizedArrayType>
272 inline const ::FEValues<dim> &
280 template <
int dim,
typename Number,
typename VectorizedArrayType>
285 return mapping_info_storage;
290 template <
int dim,
typename Number,
typename VectorizedArrayType>
295 return quadrature_1d;
bool is_initialized() const
const Quadrature< 1 > & get_quadrature() const
@ update_quadrature_points
Transformed quadrature points.
FE_Nothing< dim > fe_dummy
typename ::Triangulation< dim >::cell_iterator present_cell
static ::ExceptionBase & ExcNotImplemented()
AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > > normal_vectors
AlignedVector< Point< spacedim, VectorizedArrayType > > quadrature_points
AlignedVector< unsigned int > quadrature_point_offsets
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void resize(const size_type size_in)
MappingInfoStorage< dim, dim, Number, VectorizedArrayType > mapping_info_storage
::FEValues< dim > fe_values
const Quadrature< 1 > quadrature_1d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit(typename ::Triangulation< dim >::cell_iterator cell)
typename ::Triangulation< dim >::cell_iterator get_cell() const
AlignedVector< Tensor< 1, spacedim, VectorizedArrayType > > normals_times_jacobians[2]
std::vector< QuadratureDescriptor > descriptor
AlignedVector< VectorizedArrayType > JxW_values
#define DEAL_II_NAMESPACE_OPEN
@ update_jacobians
Volume element.
MappingDataOnTheFly(const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
const ::FEValues< dim > & get_fe_values() const
@ update_jacobian_grads
Gradient of volume element.
@ update_JxW_values
Transformed quadrature weights.
@ update_normal_vectors
Normal vectors.
#define Assert(cond, exc)
const MappingInfoStorage< dim, dim, Number, VectorizedArrayType > & get_data_storage() const
AlignedVector< unsigned int > data_index_offsets
MappingDataOnTheFly(const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
AlignedVector< Tensor< 2, spacedim, VectorizedArrayType > > jacobians[2]
#define DEAL_II_NAMESPACE_CLOSE
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)