Reference documentation for deal.II version 9.1.0-pre
evaluation_selector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2018 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/matrix_free/evaluation_kernels.h>
21 
22 DEAL_II_NAMESPACE_OPEN
23 
24 #ifndef DOXYGEN
25 namespace internal
26 {
27  namespace EvaluationSelectorImplementation
28  {
29  // The following classes serve the purpose of choosing the correct template
30  // specialization of the FEEvaluationImpl* classes in case fe_degree
31  // and n_q_points_1d are only given as runtime parameters.
32  // The logic is the following:
33  // 1. Start with fe_degree=0, n_q_points_1d=0 and DEPTH=0.
34  // 2. If the current assumption on fe_degree doesn't match the runtime
35  // parameter, increase fe_degree by one and try again.
36  // If fe_degree==10 use the class Default which serves as a fallback.
37  // 3. After fixing the fe_degree, DEPTH is increased (DEPTH=1) and we start
38  // with
39  // n_q_points=fe_degree+1.
40  // 4. If the current assumption on n_q_points_1d doesn't match the runtime
41  // parameter, increase n_q_points_1d by one and try again.
42  // If n_q_points_1d==degree+3 use the class Default which serves as a
43  // fallback.
44 
49  template <int dim, int n_components, typename Number>
50  struct Default
51  {
52  static inline void
53  evaluate(
55  Number * values_dofs_actual,
56  Number * values_quad,
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)
63  {
65  internal::MatrixFreeFunctions::tensor_general,
66  dim,
67  -1,
68  0,
70  Number>::evaluate(shape_info,
71  values_dofs_actual,
72  values_quad,
73  gradients_quad,
74  hessians_quad,
75  scratch_data,
76  evaluate_values,
77  evaluate_gradients,
78  evaluate_hessians);
79  }
80 
81  static inline void
82  integrate(
84  Number * values_dofs_actual,
85  Number * values_quad,
86  Number * gradients_quad,
87  Number * scratch_data,
88  const bool integrate_values,
89  const bool integrate_gradients)
90  {
92  internal::MatrixFreeFunctions::tensor_general,
93  dim,
94  -1,
95  0,
97  Number>::integrate(shape_info,
98  values_dofs_actual,
99  values_quad,
100  gradients_quad,
101  scratch_data,
102  integrate_values,
103  integrate_gradients,
104  false);
105  }
106  };
107 
108 
112  template <int dim,
113  int n_components,
114  typename Number,
115  int DEPTH = 0,
116  int degree = 0,
117  int n_q_points_1d = 0,
118  class Enable = void>
119  struct Factory : Default<dim, n_components, Number>
120  {};
121 
127  template <int n_q_points_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>
130  {};
131 
137  template <int degree,
138  int n_q_points_1d,
139  int dim,
140  int n_components,
141  typename Number>
142  struct Factory<dim,
143  n_components,
144  Number,
145  1,
146  degree,
147  n_q_points_1d,
148  typename std::enable_if<n_q_points_1d == degree + 3>::type>
149  : Default<dim, n_components, Number>
150  {};
151 
155  template <int degree,
156  int n_q_points_1d,
157  int dim,
158  int n_components,
159  typename Number>
160  struct Factory<dim, n_components, Number, 0, degree, n_q_points_1d>
161  {
162  static inline void
163  evaluate(
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)
173  {
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>::
178  evaluate(shape_info,
179  values_dofs_actual,
180  values_quad,
181  gradients_quad,
182  hessians_quad,
183  scratch_data,
184  evaluate_values,
185  evaluate_gradients,
186  evaluate_hessians);
187  else
188  Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
189  evaluate(shape_info,
190  values_dofs_actual,
191  values_quad,
192  gradients_quad,
193  hessians_quad,
194  scratch_data,
195  evaluate_values,
196  evaluate_gradients,
197  evaluate_hessians);
198  }
199 
200  static inline void
201  integrate(
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)
209  {
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,
215  values_dofs_actual,
216  values_quad,
217  gradients_quad,
218  scratch_data,
219  integrate_values,
220  integrate_gradients);
221  else
222  Factory<dim, n_components, Number, 0, degree + 1, n_q_points_1d>::
223  integrate(shape_info,
224  values_dofs_actual,
225  values_quad,
226  gradients_quad,
227  scratch_data,
228  integrate_values,
229  integrate_gradients);
230  }
231  };
232 
237  template <int degree,
238  int n_q_points_1d,
239  int dim,
240  int n_components,
241  typename Number>
242  struct Factory<dim,
243  n_components,
244  Number,
245  1,
246  degree,
247  n_q_points_1d,
248  typename std::enable_if<(n_q_points_1d < degree + 3)>::type>
249  {
250  static inline void
251  evaluate(
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)
261  {
262  const int runtime_n_q_points_1d = shape_info.n_q_points_1d;
263  if (runtime_n_q_points_1d == n_q_points_1d)
264  {
265  if (n_q_points_1d == degree + 1 &&
266  shape_info.element_type ==
267  internal::MatrixFreeFunctions::tensor_symmetric_collocation)
268  internal::
270  evaluate(shape_info,
271  values_dofs_actual,
272  values_quad,
273  gradients_quad,
274  hessians_quad,
275  scratch_data,
276  evaluate_values,
277  evaluate_gradients,
278  evaluate_hessians);
279  else if (degree < n_q_points_1d)
281  dim,
282  degree,
283  n_q_points_1d,
284  n_components,
285  Number>::evaluate(shape_info,
286  values_dofs_actual,
287  values_quad,
288  gradients_quad,
289  hessians_quad,
290  scratch_data,
291  evaluate_values,
292  evaluate_gradients,
293  evaluate_hessians);
294  else
296  internal::MatrixFreeFunctions::tensor_symmetric,
297  dim,
298  degree,
299  n_q_points_1d,
300  n_components,
301  Number>::evaluate(shape_info,
302  values_dofs_actual,
303  values_quad,
304  gradients_quad,
305  hessians_quad,
306  scratch_data,
307  evaluate_values,
308  evaluate_gradients,
309  evaluate_hessians);
310  }
311  else
312  Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
313  evaluate(shape_info,
314  values_dofs_actual,
315  values_quad,
316  gradients_quad,
317  hessians_quad,
318  scratch_data,
319  evaluate_values,
320  evaluate_gradients,
321  evaluate_hessians);
322  }
323 
324  static inline void
325  integrate(
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)
333  {
334  const int runtime_n_q_points_1d = shape_info.n_q_points_1d;
335  if (runtime_n_q_points_1d == n_q_points_1d)
336  {
337  if (n_q_points_1d == degree + 1 &&
338  shape_info.element_type ==
339  internal::MatrixFreeFunctions::tensor_symmetric_collocation)
340  internal::
342  integrate(shape_info,
343  values_dofs_actual,
344  values_quad,
345  gradients_quad,
346  scratch_data,
347  integrate_values,
348  integrate_gradients,
349  false);
350  else if (degree < n_q_points_1d)
352  dim,
353  degree,
354  n_q_points_1d,
355  n_components,
356  Number>::integrate(shape_info,
357  values_dofs_actual,
358  values_quad,
359  gradients_quad,
360  scratch_data,
361  integrate_values,
362  integrate_gradients,
363  false);
364  else
366  internal::MatrixFreeFunctions::tensor_symmetric,
367  dim,
368  degree,
369  n_q_points_1d,
370  n_components,
371  Number>::integrate(shape_info,
372  values_dofs_actual,
373  values_quad,
374  gradients_quad,
375  scratch_data,
376  integrate_values,
377  integrate_gradients,
378  false);
379  }
380  else
381  Factory<dim, n_components, Number, 1, degree, n_q_points_1d + 1>::
382  integrate(shape_info,
383  values_dofs_actual,
384  values_quad,
385  gradients_quad,
386  scratch_data,
387  integrate_values,
388  integrate_gradients);
389  }
390  };
391 
392 
393 
398  template <int dim, int n_components, typename Number>
399  void
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)
410  {
411  Assert(shape_info.element_type <=
412  internal::MatrixFreeFunctions::tensor_symmetric,
413  ExcInternalError());
414  Factory<dim, n_components, Number>::evaluate(shape_info,
415  values_dofs_actual,
416  values_quad,
417  gradients_quad,
418  hessians_quad,
419  scratch_data,
420  evaluate_values,
421  evaluate_gradients,
422  evaluate_hessians);
423  }
424 
425 
426 
431  template <int dim, int n_components, typename Number>
432  void
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)
441  {
442  Assert(shape_info.element_type <=
443  internal::MatrixFreeFunctions::tensor_symmetric,
444  ExcInternalError());
445  Factory<dim, n_components, Number>::integrate(shape_info,
446  values_dofs_actual,
447  values_quad,
448  gradients_quad,
449  scratch_data,
450  integrate_values,
451  integrate_gradients);
452  }
453  } // namespace EvaluationSelectorImplementation
454 } // namespace internal
455 #endif
456 
457 
469 template <int dim,
470  int fe_degree,
471  int n_q_points_1d,
472  int n_components,
473  typename Number>
475 {
483  static void
484  evaluate(const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
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);
493 
501  static void
502  integrate(const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
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);
509 };
510 
521 template <int dim, int n_q_points_1d, int n_components, typename Number>
522 struct SelectEvaluator<dim, -1, n_q_points_1d, n_components, Number>
523 {
532  static void
533  evaluate(const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
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);
542 
551  static void
552  integrate(const internal::MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
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);
559 };
560 
561 //----------------------Implementation for SelectEvaluator---------------------
562 #ifndef DOXYGEN
563 
564 template <int dim,
565  int fe_degree,
566  int n_q_points_1d,
567  int n_components,
568  typename Number>
569 inline void
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)
580 {
581  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
582 
583  if (fe_degree + 1 == n_q_points_1d &&
584  shape_info.element_type ==
585  internal::MatrixFreeFunctions::tensor_symmetric_collocation)
586  {
587  internal::
589  evaluate(shape_info,
590  values_dofs_actual,
591  values_quad,
592  gradients_quad,
593  hessians_quad,
594  scratch_data,
595  evaluate_values,
596  evaluate_gradients,
597  evaluate_hessians);
598  }
599  else if (fe_degree < n_q_points_1d &&
600  shape_info.element_type <=
601  internal::MatrixFreeFunctions::tensor_symmetric)
602  {
604  dim,
605  fe_degree,
606  n_q_points_1d,
607  n_components,
608  Number>::evaluate(shape_info,
609  values_dofs_actual,
610  values_quad,
611  gradients_quad,
612  hessians_quad,
613  scratch_data,
614  evaluate_values,
615  evaluate_gradients,
616  evaluate_hessians);
617  }
618  else if (shape_info.element_type ==
619  internal::MatrixFreeFunctions::tensor_symmetric)
620  {
622  internal::MatrixFreeFunctions::tensor_symmetric,
623  dim,
624  fe_degree,
625  n_q_points_1d,
626  n_components,
627  Number>::evaluate(shape_info,
628  values_dofs_actual,
629  values_quad,
630  gradients_quad,
631  hessians_quad,
632  scratch_data,
633  evaluate_values,
634  evaluate_gradients,
635  evaluate_hessians);
636  }
637  else if (shape_info.element_type ==
638  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
639  {
641  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
642  dim,
643  fe_degree,
644  n_q_points_1d,
645  n_components,
646  Number>::evaluate(shape_info,
647  values_dofs_actual,
648  values_quad,
649  gradients_quad,
650  hessians_quad,
651  scratch_data,
652  evaluate_values,
653  evaluate_gradients,
654  evaluate_hessians);
655  }
656  else if (shape_info.element_type ==
657  internal::MatrixFreeFunctions::truncated_tensor)
658  {
660  internal::MatrixFreeFunctions::truncated_tensor,
661  dim,
662  fe_degree,
663  n_q_points_1d,
664  n_components,
665  Number>::evaluate(shape_info,
666  values_dofs_actual,
667  values_quad,
668  gradients_quad,
669  hessians_quad,
670  scratch_data,
671  evaluate_values,
672  evaluate_gradients,
673  evaluate_hessians);
674  }
675  else if (shape_info.element_type ==
676  internal::MatrixFreeFunctions::tensor_general)
677  {
678  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
679  dim,
680  fe_degree,
681  n_q_points_1d,
682  n_components,
683  Number>::evaluate(shape_info,
684  values_dofs_actual,
685  values_quad,
686  gradients_quad,
687  hessians_quad,
688  scratch_data,
689  evaluate_values,
690  evaluate_gradients,
691  evaluate_hessians);
692  }
693  else
694  AssertThrow(false, ExcNotImplemented());
695 }
696 
697 
698 
699 template <int dim,
700  int fe_degree,
701  int n_q_points_1d,
702  int n_components,
703  typename Number>
704 inline void
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)
713 {
714  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
715 
716  if (fe_degree + 1 == n_q_points_1d &&
717  shape_info.element_type ==
718  internal::MatrixFreeFunctions::tensor_symmetric_collocation)
719  {
720  internal::
722  integrate(shape_info,
723  values_dofs_actual,
724  values_quad,
725  gradients_quad,
726  scratch_data,
727  integrate_values,
728  integrate_gradients,
729  false);
730  }
731  else if (fe_degree < n_q_points_1d &&
732  shape_info.element_type <=
733  internal::MatrixFreeFunctions::tensor_symmetric)
734  {
736  dim,
737  fe_degree,
738  n_q_points_1d,
739  n_components,
740  Number>::integrate(shape_info,
741  values_dofs_actual,
742  values_quad,
743  gradients_quad,
744  scratch_data,
745  integrate_values,
746  integrate_gradients,
747  false);
748  }
749  else if (shape_info.element_type ==
750  internal::MatrixFreeFunctions::tensor_symmetric)
751  {
753  internal::MatrixFreeFunctions::tensor_symmetric,
754  dim,
755  fe_degree,
756  n_q_points_1d,
757  n_components,
758  Number>::integrate(shape_info,
759  values_dofs_actual,
760  values_quad,
761  gradients_quad,
762  scratch_data,
763  integrate_values,
764  integrate_gradients,
765  false);
766  }
767  else if (shape_info.element_type ==
768  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
769  {
771  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
772  dim,
773  fe_degree,
774  n_q_points_1d,
775  n_components,
776  Number>::integrate(shape_info,
777  values_dofs_actual,
778  values_quad,
779  gradients_quad,
780  scratch_data,
781  integrate_values,
782  integrate_gradients,
783  false);
784  }
785  else if (shape_info.element_type ==
786  internal::MatrixFreeFunctions::truncated_tensor)
787  {
789  internal::MatrixFreeFunctions::truncated_tensor,
790  dim,
791  fe_degree,
792  n_q_points_1d,
793  n_components,
794  Number>::integrate(shape_info,
795  values_dofs_actual,
796  values_quad,
797  gradients_quad,
798  scratch_data,
799  integrate_values,
800  integrate_gradients,
801  false);
802  }
803  else if (shape_info.element_type ==
804  internal::MatrixFreeFunctions::tensor_general)
805  {
806  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
807  dim,
808  fe_degree,
809  n_q_points_1d,
810  n_components,
811  Number>::integrate(shape_info,
812  values_dofs_actual,
813  values_quad,
814  gradients_quad,
815  scratch_data,
816  integrate_values,
817  integrate_gradients,
818  false);
819  }
820  else
821  AssertThrow(false, ExcNotImplemented());
822 }
823 
824 
825 
826 template <int dim, int dummy, int n_components, typename Number>
827 inline void
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)
838 {
839  if (shape_info.element_type ==
840  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
841  {
843  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
844  dim,
845  -1,
846  0,
847  n_components,
848  Number>::evaluate(shape_info,
849  values_dofs_actual,
850  values_quad,
851  gradients_quad,
852  hessians_quad,
853  scratch_data,
854  evaluate_values,
855  evaluate_gradients,
856  evaluate_hessians);
857  }
858  else if (shape_info.element_type ==
859  internal::MatrixFreeFunctions::truncated_tensor)
860  {
862  internal::MatrixFreeFunctions::truncated_tensor,
863  dim,
864  -1,
865  0,
866  n_components,
867  Number>::evaluate(shape_info,
868  values_dofs_actual,
869  values_quad,
870  gradients_quad,
871  hessians_quad,
872  scratch_data,
873  evaluate_values,
874  evaluate_gradients,
875  evaluate_hessians);
876  }
877  else if (shape_info.element_type ==
878  internal::MatrixFreeFunctions::tensor_general)
879  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
880  dim,
881  -1,
882  0,
883  n_components,
884  Number>::evaluate(shape_info,
885  values_dofs_actual,
886  values_quad,
887  gradients_quad,
888  hessians_quad,
889  scratch_data,
890  evaluate_values,
891  evaluate_gradients,
892  evaluate_hessians);
893  else
894  internal::EvaluationSelectorImplementation::
895  symmetric_selector_evaluate<dim, n_components, Number>(shape_info,
896  values_dofs_actual,
897  values_quad,
898  gradients_quad,
899  hessians_quad,
900  scratch_data,
901  evaluate_values,
902  evaluate_gradients,
903  evaluate_hessians);
904 }
905 
906 
907 
908 template <int dim, int dummy, int n_components, typename Number>
909 inline void
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)
918 {
919  if (shape_info.element_type ==
920  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
921  {
923  internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
924  dim,
925  -1,
926  0,
927  n_components,
928  Number>::integrate(shape_info,
929  values_dofs_actual,
930  values_quad,
931  gradients_quad,
932  scratch_data,
933  integrate_values,
934  integrate_gradients,
935  false);
936  }
937  else if (shape_info.element_type ==
938  internal::MatrixFreeFunctions::truncated_tensor)
939  {
941  internal::MatrixFreeFunctions::truncated_tensor,
942  dim,
943  -1,
944  0,
945  n_components,
946  Number>::integrate(shape_info,
947  values_dofs_actual,
948  values_quad,
949  gradients_quad,
950  scratch_data,
951  integrate_values,
952  integrate_gradients,
953  false);
954  }
955  else if (shape_info.element_type ==
956  internal::MatrixFreeFunctions::tensor_general)
957  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
958  dim,
959  -1,
960  0,
961  n_components,
962  Number>::integrate(shape_info,
963  values_dofs_actual,
964  values_quad,
965  gradients_quad,
966  scratch_data,
967  integrate_values,
968  integrate_gradients,
969  false);
970  else
971  internal::EvaluationSelectorImplementation::
972  symmetric_selector_integrate<dim, n_components, Number>(
973  shape_info,
974  values_dofs_actual,
975  values_quad,
976  gradients_quad,
977  scratch_data,
978  integrate_values,
979  integrate_gradients);
980 }
981 #endif // DOXYGEN
982 
983 DEAL_II_NAMESPACE_CLOSE
984 
985 #endif
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)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
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)
Definition: exceptions.h:1227
static::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static::ExceptionBase & ExcInternalError()