#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping_q_generic.h>
-#include <deal.II/matrix_free/evaluation_kernels.h>
-#include <deal.II/matrix_free/evaluation_selector.h>
+#include <deal.II/matrix_free/evaluation_template_factory.h>
#include <deal.II/matrix_free/mapping_info.h>
DEAL_II_NAMESPACE_OPEN
start_indices,
cell_points.data());
- SelectEvaluator<dim, -1, 0, VectorizedDouble>::evaluate(
+ FEEvaluationFactory<dim, double, VectorizedDouble>::evaluate(
dim,
EvaluationFlags::values | EvaluationFlags::gradients |
(update_flags_cells & update_jacobian_grads ?
// now let the matrix-free evaluators provide us with the
// data on faces
- FEFaceEvaluationImplEvaluateSelector<dim, VectorizedDouble>::
- template run<-1, 0>(dim,
- shape_info,
- cell_points.data(),
- face_quads.data(),
- face_grads.data(),
- scratch_data.data(),
- true,
- true,
- face_no,
- GeometryInfo<dim>::max_children_per_cell,
- faces[face].face_orientation > 8 ?
- faces[face].face_orientation - 8 :
- 0,
- my_data.descriptor[0].face_orientations);
+ FEFaceEvaluationFactory<dim, double, VectorizedDouble>::evaluate(
+ dim,
+ shape_info,
+ cell_points.data(),
+ face_quads.data(),
+ face_grads.data(),
+ scratch_data.data(),
+ true,
+ true,
+ face_no,
+ GeometryInfo<dim>::max_children_per_cell,
+ faces[face].face_orientation > 8 ?
+ faces[face].face_orientation - 8 :
+ 0,
+ my_data.descriptor[0].face_orientations);
if (update_flags_faces & update_quadrature_points)
start_indices,
cell_points.data());
- FEFaceEvaluationImplEvaluateSelector<dim, VectorizedDouble>::
- template run<-1, 0>(
- dim,
- shape_info,
- cell_points.data(),
- face_quads.data(),
- face_grads.data(),
- scratch_data.data(),
- false,
- true,
- faces[face].exterior_face_no,
- faces[face].subface_index,
- faces[face].face_orientation < 8 ?
- faces[face].face_orientation :
- 0,
- my_data.descriptor[0].face_orientations);
+ FEFaceEvaluationFactory<dim, double, VectorizedDouble>::
+ evaluate(dim,
+ shape_info,
+ cell_points.data(),
+ face_quads.data(),
+ face_grads.data(),
+ scratch_data.data(),
+ false,
+ true,
+ faces[face].exterior_face_no,
+ faces[face].subface_index,
+ faces[face].face_orientation < 8 ?
+ faces[face].face_orientation :
+ 0,
+ my_data.descriptor[0].face_orientations);
for (unsigned int q = 0; q < n_points_compute; ++q)
{
/*
* MatrixFree implementation of project() for an arbitrary number of
- * components and arbitrary degree of the FiniteElement.
+ * components of the FiniteElement.
*/
- template <int components,
- int fe_degree,
- int dim,
- typename Number,
- int spacedim>
+ template <int components, int dim, typename Number, int spacedim>
void
project_matrix_free(
const Mapping<dim, spacedim> & mapping,
Assert(dof.get_fe(0).n_components() == function.n_components,
ExcDimensionMismatch(dof.get_fe(0).n_components(),
function.n_components));
- Assert(fe_degree == -1 ||
- dof.get_fe().degree == static_cast<unsigned int>(fe_degree),
- ExcDimensionMismatch(fe_degree, dof.get_fe().degree));
Assert(dof.get_fe(0).n_components() == components,
ExcDimensionMismatch(components, dof.get_fe(0).n_components()));
additional_data);
using MatrixType = MatrixFreeOperators::MassOperator<
dim,
- fe_degree,
- fe_degree + 2,
+ -1,
+ 0,
components,
LinearAlgebra::distributed::Vector<Number>>;
MatrixType mass_matrix;
// account for inhomogeneous constraints
inhomogeneities.update_ghost_values();
- FEEvaluation<dim, fe_degree, fe_degree + 2, components, Number> phi(
- *matrix_free);
+ FEEvaluation<dim, -1, 0, components, Number> phi(*matrix_free);
for (unsigned int cell = 0; cell < matrix_free->n_macro_cells(); ++cell)
{
phi.reinit(cell);
- /**
- * Helper interface. After figuring out the number of components in
- * project_matrix_free_component, we determine the degree of the
- * FiniteElement and call project_matrix_free with the appropriate
- * template arguments.
- */
- template <int components, int dim, typename Number, int spacedim>
- void
- project_matrix_free_degree(
- const Mapping<dim, spacedim> & mapping,
- const DoFHandler<dim, spacedim> &dof,
- const AffineConstraints<Number> &constraints,
- const Quadrature<dim> & quadrature,
- const Function<
- spacedim,
- typename LinearAlgebra::distributed::Vector<Number>::value_type>
- & function,
- LinearAlgebra::distributed::Vector<Number> &work_result,
- const bool enforce_zero_boundary,
- const Quadrature<dim - 1> & q_boundary,
- const bool project_to_boundary_first)
- {
- switch (dof.get_fe().degree)
- {
- case 1:
- project_matrix_free<components, 1>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
- break;
-
- case 2:
- project_matrix_free<components, 2>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
- break;
-
- case 3:
- project_matrix_free<components, 3>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
- break;
-
- default:
- project_matrix_free<components, -1>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
- }
- }
-
-
-
// Helper interface for the matrix-free implementation of project().
// Used to determine the number of components.
template <int dim, typename Number, int spacedim>
switch (dof.get_fe(0).n_components())
{
case 1:
- project_matrix_free_degree<1>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
+ project_matrix_free<1>(mapping,
+ dof,
+ constraints,
+ quadrature,
+ function,
+ work_result,
+ enforce_zero_boundary,
+ q_boundary,
+ project_to_boundary_first);
break;
case 2:
- project_matrix_free_degree<2>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
+ project_matrix_free<2>(mapping,
+ dof,
+ constraints,
+ quadrature,
+ function,
+ work_result,
+ enforce_zero_boundary,
+ q_boundary,
+ project_to_boundary_first);
break;
case 3:
- project_matrix_free_degree<3>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
+ project_matrix_free<3>(mapping,
+ dof,
+ constraints,
+ quadrature,
+ function,
+ work_result,
+ enforce_zero_boundary,
+ q_boundary,
+ project_to_boundary_first);
break;
case 4:
- project_matrix_free_degree<4>(mapping,
- dof,
- constraints,
- quadrature,
- function,
- work_result,
- enforce_zero_boundary,
- q_boundary,
- project_to_boundary_first);
+ project_matrix_free<4>(mapping,
+ dof,
+ constraints,
+ quadrature,
+ function,
+ work_result,
+ enforce_zero_boundary,
+ q_boundary,
+ project_to_boundary_first);
break;
default:
vec_result);
}
- template <int dim, typename VectorType, int spacedim, int fe_degree>
+ template <int dim, typename VectorType, int spacedim>
void
project_parallel(
const Mapping<dim, spacedim> & mapping,
ExcDimensionMismatch(dof.get_fe(0).n_components(), 1));
Assert(vec_result.size() == dof.n_dofs(),
ExcDimensionMismatch(vec_result.size(), dof.n_dofs()));
- Assert(fe_degree == -1 ||
- dof.get_fe().degree == static_cast<unsigned int>(fe_degree),
- ExcDimensionMismatch(fe_degree, dof.get_fe().degree));
// set up mass matrix and right hand side
typename MatrixFree<dim, Number>::AdditionalData additional_data;
constraints,
QGauss<1>(dof.get_fe().degree + 2),
additional_data);
- using MatrixType = MatrixFreeOperators::MassOperator<
- dim,
- fe_degree,
- fe_degree + 2,
- 1,
- LinearAlgebra::distributed::Vector<Number>>;
+ using MatrixType = MatrixFreeOperators::
+ MassOperator<dim, -1, 0, 1, LinearAlgebra::distributed::Vector<Number>>;
MatrixType mass_matrix;
mass_matrix.initialize(matrix_free);
mass_matrix.compute_diagonal();
- template <int dim,
- typename VectorType,
- int spacedim,
- int fe_degree,
- int n_q_points_1d>
+ template <int dim, typename VectorType, int spacedim>
void
project_parallel(
std::shared_ptr<const MatrixFree<dim, typename VectorType::value_type>>
ExcDimensionMismatch(dof.get_fe(0).n_components(), 1));
Assert(vec_result.size() == dof.n_dofs(),
ExcDimensionMismatch(vec_result.size(), dof.n_dofs()));
- Assert(fe_degree == -1 ||
- dof.get_fe().degree == static_cast<unsigned int>(fe_degree),
- ExcDimensionMismatch(fe_degree, dof.get_fe().degree));
- using MatrixType = MatrixFreeOperators::MassOperator<
- dim,
- fe_degree,
- n_q_points_1d,
- 1,
- LinearAlgebra::distributed::Vector<Number>>;
+ using MatrixType = MatrixFreeOperators::
+ MassOperator<dim, -1, 0, 1, LinearAlgebra::distributed::Vector<Number>>;
MatrixType mass_matrix;
mass_matrix.initialize(matrix_free, {fe_component});
mass_matrix.compute_diagonal();
// assemble right hand side:
{
- FEEvaluation<dim, fe_degree, n_q_points_1d, 1, Number> fe_eval(
- *matrix_free, fe_component);
- const unsigned int n_cells = matrix_free->n_macro_cells();
+ FEEvaluation<dim, -1, 0, 1, Number> fe_eval(*matrix_free, fe_component);
+ const unsigned int n_cells = matrix_free->n_cell_batches();
const unsigned int n_q_points = fe_eval.n_q_points;
for (unsigned int cell = 0; cell < n_cells; ++cell)
const unsigned int)> & func,
VectorType & vec_result)
{
- switch (dof.get_fe().degree)
- {
- case 1:
- internal::project_parallel<dim, VectorType, spacedim, 1>(
- mapping, dof, constraints, quadrature, func, vec_result);
- break;
- case 2:
- internal::project_parallel<dim, VectorType, spacedim, 2>(
- mapping, dof, constraints, quadrature, func, vec_result);
- break;
- case 3:
- internal::project_parallel<dim, VectorType, spacedim, 3>(
- mapping, dof, constraints, quadrature, func, vec_result);
- break;
- default:
- internal::project_parallel<dim, VectorType, spacedim, -1>(
- mapping, dof, constraints, quadrature, func, vec_result);
- }
+ internal::project_parallel<dim, VectorType, spacedim>(
+ mapping, dof, constraints, quadrature, func, vec_result);
}
typename VectorType::value_type,
VectorizedArray<typename VectorType::value_type>>> matrix_free,
const AffineConstraints<typename VectorType::value_type> &constraints,
- const unsigned int n_q_points_1d,
+ const unsigned int,
const std::function<VectorizedArray<typename VectorType::value_type>(
const unsigned int,
const unsigned int)> &func,
VectorType & vec_result,
const unsigned int fe_component)
{
- const unsigned int fe_degree =
- matrix_free->get_dof_handler(fe_component).get_fe().degree;
-
- if (fe_degree + 1 == n_q_points_1d)
- switch (fe_degree)
- {
- case 1:
- internal::project_parallel<dim, VectorType, dim, 1, 2>(
- matrix_free, constraints, func, vec_result, fe_component);
- break;
- case 2:
- internal::project_parallel<dim, VectorType, dim, 2, 3>(
- matrix_free, constraints, func, vec_result, fe_component);
- break;
- case 3:
- internal::project_parallel<dim, VectorType, dim, 3, 4>(
- matrix_free, constraints, func, vec_result, fe_component);
- break;
- default:
- internal::project_parallel<dim, VectorType, dim, -1, 0>(
- matrix_free, constraints, func, vec_result, fe_component);
- }
- else
- internal::project_parallel<dim, VectorType, dim, -1, 0>(
- matrix_free, constraints, func, vec_result, fe_component);
+ internal::project_parallel<dim, VectorType, dim>(
+ matrix_free, constraints, func, vec_result, fe_component);
}