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)
65 internal::MatrixFreeFunctions::tensor_general,
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)
92 internal::MatrixFreeFunctions::tensor_general,
97 Number>::integrate(shape_info,
117 int n_q_points_1d = 0,
119 struct Factory : Default<dim, n_components, Number>
127 template <
int n_q_po
ints_1d,
int dim,
int n_components,
typename Number>
128 struct Factory<dim, n_components, Number, 0, 10, n_q_points_1d>
129 : Default<dim, n_components, Number>
137 template <
int degree,
148 typename
std::enable_if<n_q_points_1d == degree + 3>::type>
149 : Default<dim, n_components, Number>
155 template <
int degree,
160 struct Factory<dim, n_components, Number, 0, degree, n_q_points_1d>
165 Number * values_dofs_actual,
166 Number * values_quad,
167 Number * gradients_quad,
168 Number * hessians_quad,
169 Number * scratch_data,
170 const bool evaluate_values,
171 const bool evaluate_gradients,
172 const bool evaluate_hessians)
174 const unsigned int runtime_degree = shape_info.
fe_degree;
175 constexpr
unsigned int start_n_q_points = degree + 1;
176 if (runtime_degree == degree)
177 Factory<dim, n_components, Number, 1, degree, start_n_q_points>::
188 Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
203 Number * values_dofs_actual,
204 Number * values_quad,
205 Number * gradients_quad,
206 Number * scratch_data,
207 const bool integrate_values,
208 const bool integrate_gradients)
210 const int runtime_degree = shape_info.
fe_degree;
211 constexpr
unsigned int start_n_q_points = degree + 1;
212 if (runtime_degree == degree)
213 Factory<dim, n_components, Number, 1, degree, start_n_q_points>::
214 integrate(shape_info,
220 integrate_gradients);
222 Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
223 integrate(shape_info,
229 integrate_gradients);
237 template <
int degree,
248 typename std::enable_if<(n_q_points_1d < degree + 3)>::type>
253 Number * values_dofs_actual,
254 Number * values_quad,
255 Number * gradients_quad,
256 Number * hessians_quad,
257 Number * scratch_data,
258 const bool evaluate_values,
259 const bool evaluate_gradients,
260 const bool evaluate_hessians)
263 if (runtime_n_q_points_1d == n_q_points_1d)
265 if (n_q_points_1d == degree + 1 &&
267 internal::MatrixFreeFunctions::tensor_symmetric_collocation)
279 else if (degree < n_q_points_1d)
285 Number>::evaluate(shape_info,
296 internal::MatrixFreeFunctions::tensor_symmetric,
301 Number>::evaluate(shape_info,
312 Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
327 Number * values_dofs_actual,
328 Number * values_quad,
329 Number * gradients_quad,
330 Number * scratch_data,
331 const bool integrate_values,
332 const bool integrate_gradients)
335 if (runtime_n_q_points_1d == n_q_points_1d)
337 if (n_q_points_1d == degree + 1 &&
339 internal::MatrixFreeFunctions::tensor_symmetric_collocation)
350 else if (degree < n_q_points_1d)
356 Number>::integrate(shape_info,
366 internal::MatrixFreeFunctions::tensor_symmetric,
371 Number>::integrate(shape_info,
381 Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
382 integrate(shape_info,
388 integrate_gradients);
398 template <
int dim,
int n_components,
typename Number>
400 symmetric_selector_evaluate(
402 Number * values_dofs_actual,
403 Number * values_quad,
404 Number * gradients_quad,
405 Number * hessians_quad,
406 Number * scratch_data,
407 const bool evaluate_values,
408 const bool evaluate_gradients,
409 const bool evaluate_hessians)
412 internal::MatrixFreeFunctions::tensor_symmetric,
414 Factory<dim, n_components, Number>::evaluate(shape_info,
431 template <
int dim,
int n_components,
typename Number>
433 symmetric_selector_integrate(
435 Number * values_dofs_actual,
436 Number * values_quad,
437 Number * gradients_quad,
438 Number * scratch_data,
439 const bool integrate_values,
440 const bool integrate_gradients)
443 internal::MatrixFreeFunctions::tensor_symmetric,
445 Factory<dim, n_components, Number>::integrate(shape_info,
451 integrate_gradients);
485 Number * values_dofs_actual,
486 Number * values_quad,
487 Number * gradients_quad,
488 Number * hessians_quad,
489 Number * scratch_data,
490 const bool evaluate_values,
491 const bool evaluate_gradients,
492 const bool evaluate_hessians);
503 Number * values_dofs_actual,
504 Number * values_quad,
505 Number * gradients_quad,
506 Number * scratch_data,
507 const bool integrate_values,
508 const bool integrate_gradients);
521 template <
int dim,
int n_q_po
ints_1d,
int n_components,
typename Number>
534 Number * values_dofs_actual,
535 Number * values_quad,
536 Number * gradients_quad,
537 Number * hessians_quad,
538 Number * scratch_data,
539 const bool evaluate_values,
540 const bool evaluate_gradients,
541 const bool evaluate_hessians);
553 Number * values_dofs_actual,
554 Number * values_quad,
555 Number * gradients_quad,
556 Number * scratch_data,
557 const bool integrate_values,
558 const bool integrate_gradients);
572 Number * values_dofs_actual,
573 Number * values_quad,
574 Number * gradients_quad,
575 Number * hessians_quad,
576 Number * scratch_data,
577 const bool evaluate_values,
578 const bool evaluate_gradients,
579 const bool evaluate_hessians)
583 if (fe_degree + 1 == n_q_points_1d &&
585 internal::MatrixFreeFunctions::tensor_symmetric_collocation)
599 else if (fe_degree < n_q_points_1d &&
601 internal::MatrixFreeFunctions::tensor_symmetric)
608 Number>::evaluate(shape_info,
619 internal::MatrixFreeFunctions::tensor_symmetric)
622 internal::MatrixFreeFunctions::tensor_symmetric,
627 Number>::evaluate(shape_info,
638 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
641 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
646 Number>::evaluate(shape_info,
657 internal::MatrixFreeFunctions::truncated_tensor)
660 internal::MatrixFreeFunctions::truncated_tensor,
665 Number>::evaluate(shape_info,
676 internal::MatrixFreeFunctions::tensor_general)
683 Number>::evaluate(shape_info,
707 Number * values_dofs_actual,
708 Number * values_quad,
709 Number * gradients_quad,
710 Number * scratch_data,
711 const bool integrate_values,
712 const bool integrate_gradients)
716 if (fe_degree + 1 == n_q_points_1d &&
718 internal::MatrixFreeFunctions::tensor_symmetric_collocation)
731 else if (fe_degree < n_q_points_1d &&
733 internal::MatrixFreeFunctions::tensor_symmetric)
740 Number>::integrate(shape_info,
750 internal::MatrixFreeFunctions::tensor_symmetric)
753 internal::MatrixFreeFunctions::tensor_symmetric,
758 Number>::integrate(shape_info,
768 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
771 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
776 Number>::integrate(shape_info,
786 internal::MatrixFreeFunctions::truncated_tensor)
789 internal::MatrixFreeFunctions::truncated_tensor,
794 Number>::integrate(shape_info,
804 internal::MatrixFreeFunctions::tensor_general)
811 Number>::integrate(shape_info,
826 template <
int dim,
int dummy,
int n_components,
typename Number>
830 Number * values_dofs_actual,
831 Number * values_quad,
832 Number * gradients_quad,
833 Number * hessians_quad,
834 Number * scratch_data,
835 const bool evaluate_values,
836 const bool evaluate_gradients,
837 const bool evaluate_hessians)
840 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
843 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
848 Number>::evaluate(shape_info,
859 internal::MatrixFreeFunctions::truncated_tensor)
862 internal::MatrixFreeFunctions::truncated_tensor,
867 Number>::evaluate(shape_info,
878 internal::MatrixFreeFunctions::tensor_general)
884 Number>::evaluate(shape_info,
894 internal::EvaluationSelectorImplementation::
895 symmetric_selector_evaluate<dim, n_components, Number>(shape_info,
908 template <
int dim,
int dummy,
int n_components,
typename Number>
912 Number * values_dofs_actual,
913 Number * values_quad,
914 Number * gradients_quad,
915 Number * scratch_data,
916 const bool integrate_values,
917 const bool integrate_gradients)
920 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
923 internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
928 Number>::integrate(shape_info,
938 internal::MatrixFreeFunctions::truncated_tensor)
941 internal::MatrixFreeFunctions::truncated_tensor,
946 Number>::integrate(shape_info,
956 internal::MatrixFreeFunctions::tensor_general)
962 Number>::integrate(shape_info,
971 internal::EvaluationSelectorImplementation::
972 symmetric_selector_integrate<dim, n_components, Number>(
979 integrate_gradients);
983 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)
#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()
static::ExceptionBase & ExcInternalError()