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\}}\)
evaluation_selector.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_matrix_free_evaluation_selector_h
18 #define dealii_matrix_free_evaluation_selector_h
19 
20 #include <deal.II/base/config.h>
21 
23 
25 
26 #ifndef DOXYGEN
27 namespace internal
28 {
29  namespace EvaluationSelectorImplementation
30  {
31  // The following classes serve the purpose of choosing the correct template
32  // specialization of the FEEvaluationImpl* classes in case fe_degree
33  // and n_q_points_1d are only given as runtime parameters.
34  // The logic is the following:
35  // 1. Start with fe_degree=0, n_q_points_1d=0 and DEPTH=0.
36  // 2. If the current assumption on fe_degree doesn't match the runtime
37  // parameter, increase fe_degree by one and try again.
38  // If fe_degree==10 use the class Default which serves as a fallback.
39  // 3. After fixing the fe_degree, DEPTH is increased (DEPTH=1) and we start
40  // with
41  // n_q_points=fe_degree+1.
42  // 4. If the current assumption on n_q_points_1d doesn't match the runtime
43  // parameter, increase n_q_points_1d by one and try again.
44  // If n_q_points_1d==degree+3 use the class Default which serves as a
45  // fallback.
46 
51  template <int dim, int n_components, typename Number>
52  struct Default
53  {
54  static inline void
55  evaluate(
57  Number * values_dofs_actual,
58  Number * values_quad,
59  Number * gradients_quad,
60  Number * hessians_quad,
61  Number * scratch_data,
62  const bool evaluate_values,
63  const bool evaluate_gradients,
64  const bool evaluate_hessians)
65  {
68  dim,
69  -1,
70  0,
72  Number>::evaluate(shape_info,
73  values_dofs_actual,
74  values_quad,
75  gradients_quad,
76  hessians_quad,
77  scratch_data,
78  evaluate_values,
79  evaluate_gradients,
80  evaluate_hessians);
81  }
82 
83  static inline void
84  integrate(
86  Number * values_dofs_actual,
87  Number * values_quad,
88  Number * gradients_quad,
89  Number * scratch_data,
90  const bool integrate_values,
91  const bool integrate_gradients,
92  const bool sum_into_values_array = false)
93  {
96  dim,
97  -1,
98  0,
100  Number>::integrate(shape_info,
101  values_dofs_actual,
102  values_quad,
103  gradients_quad,
104  scratch_data,
105  integrate_values,
106  integrate_gradients,
107  sum_into_values_array);
108  }
109  };
110 
111 
115  template <int dim,
116  int n_components,
117  typename Number,
118  int DEPTH = 0,
119  int degree = 0,
120  int n_q_points_1d = 0,
121  class Enable = void>
122  struct Factory : Default<dim, n_components, Number>
123  {};
124 
130  template <int n_q_points_1d, int dim, int n_components, typename Number>
131  struct Factory<dim, n_components, Number, 0, 10, n_q_points_1d>
132  : Default<dim, n_components, Number>
133  {};
134 
140  template <int degree,
141  int n_q_points_1d,
142  int dim,
143  int n_components,
144  typename Number>
145  struct Factory<dim,
146  n_components,
147  Number,
148  1,
149  degree,
150  n_q_points_1d,
151  typename std::enable_if<n_q_points_1d == degree + 3>::type>
152  : Default<dim, n_components, Number>
153  {};
154 
158  template <int degree,
159  int n_q_points_1d,
160  int dim,
161  int n_components,
162  typename Number>
163  struct Factory<dim, n_components, Number, 0, degree, n_q_points_1d>
164  {
165  static inline void
166  evaluate(
168  Number * values_dofs_actual,
169  Number * values_quad,
170  Number * gradients_quad,
171  Number * hessians_quad,
172  Number * scratch_data,
173  const bool evaluate_values,
174  const bool evaluate_gradients,
175  const bool evaluate_hessians)
176  {
177  const unsigned int runtime_degree = shape_info.data.front().fe_degree;
178  constexpr unsigned int start_n_q_points = degree + 1;
179  if (runtime_degree == degree)
181  evaluate(shape_info,
182  values_dofs_actual,
183  values_quad,
184  gradients_quad,
185  hessians_quad,
186  scratch_data,
187  evaluate_values,
188  evaluate_gradients,
189  evaluate_hessians);
190  else
192  evaluate(shape_info,
193  values_dofs_actual,
194  values_quad,
195  gradients_quad,
196  hessians_quad,
197  scratch_data,
198  evaluate_values,
199  evaluate_gradients,
200  evaluate_hessians);
201  }
202 
203  static inline void
204  integrate(
206  Number * values_dofs_actual,
207  Number * values_quad,
208  Number * gradients_quad,
209  Number * scratch_data,
210  const bool integrate_values,
211  const bool integrate_gradients,
212  const bool sum_into_values_array = false)
213  {
214  const int runtime_degree = shape_info.data.front().fe_degree;
215  constexpr unsigned int start_n_q_points = degree + 1;
216  if (runtime_degree == degree)
218  integrate(shape_info,
219  values_dofs_actual,
220  values_quad,
221  gradients_quad,
222  scratch_data,
223  integrate_values,
224  integrate_gradients,
225  sum_into_values_array);
226  else
228  integrate(shape_info,
229  values_dofs_actual,
230  values_quad,
231  gradients_quad,
232  scratch_data,
233  integrate_values,
234  integrate_gradients,
235  sum_into_values_array);
236  }
237  };
238 
243  template <int degree,
244  int n_q_points_1d,
245  int dim,
246  int n_components,
247  typename Number>
248  struct Factory<dim,
249  n_components,
250  Number,
251  1,
252  degree,
253  n_q_points_1d,
254  typename std::enable_if<(n_q_points_1d < degree + 3)>::type>
255  {
263  static constexpr bool use_collocation =
264  n_q_points_1d > degree &&n_q_points_1d <= 3 * degree / 2 + 1 &&
265  n_q_points_1d < 200;
266 
267  static inline void
268  evaluate(
270  Number * values_dofs_actual,
271  Number * values_quad,
272  Number * gradients_quad,
273  Number * hessians_quad,
274  Number * scratch_data,
275  const bool evaluate_values,
276  const bool evaluate_gradients,
277  const bool evaluate_hessians)
278  {
279  const int runtime_n_q_points_1d = shape_info.data.front().n_q_points_1d;
280  if (runtime_n_q_points_1d == n_q_points_1d)
281  {
282  if (n_q_points_1d == degree + 1 &&
283  shape_info.element_type ==
285  internal::
287  evaluate(shape_info,
288  values_dofs_actual,
289  values_quad,
290  gradients_quad,
291  hessians_quad,
292  scratch_data,
293  evaluate_values,
294  evaluate_gradients,
295  evaluate_hessians);
296  else if (use_collocation)
298  dim,
299  degree,
300  n_q_points_1d,
301  n_components,
302  Number>::evaluate(shape_info,
303  values_dofs_actual,
304  values_quad,
305  gradients_quad,
306  hessians_quad,
307  scratch_data,
308  evaluate_values,
309  evaluate_gradients,
310  evaluate_hessians);
311  else
314  dim,
315  degree,
316  n_q_points_1d,
317  n_components,
318  Number>::evaluate(shape_info,
319  values_dofs_actual,
320  values_quad,
321  gradients_quad,
322  hessians_quad,
323  scratch_data,
324  evaluate_values,
325  evaluate_gradients,
326  evaluate_hessians);
327  }
328  else
330  evaluate(shape_info,
331  values_dofs_actual,
332  values_quad,
333  gradients_quad,
334  hessians_quad,
335  scratch_data,
336  evaluate_values,
337  evaluate_gradients,
338  evaluate_hessians);
339  }
340 
341  static inline void
342  integrate(
344  Number * values_dofs_actual,
345  Number * values_quad,
346  Number * gradients_quad,
347  Number * scratch_data,
348  const bool integrate_values,
349  const bool integrate_gradients,
350  const bool sum_into_values_array)
351  {
352  const int runtime_n_q_points_1d = shape_info.data.front().n_q_points_1d;
353  if (runtime_n_q_points_1d == n_q_points_1d)
354  {
355  if (n_q_points_1d == degree + 1 &&
356  shape_info.element_type ==
358  internal::
360  integrate(shape_info,
361  values_dofs_actual,
362  values_quad,
363  gradients_quad,
364  scratch_data,
365  integrate_values,
366  integrate_gradients,
367  sum_into_values_array);
368  else if (use_collocation)
370  dim,
371  degree,
372  n_q_points_1d,
373  n_components,
374  Number>::integrate(shape_info,
375  values_dofs_actual,
376  values_quad,
377  gradients_quad,
378  scratch_data,
379  integrate_values,
380  integrate_gradients,
381  sum_into_values_array);
382  else
385  dim,
386  degree,
387  n_q_points_1d,
388  n_components,
389  Number>::integrate(shape_info,
390  values_dofs_actual,
391  values_quad,
392  gradients_quad,
393  scratch_data,
394  integrate_values,
395  integrate_gradients,
396  sum_into_values_array);
397  }
398  else
400  integrate(shape_info,
401  values_dofs_actual,
402  values_quad,
403  gradients_quad,
404  scratch_data,
405  integrate_values,
406  integrate_gradients,
407  sum_into_values_array);
408  }
409  };
410 
411 
412 
417  template <int dim, int n_components, typename Number>
418  void
419  symmetric_selector_evaluate(
421  Number * values_dofs_actual,
422  Number * values_quad,
423  Number * gradients_quad,
424  Number * hessians_quad,
425  Number * scratch_data,
426  const bool evaluate_values,
427  const bool evaluate_gradients,
428  const bool evaluate_hessians)
429  {
430  Assert(shape_info.element_type <=
432  ExcInternalError());
434  values_dofs_actual,
435  values_quad,
436  gradients_quad,
437  hessians_quad,
438  scratch_data,
439  evaluate_values,
440  evaluate_gradients,
441  evaluate_hessians);
442  }
443 
444 
445 
450  template <int dim, int n_components, typename Number>
451  void
452  symmetric_selector_integrate(
454  Number * values_dofs_actual,
455  Number * values_quad,
456  Number * gradients_quad,
457  Number * scratch_data,
458  const bool integrate_values,
459  const bool integrate_gradients,
460  const bool sum_into_values_array = false)
461  {
462  Assert(shape_info.element_type <=
464  ExcInternalError());
466  values_dofs_actual,
467  values_quad,
468  gradients_quad,
469  scratch_data,
470  integrate_values,
471  integrate_gradients,
472  sum_into_values_array);
473  }
474  } // namespace EvaluationSelectorImplementation
475 } // namespace internal
476 #endif
477 
478 
490 template <int dim,
491  int fe_degree,
492  int n_q_points_1d,
493  int n_components,
494  typename Number>
496 {
504  static constexpr bool use_collocation =
505  n_q_points_1d > fe_degree &&n_q_points_1d <= 3 * fe_degree / 2 + 1 &&
506  n_q_points_1d < 200;
507 
515  static void
517  Number * values_dofs_actual,
518  Number * values_quad,
519  Number * gradients_quad,
520  Number * hessians_quad,
521  Number * scratch_data,
522  const bool evaluate_values,
523  const bool evaluate_gradients,
524  const bool evaluate_hessians);
525 
533  static void
535  Number * values_dofs_actual,
536  Number * values_quad,
537  Number * gradients_quad,
538  Number * scratch_data,
539  const bool integrate_values,
540  const bool integrate_gradients,
541  const bool sum_into_values_array = false);
542 };
543 
554 template <int dim, int n_q_points_1d, int n_components, typename Number>
555 struct SelectEvaluator<dim, -1, n_q_points_1d, n_components, Number>
556 {
565  static void
567  Number * values_dofs_actual,
568  Number * values_quad,
569  Number * gradients_quad,
570  Number * hessians_quad,
571  Number * scratch_data,
572  const bool evaluate_values,
573  const bool evaluate_gradients,
574  const bool evaluate_hessians);
575 
584  static void
586  Number * values_dofs_actual,
587  Number * values_quad,
588  Number * gradients_quad,
589  Number * scratch_data,
590  const bool integrate_values,
591  const bool integrate_gradients,
592  const bool sum_into_values_array = false);
593 };
594 
595 //----------------------Implementation for SelectEvaluator---------------------
596 #ifndef DOXYGEN
597 
598 template <int dim,
599  int fe_degree,
600  int n_q_points_1d,
601  int n_components,
602  typename Number>
603 inline void
606  Number * values_dofs_actual,
607  Number * values_quad,
608  Number * gradients_quad,
609  Number * hessians_quad,
610  Number * scratch_data,
611  const bool evaluate_values,
612  const bool evaluate_gradients,
613  const bool evaluate_hessians)
614 {
615  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
616 
617  if (fe_degree + 1 == n_q_points_1d &&
618  shape_info.element_type ==
620  {
621  internal::
623  evaluate(shape_info,
624  values_dofs_actual,
625  values_quad,
626  gradients_quad,
627  hessians_quad,
628  scratch_data,
629  evaluate_values,
630  evaluate_gradients,
631  evaluate_hessians);
632  }
633  // '<=' on type means tensor_symmetric or tensor_symmetric_hermite, see
634  // shape_info.h for more details
635  else if (use_collocation && shape_info.element_type <=
637  {
639  dim,
640  fe_degree,
641  n_q_points_1d,
642  n_components,
643  Number>::evaluate(shape_info,
644  values_dofs_actual,
645  values_quad,
646  gradients_quad,
647  hessians_quad,
648  scratch_data,
649  evaluate_values,
650  evaluate_gradients,
651  evaluate_hessians);
652  }
653  else if (shape_info.element_type <=
655  {
658  dim,
659  fe_degree,
660  n_q_points_1d,
661  n_components,
662  Number>::evaluate(shape_info,
663  values_dofs_actual,
664  values_quad,
665  gradients_quad,
666  hessians_quad,
667  scratch_data,
668  evaluate_values,
669  evaluate_gradients,
670  evaluate_hessians);
671  }
672  else if (shape_info.element_type ==
674  {
677  dim,
678  fe_degree,
679  n_q_points_1d,
680  n_components,
681  Number>::evaluate(shape_info,
682  values_dofs_actual,
683  values_quad,
684  gradients_quad,
685  hessians_quad,
686  scratch_data,
687  evaluate_values,
688  evaluate_gradients,
689  evaluate_hessians);
690  }
691  else if (shape_info.element_type ==
693  {
696  dim,
697  fe_degree,
698  n_q_points_1d,
699  n_components,
700  Number>::evaluate(shape_info,
701  values_dofs_actual,
702  values_quad,
703  gradients_quad,
704  hessians_quad,
705  scratch_data,
706  evaluate_values,
707  evaluate_gradients,
708  evaluate_hessians);
709  }
710  else if (shape_info.element_type ==
712  {
714  dim,
715  fe_degree,
716  n_q_points_1d,
717  n_components,
718  Number>::evaluate(shape_info,
719  values_dofs_actual,
720  values_quad,
721  gradients_quad,
722  hessians_quad,
723  scratch_data,
724  evaluate_values,
725  evaluate_gradients,
726  evaluate_hessians);
727  }
728  else
729  AssertThrow(false, ExcNotImplemented());
730 }
731 
732 
733 
734 template <int dim,
735  int fe_degree,
736  int n_q_points_1d,
737  int n_components,
738  typename Number>
739 inline void
742  Number * values_dofs_actual,
743  Number * values_quad,
744  Number * gradients_quad,
745  Number * scratch_data,
746  const bool integrate_values,
747  const bool integrate_gradients,
748  const bool sum_into_values_array)
749 {
750  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
751 
752  if (fe_degree + 1 == n_q_points_1d &&
753  shape_info.element_type ==
755  {
756  internal::
758  integrate(shape_info,
759  values_dofs_actual,
760  values_quad,
761  gradients_quad,
762  scratch_data,
763  integrate_values,
764  integrate_gradients,
765  sum_into_values_array);
766  }
767  // '<=' on type means tensor_symmetric or tensor_symmetric_hermite, see
768  // shape_info.h for more details
769  else if (use_collocation && shape_info.element_type <=
771  {
773  dim,
774  fe_degree,
775  n_q_points_1d,
776  n_components,
777  Number>::integrate(shape_info,
778  values_dofs_actual,
779  values_quad,
780  gradients_quad,
781  scratch_data,
782  integrate_values,
783  integrate_gradients,
784  sum_into_values_array);
785  }
786  else if (shape_info.element_type <=
788  {
791  dim,
792  fe_degree,
793  n_q_points_1d,
794  n_components,
795  Number>::integrate(shape_info,
796  values_dofs_actual,
797  values_quad,
798  gradients_quad,
799  scratch_data,
800  integrate_values,
801  integrate_gradients,
802  sum_into_values_array);
803  }
804  else if (shape_info.element_type ==
806  {
809  dim,
810  fe_degree,
811  n_q_points_1d,
812  n_components,
813  Number>::integrate(shape_info,
814  values_dofs_actual,
815  values_quad,
816  gradients_quad,
817  scratch_data,
818  integrate_values,
819  integrate_gradients,
820  sum_into_values_array);
821  }
822  else if (shape_info.element_type ==
824  {
827  dim,
828  fe_degree,
829  n_q_points_1d,
830  n_components,
831  Number>::integrate(shape_info,
832  values_dofs_actual,
833  values_quad,
834  gradients_quad,
835  scratch_data,
836  integrate_values,
837  integrate_gradients,
838  sum_into_values_array);
839  }
840  else if (shape_info.element_type ==
842  {
844  dim,
845  fe_degree,
846  n_q_points_1d,
847  n_components,
848  Number>::integrate(shape_info,
849  values_dofs_actual,
850  values_quad,
851  gradients_quad,
852  scratch_data,
853  integrate_values,
854  integrate_gradients,
855  sum_into_values_array);
856  }
857  else
858  AssertThrow(false, ExcNotImplemented());
859 }
860 
861 
862 
863 template <int dim, int dummy, int n_components, typename Number>
864 inline void
867  Number * values_dofs_actual,
868  Number * values_quad,
869  Number * gradients_quad,
870  Number * hessians_quad,
871  Number * scratch_data,
872  const bool evaluate_values,
873  const bool evaluate_gradients,
874  const bool evaluate_hessians)
875 {
876  if (shape_info.element_type ==
878  {
881  dim,
882  -1,
883  0,
884  n_components,
885  Number>::evaluate(shape_info,
886  values_dofs_actual,
887  values_quad,
888  gradients_quad,
889  hessians_quad,
890  scratch_data,
891  evaluate_values,
892  evaluate_gradients,
893  evaluate_hessians);
894  }
895  else if (shape_info.element_type ==
897  {
900  dim,
901  -1,
902  0,
903  n_components,
904  Number>::evaluate(shape_info,
905  values_dofs_actual,
906  values_quad,
907  gradients_quad,
908  hessians_quad,
909  scratch_data,
910  evaluate_values,
911  evaluate_gradients,
912  evaluate_hessians);
913  }
914  else if (shape_info.element_type ==
917  dim,
918  -1,
919  0,
920  n_components,
921  Number>::evaluate(shape_info,
922  values_dofs_actual,
923  values_quad,
924  gradients_quad,
925  hessians_quad,
926  scratch_data,
927  evaluate_values,
928  evaluate_gradients,
929  evaluate_hessians);
930  else
931  internal::EvaluationSelectorImplementation::
932  symmetric_selector_evaluate<dim, n_components, Number>(shape_info,
933  values_dofs_actual,
934  values_quad,
935  gradients_quad,
936  hessians_quad,
937  scratch_data,
938  evaluate_values,
939  evaluate_gradients,
940  evaluate_hessians);
941 }
942 
943 
944 
945 template <int dim, int dummy, int n_components, typename Number>
946 inline void
949  Number * values_dofs_actual,
950  Number * values_quad,
951  Number * gradients_quad,
952  Number * scratch_data,
953  const bool integrate_values,
954  const bool integrate_gradients,
955  const bool sum_into_values_array)
956 {
957  if (shape_info.element_type ==
959  {
962  dim,
963  -1,
964  0,
965  n_components,
966  Number>::integrate(shape_info,
967  values_dofs_actual,
968  values_quad,
969  gradients_quad,
970  scratch_data,
971  integrate_values,
972  integrate_gradients,
973  sum_into_values_array);
974  }
975  else if (shape_info.element_type ==
977  {
980  dim,
981  -1,
982  0,
983  n_components,
984  Number>::integrate(shape_info,
985  values_dofs_actual,
986  values_quad,
987  gradients_quad,
988  scratch_data,
989  integrate_values,
990  integrate_gradients,
991  sum_into_values_array);
992  }
993  else if (shape_info.element_type ==
996  dim,
997  -1,
998  0,
999  n_components,
1000  Number>::integrate(shape_info,
1001  values_dofs_actual,
1002  values_quad,
1003  gradients_quad,
1004  scratch_data,
1005  integrate_values,
1006  integrate_gradients,
1007  sum_into_values_array);
1008  else
1009  internal::EvaluationSelectorImplementation::
1010  symmetric_selector_integrate<dim, n_components, Number>(
1011  shape_info,
1012  values_dofs_actual,
1013  values_quad,
1014  gradients_quad,
1015  scratch_data,
1016  integrate_values,
1017  integrate_gradients,
1018  sum_into_values_array);
1019 }
1020 #endif // DOXYGEN
1021 
1023 
1024 #endif
SelectEvaluator::use_collocation
static constexpr bool use_collocation
Definition: evaluation_selector.h:504
internal::MatrixFreeFunctions::ShapeInfo::data
std::vector< UnivariateShapeData< Number > > data
Definition: shape_info.h:355
internal::MatrixFreeFunctions::ShapeInfo::element_type
ElementType element_type
Definition: shape_info.h:297
SelectEvaluator
Definition: evaluation_selector.h:495
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
internal::MatrixFreeFunctions::ShapeInfo
Definition: shape_info.h:290
internal::MatrixFreeFunctions::tensor_symmetric
@ tensor_symmetric
Definition: shape_info.h:76
SelectEvaluator::integrate
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)
DoFHandler::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
internal::MatrixFreeFunctions::truncated_tensor
@ truncated_tensor
Definition: shape_info.h:87
DoFTools::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
SelectEvaluator::evaluate
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)
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
internal::FEEvaluationImplTransformToCollocation
Definition: evaluation_kernels.h:1194
internal::MatrixFreeFunctions::tensor_general
@ tensor_general
Definition: shape_info.h:81
evaluation_kernels.h
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
internal::MatrixFreeFunctions::tensor_symmetric_collocation
@ tensor_symmetric_collocation
Definition: shape_info.h:62
internal::FEEvaluationImpl
Definition: evaluation_kernels.h:111
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
internal::FEEvaluationImplCollocation::evaluate
static void evaluate(const MatrixFreeFunctions::ShapeInfo< Number > &shape_info, const Number *values_dofs, Number *values_quad, Number *gradients_quad, Number *hessians_quad, Number *scratch_data, const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians)
Definition: evaluation_kernels.h:1042
config.h
internal
Definition: aligned_vector.h:369
internal::FEEvaluationImplCollocation::integrate
static void integrate(const MatrixFreeFunctions::ShapeInfo< Number > &shape_info, Number *values_dofs, Number *values_quad, Number *gradients_quad, Number *scratch_data, const bool integrate_values, const bool integrate_gradients, const bool add_into_values_array)
Definition: evaluation_kernels.h:1120
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0
@ tensor_symmetric_plus_dg0
Definition: shape_info.h:94
Factory
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531