17 #ifndef dealii_matrix_free_evaluation_selector_h 18 #define dealii_matrix_free_evaluation_selector_h 20 #include <deal.II/matrix_free/evaluation_kernels.h> 22 DEAL_II_NAMESPACE_OPEN
27 namespace EvaluationSelectorImplementation
49 template <
int dim,
int n_components,
typename Number>
55 Number * values_dofs_actual,
57 Number * gradients_quad,
58 Number * hessians_quad,
59 Number * scratch_data,
60 const bool evaluate_values,
61 const bool evaluate_gradients,
62 const bool evaluate_hessians)
70 Number>::evaluate(shape_info,
84 Number * values_dofs_actual,
86 Number * gradients_quad,
87 Number * scratch_data,
88 const bool integrate_values,
89 const bool integrate_gradients,
90 const bool sum_into_values_array =
false)
98 Number>::integrate(shape_info,
105 sum_into_values_array);
118 int n_q_points_1d = 0,
120 struct Factory : Default<dim, n_components, Number>
128 template <
int n_q_po
ints_1d,
int dim,
int n_components,
typename Number>
129 struct Factory<dim,
n_components, Number, 0, 10, n_q_points_1d>
130 : Default<dim, n_components, Number>
138 template <
int degree,
149 typename
std::enable_if<n_q_points_1d == degree + 3>::type>
150 : Default<dim, n_components, Number>
156 template <
int degree,
161 struct Factory<dim,
n_components, Number, 0, degree, n_q_points_1d>
166 Number * values_dofs_actual,
167 Number * values_quad,
168 Number * gradients_quad,
169 Number * hessians_quad,
170 Number * scratch_data,
171 const bool evaluate_values,
172 const bool evaluate_gradients,
173 const bool evaluate_hessians)
175 const unsigned int runtime_degree = shape_info.
fe_degree;
176 constexpr
unsigned int start_n_q_points = degree + 1;
177 if (runtime_degree == degree)
178 Factory<dim, n_components, Number, 1, degree, start_n_q_points>::
189 Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
204 Number * values_dofs_actual,
205 Number * values_quad,
206 Number * gradients_quad,
207 Number * scratch_data,
208 const bool integrate_values,
209 const bool integrate_gradients,
210 const bool sum_into_values_array =
false)
212 const int runtime_degree = shape_info.
fe_degree;
213 constexpr
unsigned int start_n_q_points = degree + 1;
214 if (runtime_degree == degree)
215 Factory<dim, n_components, Number, 1, degree, start_n_q_points>::
216 integrate(shape_info,
223 sum_into_values_array);
225 Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
226 integrate(shape_info,
233 sum_into_values_array);
241 template <
int degree,
252 typename
std::enable_if<(n_q_points_1d < degree + 3)>::type>
261 static constexpr
bool use_collocation =
262 n_q_points_1d > degree &&n_q_points_1d <= 3 * degree / 2 + 1 &&
268 Number * values_dofs_actual,
269 Number * values_quad,
270 Number * gradients_quad,
271 Number * hessians_quad,
272 Number * scratch_data,
273 const bool evaluate_values,
274 const bool evaluate_gradients,
275 const bool evaluate_hessians)
278 if (runtime_n_q_points_1d == n_q_points_1d)
280 if (n_q_points_1d == degree + 1 &&
284 FEEvaluationImplCollocation<dim, degree, n_components, Number>::
294 else if (use_collocation)
300 Number>::evaluate(shape_info,
316 Number>::evaluate(shape_info,
327 Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
342 Number * values_dofs_actual,
343 Number * values_quad,
344 Number * gradients_quad,
345 Number * scratch_data,
346 const bool integrate_values,
347 const bool integrate_gradients,
348 const bool sum_into_values_array)
351 if (runtime_n_q_points_1d == n_q_points_1d)
353 if (n_q_points_1d == degree + 1 &&
357 FEEvaluationImplCollocation<dim, degree, n_components, Number>::
358 integrate(shape_info,
365 sum_into_values_array);
366 else if (use_collocation)
372 Number>::integrate(shape_info,
379 sum_into_values_array);
387 Number>::integrate(shape_info,
394 sum_into_values_array);
397 Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
398 integrate(shape_info,
405 sum_into_values_array);
415 template <
int dim,
int n_components,
typename Number>
417 symmetric_selector_evaluate(
419 Number * values_dofs_actual,
420 Number * values_quad,
421 Number * gradients_quad,
422 Number * hessians_quad,
423 Number * scratch_data,
424 const bool evaluate_values,
425 const bool evaluate_gradients,
426 const bool evaluate_hessians)
431 Factory<dim, n_components, Number>::evaluate(shape_info,
448 template <
int dim,
int n_components,
typename Number>
450 symmetric_selector_integrate(
452 Number * values_dofs_actual,
453 Number * values_quad,
454 Number * gradients_quad,
455 Number * scratch_data,
456 const bool integrate_values,
457 const bool integrate_gradients,
458 const bool sum_into_values_array =
false)
463 Factory<dim, n_components, Number>::integrate(shape_info,
470 sum_into_values_array);
503 n_q_points_1d > fe_degree &&n_q_points_1d <= 3 * fe_degree / 2 + 1 &&
515 Number * values_dofs_actual,
516 Number * values_quad,
517 Number * gradients_quad,
518 Number * hessians_quad,
519 Number * scratch_data,
520 const bool evaluate_values,
521 const bool evaluate_gradients,
522 const bool evaluate_hessians);
533 Number * values_dofs_actual,
534 Number * values_quad,
535 Number * gradients_quad,
536 Number * scratch_data,
537 const bool integrate_values,
538 const bool integrate_gradients,
539 const bool sum_into_values_array =
false);
552 template <
int dim,
int n_q_po
ints_1d,
int n_components,
typename Number>
565 Number * values_dofs_actual,
566 Number * values_quad,
567 Number * gradients_quad,
568 Number * hessians_quad,
569 Number * scratch_data,
570 const bool evaluate_values,
571 const bool evaluate_gradients,
572 const bool evaluate_hessians);
584 Number * values_dofs_actual,
585 Number * values_quad,
586 Number * gradients_quad,
587 Number * scratch_data,
588 const bool integrate_values,
589 const bool integrate_gradients,
590 const bool sum_into_values_array =
false);
604 Number * values_dofs_actual,
605 Number * values_quad,
606 Number * gradients_quad,
607 Number * hessians_quad,
608 Number * scratch_data,
609 const bool evaluate_values,
610 const bool evaluate_gradients,
611 const bool evaluate_hessians)
615 if (fe_degree + 1 == n_q_points_1d &&
620 FEEvaluationImplCollocation<dim, fe_degree, n_components, Number>::
641 Number>::evaluate(shape_info,
660 Number>::evaluate(shape_info,
679 Number>::evaluate(shape_info,
698 Number>::evaluate(shape_info,
716 Number>::evaluate(shape_info,
740 Number * values_dofs_actual,
741 Number * values_quad,
742 Number * gradients_quad,
743 Number * scratch_data,
744 const bool integrate_values,
745 const bool integrate_gradients,
746 const bool sum_into_values_array)
750 if (fe_degree + 1 == n_q_points_1d &&
755 FEEvaluationImplCollocation<dim, fe_degree, n_components, Number>::
756 integrate(shape_info,
763 sum_into_values_array);
775 Number>::integrate(shape_info,
782 sum_into_values_array);
793 Number>::integrate(shape_info,
800 sum_into_values_array);
811 Number>::integrate(shape_info,
818 sum_into_values_array);
829 Number>::integrate(shape_info,
836 sum_into_values_array);
846 Number>::integrate(shape_info,
853 sum_into_values_array);
861 template <
int dim,
int dummy,
int n_components,
typename Number>
865 Number * values_dofs_actual,
866 Number * values_quad,
867 Number * gradients_quad,
868 Number * hessians_quad,
869 Number * scratch_data,
870 const bool evaluate_values,
871 const bool evaluate_gradients,
872 const bool evaluate_hessians)
883 Number>::evaluate(shape_info,
902 Number>::evaluate(shape_info,
919 Number>::evaluate(shape_info,
929 internal::EvaluationSelectorImplementation::
930 symmetric_selector_evaluate<dim, n_components, Number>(shape_info,
943 template <
int dim,
int dummy,
int n_components,
typename Number>
947 Number * values_dofs_actual,
948 Number * values_quad,
949 Number * gradients_quad,
950 Number * scratch_data,
951 const bool integrate_values,
952 const bool integrate_gradients,
953 const bool sum_into_values_array)
964 Number>::integrate(shape_info,
971 sum_into_values_array);
982 Number>::integrate(shape_info,
989 sum_into_values_array);
998 Number>::integrate(shape_info,
1004 integrate_gradients,
1005 sum_into_values_array);
1007 internal::EvaluationSelectorImplementation::
1008 symmetric_selector_integrate<dim, n_components, Number>(
1015 integrate_gradients,
1016 sum_into_values_array);
1020 DEAL_II_NAMESPACE_CLOSE
static void integrate(const internal::MatrixFreeFunctions::ShapeInfo< Number > &shape_info, Number *values_dofs_actual, Number *values_quad, Number *gradients_quad, Number *scratch_data, const bool integrate_values, const bool integrate_gradients, const bool sum_into_values_array=false)
static constexpr bool use_collocation
#define AssertThrow(cond, exc)
static void evaluate(const internal::MatrixFreeFunctions::ShapeInfo< Number > &shape_info, Number *values_dofs_actual, Number *values_quad, Number *gradients_quad, Number *hessians_quad, Number *scratch_data, const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians)
#define Assert(cond, exc)
unsigned int n_q_points_1d
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcInternalError()