16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/multithread_info.h> 18 #include <deal.II/base/quadrature.h> 19 #include <deal.II/base/signaling_nan.h> 20 #include <deal.II/base/std_cxx11/unique_ptr.h> 21 #include <deal.II/lac/vector.h> 22 #include <deal.II/lac/block_vector.h> 23 #include <deal.II/lac/la_vector.h> 24 #include <deal.II/lac/la_parallel_vector.h> 25 #include <deal.II/lac/la_parallel_block_vector.h> 26 #include <deal.II/lac/petsc_vector.h> 27 #include <deal.II/lac/petsc_block_vector.h> 28 #include <deal.II/lac/trilinos_vector.h> 29 #include <deal.II/lac/trilinos_block_vector.h> 30 #include <deal.II/grid/tria_iterator.h> 31 #include <deal.II/grid/tria_accessor.h> 32 #include <deal.II/grid/tria_boundary.h> 33 #include <deal.II/dofs/dof_accessor.h> 34 #include <deal.II/fe/mapping_q1.h> 35 #include <deal.II/fe/fe_values.h> 36 #include <deal.II/fe/fe.h> 40 DEAL_II_NAMESPACE_OPEN
45 template <
class VectorType>
46 typename VectorType::value_type
47 get_vector_element (
const VectorType &vector,
50 return vector[cell_number];
55 get_vector_element (
const IndexSet &is,
65 template <
int dim,
int spacedim>
67 std::vector<unsigned int>
79 unsigned int nth_nonzero_component = 0;
83 shape_function_to_row_table[i*fe.
n_components()+c] = row + nth_nonzero_component;
84 ++nth_nonzero_component;
89 return shape_function_to_row_table;
97 template <
int dim,
int spacedim>
99 const unsigned int component)
101 fe_values (&fe_values),
102 component (component),
103 shape_function_data (this->fe_values->fe->dofs_per_cell)
112 const std::vector<unsigned int> shape_function_to_row_table
113 = make_shape_function_to_row_table (fe);
119 if (is_primitive ==
true)
138 template <
int dim,
int spacedim>
142 component (
numbers::invalid_unsigned_int)
146 template <
int dim,
int spacedim>
157 template <
int dim,
int spacedim>
159 const unsigned int first_vector_component)
161 fe_values (&fe_values),
162 first_vector_component (first_vector_component),
163 shape_function_data (this->fe_values->fe->dofs_per_cell)
173 const std::vector<unsigned int> shape_function_to_row_table
174 = make_shape_function_to_row_table (fe);
176 for (
unsigned int d=0; d<spacedim; ++d)
184 if (is_primitive ==
true)
196 = shape_function_to_row_table[i*fe.
n_components()+component];
205 unsigned int n_nonzero_components = 0;
206 for (
unsigned int d=0; d<spacedim; ++d)
209 ++n_nonzero_components;
211 if (n_nonzero_components == 0)
213 else if (n_nonzero_components > 1)
217 for (
unsigned int d=0; d<spacedim; ++d)
232 template <
int dim,
int spacedim>
236 first_vector_component (
numbers::invalid_unsigned_int)
241 template <
int dim,
int spacedim>
251 template <
int dim,
int spacedim>
254 const unsigned int first_tensor_component)
256 fe_values(&fe_values),
257 first_tensor_component(first_tensor_component),
258 shape_function_data(this->fe_values->fe->dofs_per_cell)
261 Assert(first_tensor_component + (dim*dim+dim)/2 - 1
271 const std::vector<unsigned int> shape_function_to_row_table
272 = make_shape_function_to_row_table (fe);
274 for (
unsigned int d = 0; d < ::SymmetricTensor<2,dim>::n_independent_components; ++d)
276 const unsigned int component = first_tensor_component + d;
282 if (is_primitive ==
true)
283 shape_function_data[i].is_nonzero_shape_function_component[d]
287 shape_function_data[i].is_nonzero_shape_function_component[d]
291 if (shape_function_data[i].is_nonzero_shape_function_component[d]
293 shape_function_data[i].row_index[d]
294 = shape_function_to_row_table[i*fe.
n_components()+component];
296 shape_function_data[i].row_index[d]
303 unsigned int n_nonzero_components = 0;
304 for (
unsigned int d = 0; d < ::SymmetricTensor<2,dim>::n_independent_components; ++d)
305 if (shape_function_data[i].is_nonzero_shape_function_component[d]
307 ++n_nonzero_components;
309 if (n_nonzero_components == 0)
310 shape_function_data[i].single_nonzero_component = -2;
311 else if (n_nonzero_components > 1)
312 shape_function_data[i].single_nonzero_component = -1;
315 for (
unsigned int d = 0; d < ::SymmetricTensor<2,dim>::n_independent_components; ++d)
316 if (shape_function_data[i].is_nonzero_shape_function_component[d]
319 shape_function_data[i].single_nonzero_component
320 = shape_function_data[i].row_index[d];
321 shape_function_data[i].single_nonzero_component_index
331 template <
int dim,
int spacedim>
335 first_tensor_component(
numbers::invalid_unsigned_int)
340 template <
int dim,
int spacedim>
350 template <
int dim,
int spacedim>
353 const unsigned int first_tensor_component)
355 fe_values(&fe_values),
356 first_tensor_component(first_tensor_component),
357 shape_function_data(this->fe_values->fe->dofs_per_cell)
360 Assert(first_tensor_component + dim*dim - 1
370 const std::vector<unsigned int> shape_function_to_row_table
371 = make_shape_function_to_row_table (fe);
373 for (
unsigned int d = 0; d < dim*dim; ++d)
375 const unsigned int component = first_tensor_component + d;
381 if (is_primitive ==
true)
382 shape_function_data[i].is_nonzero_shape_function_component[d]
386 shape_function_data[i].is_nonzero_shape_function_component[d]
390 if (shape_function_data[i].is_nonzero_shape_function_component[d]
392 shape_function_data[i].row_index[d]
393 = shape_function_to_row_table[i*fe.
n_components()+component];
395 shape_function_data[i].row_index[d]
402 unsigned int n_nonzero_components = 0;
403 for (
unsigned int d = 0; d < dim*dim; ++d)
404 if (shape_function_data[i].is_nonzero_shape_function_component[d]
406 ++n_nonzero_components;
408 if (n_nonzero_components == 0)
409 shape_function_data[i].single_nonzero_component = -2;
410 else if (n_nonzero_components > 1)
411 shape_function_data[i].single_nonzero_component = -1;
414 for (
unsigned int d = 0; d < dim*dim; ++d)
415 if (shape_function_data[i].is_nonzero_shape_function_component[d]
418 shape_function_data[i].single_nonzero_component
419 = shape_function_data[i].row_index[d];
420 shape_function_data[i].single_nonzero_component_index
430 template <
int dim,
int spacedim>
434 first_tensor_component(
numbers::invalid_unsigned_int)
439 template <
int dim,
int spacedim>
455 template <
int dim,
int spacedim,
typename Number>
457 do_function_values (const ::Vector<Number> &dof_values,
462 const unsigned int dofs_per_cell = dof_values.size();
463 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
464 shape_values.n_cols() : values.
size();
467 std::fill (values.begin(), values.end(), Number());
469 for (
unsigned int shape_function=0;
470 shape_function<dofs_per_cell; ++shape_function)
471 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
473 const Number value = dof_values(shape_function);
474 if (value == Number() )
477 const double *shape_value_ptr =
478 &shape_values(shape_function_data[shape_function].row_index, 0);
479 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
480 values[q_point] += value **shape_value_ptr++;
488 template <
int order,
int dim,
int spacedim,
typename Number>
490 do_function_derivatives (const ::Vector<Number> &dof_values,
492 const std::vector<
typename Scalar<dim,spacedim>::ShapeFunctionData> &shape_function_data,
495 const unsigned int dofs_per_cell = dof_values.size();
496 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
497 shape_derivatives[0].size() : derivatives.size();
500 std::fill (derivatives.begin(), derivatives.end(),
503 for (
unsigned int shape_function=0;
504 shape_function<dofs_per_cell; ++shape_function)
505 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
507 const Number value = dof_values(shape_function);
508 if (value == Number() )
511 const ::Tensor<order,spacedim> *shape_derivative_ptr =
512 &shape_derivatives[shape_function_data[shape_function].row_index][0];
513 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
514 derivatives[q_point] += value *
521 template <
int dim,
int spacedim,
typename Number>
523 do_function_laplacians (const ::Vector<Number> &dof_values,
525 const std::vector<
typename Scalar<dim,spacedim>::ShapeFunctionData> &shape_function_data,
528 const unsigned int dofs_per_cell = dof_values.size();
529 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
530 shape_hessians[0].size() : laplacians.size();
535 for (
unsigned int shape_function=0;
536 shape_function<dofs_per_cell; ++shape_function)
537 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
539 const Number value = dof_values(shape_function);
540 if (value == Number())
543 const ::Tensor<2,spacedim> *shape_hessian_ptr =
544 &shape_hessians[shape_function_data[shape_function].row_index][0];
545 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
546 laplacians[q_point] += value * trace(*shape_hessian_ptr++);
554 template <
int dim,
int spacedim,
typename Number>
555 void do_function_values (const ::Vector<Number> &dof_values,
560 const unsigned int dofs_per_cell = dof_values.size();
561 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
562 shape_values.n_cols() : values.
size();
567 for (
unsigned int shape_function=0;
568 shape_function<dofs_per_cell; ++shape_function)
570 const int snc = shape_function_data[shape_function].single_nonzero_component;
576 const Number value = dof_values(shape_function);
577 if (value == Number())
582 const unsigned int comp =
583 shape_function_data[shape_function].single_nonzero_component_index;
584 const double *shape_value_ptr = &shape_values(snc,0);
585 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
586 values[q_point][comp] += value **shape_value_ptr++;
589 for (
unsigned int d=0;
d<spacedim; ++
d)
590 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
592 const double *shape_value_ptr =
593 &shape_values(shape_function_data[shape_function].row_index[d],0);
594 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
595 values[q_point][d] += value **shape_value_ptr++;
602 template <
int order,
int dim,
int spacedim,
typename Number>
604 do_function_derivatives (const ::Vector<Number> &dof_values,
609 const unsigned int dofs_per_cell = dof_values.size();
610 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
611 shape_derivatives[0].size() : derivatives.size();
614 std::fill (derivatives.begin(), derivatives.end(),
617 for (
unsigned int shape_function=0;
618 shape_function<dofs_per_cell; ++shape_function)
620 const int snc = shape_function_data[shape_function].single_nonzero_component;
626 const Number value = dof_values(shape_function);
627 if (value == Number())
632 const unsigned int comp =
633 shape_function_data[shape_function].single_nonzero_component_index;
634 const ::Tensor<order,spacedim> *shape_derivative_ptr =
635 &shape_derivatives[snc][0];
636 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
637 derivatives[q_point][comp] += value *
641 for (
unsigned int d=0;
d<spacedim; ++
d)
642 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
644 const ::Tensor<order,spacedim> *shape_derivative_ptr =
645 &shape_derivatives[shape_function_data[shape_function].
647 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
648 derivatives[q_point][d] += value *
656 template <
int dim,
int spacedim,
typename Number>
658 do_function_symmetric_gradients (const ::Vector<Number> &dof_values,
663 const unsigned int dofs_per_cell = dof_values.size();
664 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
665 shape_gradients[0].size() : symmetric_gradients.size();
668 std::fill (symmetric_gradients.begin(), symmetric_gradients.end(),
671 for (
unsigned int shape_function=0;
672 shape_function<dofs_per_cell; ++shape_function)
674 const int snc = shape_function_data[shape_function].single_nonzero_component;
680 const Number value = dof_values(shape_function);
681 if (value == Number())
686 const unsigned int comp =
687 shape_function_data[shape_function].single_nonzero_component_index;
688 const ::Tensor<1,spacedim> *shape_gradient_ptr =
689 &shape_gradients[snc][0];
690 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
691 symmetric_gradients[q_point] += value *
695 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
698 for (
unsigned int d=0;
d<spacedim; ++
d)
699 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
701 shape_gradients[shape_function_data[shape_function].row_index[
d]][q_point];
702 symmetric_gradients[q_point] += symmetrize(grad);
709 template <
int dim,
int spacedim,
typename Number>
711 do_function_divergences (const ::Vector<Number> &dof_values,
716 const unsigned int dofs_per_cell = dof_values.size();
717 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
718 shape_gradients[0].size() : divergences.size();
723 for (
unsigned int shape_function=0;
724 shape_function<dofs_per_cell; ++shape_function)
726 const int snc = shape_function_data[shape_function].single_nonzero_component;
732 const Number value = dof_values(shape_function);
733 if (value == Number())
738 const unsigned int comp =
739 shape_function_data[shape_function].single_nonzero_component_index;
740 const ::Tensor<1,spacedim> *shape_gradient_ptr = &shape_gradients[snc][0];
741 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
742 divergences[q_point] += value * (*shape_gradient_ptr++)[comp];
745 for (
unsigned int d=0;
d<spacedim; ++
d)
746 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
748 const ::Tensor<1,spacedim> *shape_gradient_ptr =
749 &shape_gradients[shape_function_data[shape_function].
751 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
752 divergences[q_point] += value * (*shape_gradient_ptr++)[
d];
759 template <
int dim,
int spacedim,
typename Number>
761 do_function_curls (const ::Vector<Number> &dof_values,
764 std::vector<
typename ProductType<Number,typename ::internal::CurlType<spacedim>::type>::type> &curls)
766 const unsigned int dofs_per_cell = dof_values.size();
767 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
768 shape_gradients[0].size() : curls.size();
777 Assert (
false,
ExcMessage(
"Computing the curl in 1d is not a useful operation"));
783 for (
unsigned int shape_function = 0;
784 shape_function < dofs_per_cell; ++shape_function)
786 const int snc = shape_function_data[shape_function].single_nonzero_component;
792 const Number value = dof_values (shape_function);
794 if (value == Number())
799 const ::Tensor<1, spacedim> *shape_gradient_ptr =
800 &shape_gradients[snc][0];
802 Assert (shape_function_data[shape_function].single_nonzero_component >= 0,
805 if (shape_function_data[shape_function].single_nonzero_component_index == 0)
806 for (
unsigned int q_point = 0;
807 q_point < n_quadrature_points; ++q_point)
808 curls[q_point][0] -= value * (*shape_gradient_ptr++)[1];
810 for (
unsigned int q_point = 0;
811 q_point < n_quadrature_points; ++q_point)
812 curls[q_point][0] += value * (*shape_gradient_ptr++)[0];
819 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
821 const ::Tensor<1,spacedim> *shape_gradient_ptr =
822 &shape_gradients[shape_function_data[shape_function].row_index[0]][0];
824 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
825 curls[q_point][0] -= value * (*shape_gradient_ptr++)[1];
828 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
830 const ::Tensor<1,spacedim> *shape_gradient_ptr =
831 &shape_gradients[shape_function_data[shape_function].row_index[1]][0];
833 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
834 curls[q_point][0] += value * (*shape_gradient_ptr++)[0];
843 for (
unsigned int shape_function = 0;
844 shape_function < dofs_per_cell; ++shape_function)
846 const int snc = shape_function_data[shape_function].single_nonzero_component;
852 const Number value = dof_values (shape_function);
854 if (value == Number())
859 const ::Tensor<1, spacedim> *shape_gradient_ptr = &shape_gradients[snc][0];
861 switch (shape_function_data[shape_function].single_nonzero_component_index)
865 for (
unsigned int q_point = 0;
866 q_point < n_quadrature_points; ++q_point)
868 curls[q_point][1] += value * (*shape_gradient_ptr)[2];
869 curls[q_point][2] -= value * (*shape_gradient_ptr++)[1];
877 for (
unsigned int q_point = 0;
878 q_point < n_quadrature_points; ++q_point)
880 curls[q_point][0] -= value * (*shape_gradient_ptr)[2];
881 curls[q_point][2] += value * (*shape_gradient_ptr++)[0];
889 for (
unsigned int q_point = 0;
890 q_point < n_quadrature_points; ++q_point)
892 curls[q_point][0] += value * (*shape_gradient_ptr)[1];
893 curls[q_point][1] -= value * (*shape_gradient_ptr++)[0];
908 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
910 const ::Tensor<1,spacedim> *shape_gradient_ptr =
911 &shape_gradients[shape_function_data[shape_function].row_index[0]][0];
913 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
915 curls[q_point][1] += value * (*shape_gradient_ptr)[2];
916 curls[q_point][2] -= value * (*shape_gradient_ptr++)[1];
920 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
922 const ::Tensor<1,spacedim> *shape_gradient_ptr =
923 &shape_gradients[shape_function_data[shape_function].row_index[1]][0];
925 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
927 curls[q_point][0] -= value * (*shape_gradient_ptr)[2];
928 curls[q_point][2] += value * (*shape_gradient_ptr++)[0];
932 if (shape_function_data[shape_function].is_nonzero_shape_function_component[2])
934 const ::Tensor<1,spacedim> *shape_gradient_ptr =
935 &shape_gradients[shape_function_data[shape_function].row_index[2]][0];
937 for (
unsigned int q_point = 0; q_point < n_quadrature_points; ++q_point)
939 curls[q_point][0] += value * (*shape_gradient_ptr)[1];
940 curls[q_point][1] -= value * (*shape_gradient_ptr++)[0];
951 template <
int dim,
int spacedim,
typename Number>
953 do_function_laplacians (const ::Vector<Number> &dof_values,
958 const unsigned int dofs_per_cell = dof_values.size();
959 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
960 shape_hessians[0].size() : laplacians.size();
963 std::fill (laplacians.begin(), laplacians.end(),
966 for (
unsigned int shape_function=0;
967 shape_function<dofs_per_cell; ++shape_function)
969 const int snc = shape_function_data[shape_function].single_nonzero_component;
975 const Number value = dof_values(shape_function);
976 if (value == Number())
981 const unsigned int comp =
982 shape_function_data[shape_function].single_nonzero_component_index;
983 const ::Tensor<2,spacedim> *shape_hessian_ptr =
984 &shape_hessians[snc][0];
985 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
986 laplacians[q_point][comp] += value * trace(*shape_hessian_ptr++);
989 for (
unsigned int d=0;
d<spacedim; ++
d)
990 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
992 const ::Tensor<2,spacedim> *shape_hessian_ptr =
993 &shape_hessians[shape_function_data[shape_function].
995 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
996 laplacians[q_point][d] += value * trace(*shape_hessian_ptr++);
1005 template <
int dim,
int spacedim,
typename Number>
1007 do_function_values (const ::Vector<Number> &dof_values,
1008 const ::Table<2,double> &shape_values,
1012 const unsigned int dofs_per_cell = dof_values.size();
1013 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
1014 shape_values.n_cols() : values.size();
1017 std::fill (values.begin(), values.end(),
1020 for (
unsigned int shape_function=0;
1021 shape_function<dofs_per_cell; ++shape_function)
1023 const int snc = shape_function_data[shape_function].single_nonzero_component;
1029 const Number value = dof_values(shape_function);
1030 if (value == Number())
1037 (shape_function_data[shape_function].single_nonzero_component_index);
1038 const double *shape_value_ptr = &shape_values(snc,0);
1039 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
1040 values[q_point][comp] += value **shape_value_ptr++;
1043 for (
unsigned int d=0;
1044 d<::SymmetricTensor<2,spacedim>::n_independent_components; ++
d)
1045 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
1049 const double *shape_value_ptr =
1050 &shape_values(shape_function_data[shape_function].row_index[d],0);
1051 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
1052 values[q_point][comp] += value **shape_value_ptr++;
1059 template <
int dim,
int spacedim,
typename Number>
1061 do_function_divergences (const ::Vector<Number> &dof_values,
1066 const unsigned int dofs_per_cell = dof_values.size();
1067 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
1068 shape_gradients[0].size() : divergences.size();
1071 std::fill (divergences.begin(), divergences.end(),
1074 for (
unsigned int shape_function=0;
1075 shape_function<dofs_per_cell; ++shape_function)
1077 const int snc = shape_function_data[shape_function].single_nonzero_component;
1083 const Number value = dof_values(shape_function);
1084 if (value == Number())
1089 const unsigned int comp =
1090 shape_function_data[shape_function].single_nonzero_component_index;
1092 const ::Tensor < 1, spacedim> *shape_gradient_ptr =
1093 &shape_gradients[snc][0];
1100 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1101 ++q_point, ++shape_gradient_ptr)
1103 divergences[q_point][ii] += value * (*shape_gradient_ptr)[jj];
1106 divergences[q_point][jj] += value * (*shape_gradient_ptr)[ii];
1111 for (
unsigned int d = 0;
1112 d < ::SymmetricTensor<2,spacedim>::n_independent_components; ++
d)
1113 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
1125 const unsigned int comp =
1126 shape_function_data[shape_function].single_nonzero_component_index;
1128 const ::Tensor < 1, spacedim> *shape_gradient_ptr =
1129 &shape_gradients[shape_function_data[shape_function].
1131 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1132 ++q_point, ++shape_gradient_ptr)
1134 for (
unsigned int j = 0; j < spacedim; ++j)
1137 divergences[q_point][vector_component] += value * (*shape_gradient_ptr++)[j];
1147 template <
int dim,
int spacedim,
typename Number>
1149 do_function_values (const ::Vector<Number> &dof_values,
1150 const ::Table<2,double> &shape_values,
1154 const unsigned int dofs_per_cell = dof_values.size();
1155 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
1156 shape_values.n_cols() : values.size();
1159 std::fill (values.begin(), values.end(),
1162 for (
unsigned int shape_function=0;
1163 shape_function<dofs_per_cell; ++shape_function)
1165 const int snc = shape_function_data[shape_function].single_nonzero_component;
1171 const Number value = dof_values(shape_function);
1172 if (value == Number())
1177 const unsigned int comp =
1178 shape_function_data[shape_function].single_nonzero_component_index;
1182 const double *shape_value_ptr = &shape_values(snc,0);
1183 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
1184 values[q_point][indices] += value **shape_value_ptr++;
1187 for (
unsigned int d=0;
1189 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
1193 const double *shape_value_ptr =
1194 &shape_values(shape_function_data[shape_function].row_index[d],0);
1195 for (
unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
1196 values[q_point][indices] += value **shape_value_ptr++;
1203 template <
int dim,
int spacedim,
typename Number>
1205 do_function_divergences (const ::Vector<Number> &dof_values,
1210 const unsigned int dofs_per_cell = dof_values.size();
1211 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
1212 shape_gradients[0].size() : divergences.size();
1215 std::fill (divergences.begin(), divergences.end(),
1218 for (
unsigned int shape_function=0;
1219 shape_function<dofs_per_cell; ++shape_function)
1221 const int snc = shape_function_data[shape_function].single_nonzero_component;
1227 const Number value = dof_values(shape_function);
1228 if (value == Number())
1233 const unsigned int comp =
1234 shape_function_data[shape_function].single_nonzero_component_index;
1236 const ::Tensor < 1, spacedim> *shape_gradient_ptr =
1237 &shape_gradients[snc][0];
1240 const unsigned int ii = indices[0];
1241 const unsigned int jj = indices[1];
1243 for (
unsigned int q_point = 0; q_point < n_quadrature_points;
1244 ++q_point, ++shape_gradient_ptr)
1246 divergences[q_point][jj] += value * (*shape_gradient_ptr)[ii];
1251 for (
unsigned int d = 0;
1253 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
1265 template <
int dim,
int spacedim>
1266 template <
class InputVector>
1274 Assert (fe_values->present_cell.get() != 0,
1275 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1277 fe_values->present_cell->n_dofs_for_dof_handler());
1281 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1282 internal::do_function_values<dim,spacedim>
1283 (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
1288 template <
int dim,
int spacedim>
1289 template <
class InputVector>
1297 Assert (fe_values->present_cell.get() != 0,
1298 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1300 fe_values->present_cell->n_dofs_for_dof_handler());
1304 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1305 internal::do_function_derivatives<1,dim,spacedim>
1306 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, gradients);
1311 template <
int dim,
int spacedim>
1312 template <
class InputVector>
1320 Assert (fe_values->present_cell.get() != 0,
1321 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1323 fe_values->present_cell->n_dofs_for_dof_handler());
1327 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1328 internal::do_function_derivatives<2,dim,spacedim>
1329 (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, hessians);
1334 template <
int dim,
int spacedim>
1335 template <
class InputVector>
1343 Assert (fe_values->present_cell.get() != 0,
1344 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1346 fe_values->present_cell->n_dofs_for_dof_handler());
1350 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1351 internal::do_function_laplacians<dim,spacedim>
1352 (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, laplacians);
1357 template <
int dim,
int spacedim>
1358 template <
class InputVector>
1366 Assert (fe_values->present_cell.get() != 0,
1367 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1369 fe_values->present_cell->n_dofs_for_dof_handler());
1373 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1374 internal::do_function_derivatives<3,dim,spacedim>
1375 (dof_values, fe_values->finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
1380 template <
int dim,
int spacedim>
1381 template <
class InputVector>
1389 Assert (fe_values->present_cell.get() != 0,
1390 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1392 fe_values->present_cell->n_dofs_for_dof_handler());
1396 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1397 internal::do_function_values<dim,spacedim>
1398 (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
1404 template <
int dim,
int spacedim>
1405 template <
class InputVector>
1413 Assert (fe_values->present_cell.get() != 0,
1414 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1416 fe_values->present_cell->n_dofs_for_dof_handler());
1420 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1421 internal::do_function_derivatives<1,dim,spacedim>
1422 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, gradients);
1427 template <
int dim,
int spacedim>
1428 template <
class InputVector>
1436 Assert (fe_values->present_cell.get() != 0,
1437 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1439 fe_values->present_cell->n_dofs_for_dof_handler());
1443 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1444 internal::do_function_symmetric_gradients<dim,spacedim>
1445 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data,
1446 symmetric_gradients);
1451 template <
int dim,
int spacedim>
1452 template <
class InputVector>
1460 Assert (fe_values->present_cell.get() != 0,
1461 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1463 fe_values->present_cell->n_dofs_for_dof_handler());
1468 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1469 internal::do_function_divergences<dim,spacedim>
1470 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
1473 template <
int dim,
int spacedim>
1474 template <
class InputVector>
1482 Assert (fe_values->present_cell.get () != 0,
1483 ExcMessage (
"FEValues object is not reinited to any cell"));
1485 fe_values->present_cell->n_dofs_for_dof_handler ());
1489 fe_values->present_cell->get_interpolated_dof_values (fe_function, dof_values);
1490 internal::do_function_curls<dim,spacedim>
1491 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, curls);
1495 template <
int dim,
int spacedim>
1496 template <
class InputVector>
1504 Assert (fe_values->present_cell.get() != 0,
1505 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1507 fe_values->present_cell->n_dofs_for_dof_handler());
1511 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1512 internal::do_function_derivatives<2,dim,spacedim>
1513 (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, hessians);
1518 template <
int dim,
int spacedim>
1519 template <
class InputVector>
1527 Assert (laplacians.size() == fe_values->n_quadrature_points,
1529 Assert (fe_values->present_cell.get() != 0,
1530 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1531 Assert (fe_function.size() == fe_values->present_cell->n_dofs_for_dof_handler(),
1533 fe_values->present_cell->n_dofs_for_dof_handler()));
1537 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1538 internal::do_function_laplacians<dim,spacedim>
1539 (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, laplacians);
1543 template <
int dim,
int spacedim>
1544 template <
class InputVector>
1552 Assert (fe_values->present_cell.get() != 0,
1553 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
1555 fe_values->present_cell->n_dofs_for_dof_handler());
1559 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1560 internal::do_function_derivatives<3,dim,spacedim>
1561 (dof_values, fe_values->finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
1566 template <
int dim,
int spacedim>
1567 template <
class InputVector>
1569 SymmetricTensor<2, dim, spacedim>::
1570 get_function_values(
const InputVector &fe_function,
1575 Assert(fe_values->present_cell.get() != 0,
1576 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
1578 fe_values->present_cell->n_dofs_for_dof_handler());
1582 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1583 internal::do_function_values<dim,spacedim>
1584 (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
1589 template <
int dim,
int spacedim>
1590 template <
class InputVector>
1592 SymmetricTensor<2, dim, spacedim>::
1593 get_function_divergences(
const InputVector &fe_function,
1598 Assert(fe_values->present_cell.get() != 0,
1599 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
1601 fe_values->present_cell->n_dofs_for_dof_handler());
1606 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1607 internal::do_function_divergences<dim,spacedim>
1608 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
1611 template <
int dim,
int spacedim>
1612 template <
class InputVector>
1614 Tensor<2, dim, spacedim>::
1615 get_function_values(
const InputVector &fe_function,
1620 Assert(fe_values->present_cell.get() != 0,
1621 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
1623 fe_values->present_cell->n_dofs_for_dof_handler());
1627 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1628 internal::do_function_values<dim,spacedim>
1629 (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
1634 template <
int dim,
int spacedim>
1635 template <
class InputVector>
1637 Tensor<2, dim, spacedim>::
1638 get_function_divergences(
const InputVector &fe_function,
1643 Assert(fe_values->present_cell.get() != 0,
1644 ExcMessage(
"FEValues object is not reinit'ed to any cell"));
1646 fe_values->present_cell->n_dofs_for_dof_handler());
1651 fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
1652 internal::do_function_divergences<dim,spacedim>
1653 (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
1662 template <
int dim,
int spacedim>
1671 scalars.resize (n_scalars);
1672 for (
unsigned int component=0; component<n_scalars; ++component)
1675 typedef ::FEValuesViews::Scalar<dim,spacedim> ScalarView;
1676 scalars[component].ScalarView::~ScalarView ();
1678 new (&scalars[component])
1690 const unsigned int n_vectors = (fe.
n_components() >= spacedim ?
1693 vectors.resize (n_vectors);
1694 for (
unsigned int component=0; component<n_vectors; ++component)
1697 typedef ::FEValuesViews::Vector<dim,spacedim>
VectorView;
1698 vectors[component].VectorView::~VectorView ();
1700 new (&vectors[component])
1707 const unsigned int n_symmetric_second_order_tensors
1711 symmetric_second_order_tensors.resize(n_symmetric_second_order_tensors);
1712 for (
unsigned int component = 0; component < n_symmetric_second_order_tensors; ++component)
1715 typedef ::FEValuesViews::SymmetricTensor<2, dim, spacedim> SymmetricTensorView;
1716 symmetric_second_order_tensors[component].SymmetricTensorView::~SymmetricTensorView();
1718 new (&symmetric_second_order_tensors[component])
1726 const unsigned int n_second_order_tensors
1730 second_order_tensors.resize(n_second_order_tensors);
1731 for (
unsigned int component = 0; component < n_second_order_tensors; ++component)
1734 typedef ::FEValuesViews::Tensor<2, dim, spacedim> TensorView;
1735 second_order_tensors[component].TensorView::~TensorView();
1737 new (&second_order_tensors[component])
1748 template <
int dim,
int spacedim>
1758 virtual ~CellIteratorBase ();
1785 n_dofs_for_dof_handler ()
const = 0;
1787 #include "fe_values.decl.1.inst" 1795 get_interpolated_dof_values (
const IndexSet &in,
1800 template <
int dim,
int spacedim>
1817 template <
int dim,
int spacedim>
1818 template <
typename CI>
1827 CellIterator (
const CI &cell);
1854 n_dofs_for_dof_handler ()
const;
1856 #include "fe_values.decl.2.inst" 1864 get_interpolated_dof_values (
const IndexSet &in,
1919 template <
int dim,
int spacedim>
1957 n_dofs_for_dof_handler ()
const;
1959 #include "fe_values.decl.2.inst" 1967 get_interpolated_dof_values (
const IndexSet &in,
1995 template <
int dim,
int spacedim>
1996 template <
typename CI>
2004 template <
int dim,
int spacedim>
2005 template <
typename CI>
2014 template <
int dim,
int spacedim>
2015 template <
typename CI>
2019 return cell->get_dof_handler().n_dofs();
2024 #include "fe_values.impl.1.inst" 2027 template <
int dim,
int spacedim>
2028 template <
typename CI>
2036 std::vector<types::global_dof_index> dof_indices (cell->get_fe().dofs_per_cell);
2037 cell->get_dof_indices (dof_indices);
2039 for (
unsigned int i=0; i<cell->get_fe().dofs_per_cell; ++i)
2040 out[i] = (in.
is_element (dof_indices[i]) ? 1 : 0);
2046 template <
int dim,
int spacedim>
2049 = (
"You have previously called the FEValues::reinit function with a\n" 2050 "cell iterator of type Triangulation<dim,spacedim>::cell_iterator. However,\n" 2051 "when you do this, you cannot call some functions in the FEValues\n" 2052 "class, such as the get_function_values/gradients/hessians/third_derivatives\n" 2053 "functions. If you need these functions, then you need to call\n" 2054 "FEValues::reinit with an iterator type that allows to extract\n" 2055 "degrees of freedom, such as DoFHandler<dim,spacedim>::cell_iterator.");
2058 template <
int dim,
int spacedim>
2067 template <
int dim,
int spacedim>
2076 template <
int dim,
int spacedim>
2085 #include "fe_values.impl.2.inst" 2088 template <
int dim,
int spacedim>
2103 template <
int dim,
int spacedim>
2114 numbers::signaling_nan<double>());
2158 template <
int dim,
int spacedim>
2179 template <
int dim,
int spacedim>
2188 this->shape_function_to_row_table
2189 = make_shape_function_to_row_table (
fe);
2193 unsigned int n_nonzero_shape_components = 0;
2194 for (
unsigned int i=0; i<
fe.dofs_per_cell; ++i)
2195 n_nonzero_shape_components +=
fe.n_nonzero_components (i);
2196 Assert (n_nonzero_shape_components >=
fe.dofs_per_cell,
2205 this->shape_values.reinit(n_nonzero_shape_components,
2207 this->shape_values.fill(numbers::signaling_nan<double>());
2212 this->shape_gradients.reinit(n_nonzero_shape_components,
2219 this->shape_hessians.reinit(n_nonzero_shape_components,
2226 this->shape_3rd_derivatives.reinit(n_nonzero_shape_components,
2235 template <
int dim,
int spacedim>
2253 template <
int dim,
int spacedim>
2263 fe(&
fe, typeid(*this).name()),
2267 ExcMessage (
"There is nothing useful you can do with an FEValues " 2268 "object when using a quadrature formula with zero " 2269 "quadrature points!"));
2275 template <
int dim,
int spacedim>
2278 tria_listener_refinement.disconnect ();
2279 tria_listener_mesh_transform.disconnect ();
2293 template <
typename Number,
typename Number2>
2295 do_function_values (
const Number2 *dof_values_ptr,
2296 const ::Table<2,double> &shape_values,
2297 std::vector<Number> &values)
2300 const unsigned int dofs_per_cell = shape_values.n_rows();
2301 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
2302 shape_values.n_cols() : values.size();
2306 std::fill_n (values.begin(), n_quadrature_points, Number());
2315 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2317 const Number2 value = dof_values_ptr[shape_func];
2318 if (value == Number2())
2321 const double *shape_value_ptr = &shape_values(shape_func, 0);
2322 for (
unsigned int point=0; point<n_quadrature_points; ++point)
2323 values[point] += value **shape_value_ptr++;
2327 template <
int dim,
int spacedim,
typename VectorType,
typename Number>
2329 do_function_values (
const Number *dof_values_ptr,
2330 const ::Table<2,double> &shape_values,
2332 const std::vector<unsigned int> &shape_function_to_row_table,
2334 const bool quadrature_points_fastest =
false,
2335 const unsigned int component_multiple = 1)
2338 for (
unsigned int i=0; i<values.size(); ++i)
2339 std::fill_n (values[i].begin(), values[i].size(),
2340 typename VectorType::value_type());
2345 if (dofs_per_cell == 0)
2348 const unsigned int n_quadrature_points = shape_values.n_cols();
2352 const unsigned result_components = n_components * component_multiple;
2353 (void)result_components;
2354 if (quadrature_points_fastest)
2357 for (
unsigned int i=0; i<values.size(); ++i)
2363 for (
unsigned int i=0; i<values.size(); ++i)
2370 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
2371 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2373 const Number value = dof_values_ptr[shape_func+mc*dofs_per_cell];
2374 if (value == Number())
2379 const unsigned int comp =
2381 + mc * n_components;
2383 row = shape_function_to_row_table[shape_func*n_components+comp];
2385 const double *shape_value_ptr = &shape_values(row, 0);
2387 if (quadrature_points_fastest)
2389 VectorType &values_comp = values[comp];
2390 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2391 values_comp[point] += value **shape_value_ptr++;
2394 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2395 values[point][comp] += value **shape_value_ptr++;
2398 for (
unsigned int c=0; c<n_components; ++c)
2404 row = shape_function_to_row_table[shape_func*n_components+c];
2406 const double *shape_value_ptr = &shape_values(row, 0);
2407 const unsigned int comp = c + mc * n_components;
2409 if (quadrature_points_fastest)
2411 VectorType &values_comp = values[comp];
2412 for (
unsigned int point=0;
point<n_quadrature_points;
2414 values_comp[point] += value **shape_value_ptr++;
2417 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2418 values[point][comp] += value **shape_value_ptr++;
2425 template <
int order,
int spacedim,
typename Number>
2427 do_function_derivatives (
const Number *dof_values_ptr,
2431 const unsigned int dofs_per_cell = shape_derivatives.size()[0];
2432 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
2433 shape_derivatives[0].size() : derivatives.size();
2446 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2448 const Number value = dof_values_ptr[shape_func];
2449 if (value == Number())
2453 = &shape_derivatives[shape_func][0];
2454 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2455 derivatives[point] += value *
2460 template <
int order,
int dim,
int spacedim,
typename Number>
2462 do_function_derivatives (
const Number *dof_values_ptr,
2465 const std::vector<unsigned int> &shape_function_to_row_table,
2467 const bool quadrature_points_fastest =
false,
2468 const unsigned int component_multiple = 1)
2471 for (
unsigned int i=0; i<derivatives.size(); ++i)
2472 std::fill_n (derivatives[i].begin(), derivatives[i].size(),
2478 if (dofs_per_cell == 0)
2482 const unsigned int n_quadrature_points = shape_derivatives[0].size();
2486 const unsigned result_components = n_components * component_multiple;
2487 (void)result_components;
2488 if (quadrature_points_fastest)
2491 for (
unsigned int i=0; i<derivatives.size(); ++i)
2497 for (
unsigned int i=0; i<derivatives.size(); ++i)
2504 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
2505 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2507 const Number value = dof_values_ptr[shape_func+mc*dofs_per_cell];
2508 if (value == Number())
2513 const unsigned int comp =
2515 + mc * n_components;
2517 row = shape_function_to_row_table[shape_func*n_components+comp];
2520 &shape_derivatives[row][0];
2522 if (quadrature_points_fastest)
2523 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2524 derivatives[comp][point] += value *
2527 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2528 derivatives[point][comp] += value *
2532 for (
unsigned int c=0; c<n_components; ++c)
2538 row = shape_function_to_row_table[shape_func*n_components+c];
2541 &shape_derivatives[row][0];
2542 const unsigned int comp = c + mc * n_components;
2544 if (quadrature_points_fastest)
2545 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2546 derivatives[comp][point] += value *
2549 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2550 derivatives[point][comp] += value *
2556 template <
int spacedim,
typename Number,
typename Number2>
2558 do_function_laplacians (
const Number2 *dof_values_ptr,
2560 std::vector<Number> &laplacians)
2562 const unsigned int dofs_per_cell = shape_hessians.size()[0];
2563 const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
2564 shape_hessians[0].size() : laplacians.size();
2568 std::fill_n (laplacians.begin(), n_quadrature_points, Number());
2573 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2575 const Number2 value = dof_values_ptr[shape_func];
2576 if (value == Number2())
2580 = &shape_hessians[shape_func][0];
2581 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2582 laplacians[point] += value * trace(*shape_hessian_ptr++);
2586 template <
int dim,
int spacedim,
typename VectorType,
typename Number>
2588 do_function_laplacians (
const Number *dof_values_ptr,
2591 const std::vector<unsigned int> &shape_function_to_row_table,
2592 std::vector<VectorType> &laplacians,
2593 const bool quadrature_points_fastest =
false,
2594 const unsigned int component_multiple = 1)
2597 for (
unsigned int i=0; i<laplacians.size(); ++i)
2598 std::fill_n (laplacians[i].begin(), laplacians[i].size(),
2599 typename VectorType::value_type());
2604 if (dofs_per_cell == 0)
2608 const unsigned int n_quadrature_points = shape_hessians[0].size();
2612 const unsigned result_components = n_components * component_multiple;
2613 (void)result_components;
2614 if (quadrature_points_fastest)
2617 for (
unsigned int i=0; i<laplacians.size(); ++i)
2623 for (
unsigned int i=0; i<laplacians.size(); ++i)
2630 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
2631 for (
unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
2633 const Number value = dof_values_ptr[shape_func+mc*dofs_per_cell];
2634 if (value == Number())
2639 const unsigned int comp =
2641 + mc * n_components;
2643 row = shape_function_to_row_table[shape_func*n_components+comp];
2646 &shape_hessians[row][0];
2647 if (quadrature_points_fastest)
2649 VectorType &laplacians_comp = laplacians[comp];
2650 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2651 laplacians_comp[point] += value * trace(*shape_hessian_ptr++);
2654 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2655 laplacians[point][comp] += value * trace(*shape_hessian_ptr++);
2658 for (
unsigned int c=0; c<n_components; ++c)
2664 row = shape_function_to_row_table[shape_func*n_components+c];
2667 &shape_hessians[row][0];
2668 const unsigned int comp = c + mc * n_components;
2670 if (quadrature_points_fastest)
2672 VectorType &laplacians_comp = laplacians[comp];
2673 for (
unsigned int point=0;
point<n_quadrature_points;
2675 laplacians_comp[point] += value * trace(*shape_hessian_ptr++);
2678 for (
unsigned int point=0;
point<n_quadrature_points; ++
point)
2679 laplacians[point][comp] += value * trace(*shape_hessian_ptr++);
2687 template <
int dim,
int spacedim>
2688 template <
class InputVector>
2690 const InputVector &fe_function,
2691 std::vector<typename InputVector::value_type> &values)
const 2693 typedef typename InputVector::value_type Number;
2697 Assert (present_cell.get() != 0,
2698 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
2700 present_cell->n_dofs_for_dof_handler());
2704 present_cell->get_interpolated_dof_values(fe_function, dof_values);
2705 internal::do_function_values (dof_values.begin(), this->finite_element_output.shape_values,
2711 template <
int dim,
int spacedim>
2712 template <
class InputVector>
2714 const InputVector &fe_function,
2715 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2716 std::vector<typename InputVector::value_type> &values)
const 2718 typedef typename InputVector::value_type Number;
2725 if (dofs_per_cell <= 100)
2727 Number dof_values[100];
2728 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2729 dof_values[i] = get_vector_element (fe_function, indices[i]);
2730 internal::do_function_values(&dof_values[0], this->finite_element_output.shape_values, values);
2735 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2736 dof_values[i] = get_vector_element (fe_function, indices[i]);
2737 internal::do_function_values(dof_values.
begin(), this->finite_element_output.shape_values,
2744 template <
int dim,
int spacedim>
2745 template <
class InputVector>
2747 const InputVector &fe_function,
2750 typedef typename InputVector::value_type Number;
2751 Assert (present_cell.get() != 0,
2752 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
2756 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
2760 present_cell->get_interpolated_dof_values(fe_function, dof_values);
2762 internal::do_function_values(dof_values.begin(), this->finite_element_output.shape_values, *fe,
2763 this->finite_element_output.shape_function_to_row_table, val);
2768 template <
int dim,
int spacedim>
2769 template <
class InputVector>
2771 const InputVector &fe_function,
2772 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2775 typedef typename InputVector::value_type Number;
2778 Assert (indices.size() % dofs_per_cell == 0,
2784 if (indices.size() <= 100)
2786 Number dof_values[100];
2787 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2788 dof_values[i] = get_vector_element (fe_function, indices[i]);
2789 internal::do_function_values(&dof_values[0], this->finite_element_output.shape_values, *fe,
2790 this->finite_element_output.shape_function_to_row_table, val,
2791 false, indices.size()/dofs_per_cell);
2796 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2797 dof_values[i] = get_vector_element (fe_function, indices[i]);
2798 internal::do_function_values(dof_values.
begin(), this->finite_element_output.shape_values, *fe,
2799 this->finite_element_output.shape_function_to_row_table, val,
2800 false, indices.size()/dofs_per_cell);
2806 template <
int dim,
int spacedim>
2807 template <
class InputVector>
2809 const InputVector &fe_function,
2810 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2811 VectorSlice<std::vector<std::vector<typename InputVector::value_type> > > values,
2812 bool quadrature_points_fastest)
const 2814 typedef typename InputVector::value_type Number;
2820 Assert (indices.size() % dofs_per_cell == 0,
2823 if (indices.size() <= 100)
2825 Number dof_values[100];
2826 for (
unsigned int i=0; i<indices.size(); ++i)
2827 dof_values[i] = get_vector_element (fe_function, indices[i]);
2828 internal::do_function_values(&dof_values[0], this->finite_element_output.shape_values, *fe,
2829 this->finite_element_output.shape_function_to_row_table, values,
2830 quadrature_points_fastest,
2831 indices.size()/dofs_per_cell);
2836 for (
unsigned int i=0; i<indices.size(); ++i)
2837 dof_values[i] = get_vector_element (fe_function, indices[i]);
2838 internal::do_function_values(dof_values.begin(), this->finite_element_output.shape_values, *fe,
2839 this->finite_element_output.shape_function_to_row_table, values,
2840 quadrature_points_fastest,
2841 indices.size()/dofs_per_cell);
2847 template <
int dim,
int spacedim>
2848 template <
class InputVector>
2851 const InputVector &fe_function,
2854 typedef typename InputVector::value_type Number;
2858 Assert (present_cell.get() != 0,
2859 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
2860 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
2864 present_cell->get_interpolated_dof_values(fe_function, dof_values);
2865 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_gradients,
2871 template <
int dim,
int spacedim>
2872 template <
class InputVector>
2874 const InputVector &fe_function,
2875 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2878 typedef typename InputVector::value_type Number;
2883 if (dofs_per_cell <= 100)
2885 Number dof_values[100];
2886 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2887 dof_values[i] = get_vector_element (fe_function, indices[i]);
2888 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_gradients,
2894 for (
unsigned int i=0; i<dofs_per_cell; ++i)
2895 dof_values[i] = get_vector_element (fe_function, indices[i]);
2896 internal::do_function_derivatives(dof_values.
begin(), this->finite_element_output.shape_gradients,
2904 template <
int dim,
int spacedim>
2905 template <
class InputVector>
2908 const InputVector &fe_function,
2911 typedef typename InputVector::value_type Number;
2914 Assert (present_cell.get() != 0,
2915 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
2916 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
2920 present_cell->get_interpolated_dof_values(fe_function, dof_values);
2922 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_gradients,
2923 *fe, this->finite_element_output.shape_function_to_row_table,
2929 template <
int dim,
int spacedim>
2930 template <
class InputVector>
2932 const InputVector &fe_function,
2933 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2935 bool quadrature_points_fastest)
const 2937 typedef typename InputVector::value_type Number;
2940 Assert (indices.size() % dofs_per_cell == 0,
2945 if (indices.size() <= 100)
2947 Number dof_values[100];
2948 for (
unsigned int i=0; i<indices.size(); ++i)
2949 dof_values[i] = get_vector_element (fe_function, indices[i]);
2950 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_gradients,
2951 *fe, this->finite_element_output.shape_function_to_row_table,
2952 gradients, quadrature_points_fastest,
2953 indices.size()/dofs_per_cell);
2958 for (
unsigned int i=0; i<indices.size(); ++i)
2959 dof_values[i] = get_vector_element (fe_function, indices[i]);
2960 internal::do_function_derivatives(dof_values.begin(),this->finite_element_output.shape_gradients,
2961 *fe, this->finite_element_output.shape_function_to_row_table,
2962 gradients, quadrature_points_fastest,
2963 indices.size()/dofs_per_cell);
2969 template <
int dim,
int spacedim>
2970 template <
class InputVector>
2976 typedef typename InputVector::value_type Number;
2980 Assert (present_cell.get() != 0,
2981 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
2982 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
2986 present_cell->get_interpolated_dof_values(fe_function, dof_values);
2987 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_hessians,
2993 template <
int dim,
int spacedim>
2994 template <
class InputVector>
2996 const InputVector &fe_function,
2997 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3000 typedef typename InputVector::value_type Number;
3003 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3005 if (dofs_per_cell <= 100)
3007 Number dof_values[100];
3008 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3009 dof_values[i] = get_vector_element (fe_function, indices[i]);
3010 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_hessians,
3016 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3017 dof_values[i] = get_vector_element (fe_function, indices[i]);
3018 internal::do_function_derivatives(dof_values.
begin(), this->finite_element_output.shape_hessians,
3026 template <
int dim,
int spacedim>
3027 template <
class InputVector>
3032 bool quadrature_points_fastest)
const 3034 typedef typename InputVector::value_type Number;
3037 Assert (present_cell.get() != 0,
3038 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
3039 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3043 present_cell->get_interpolated_dof_values(fe_function, dof_values);
3045 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_hessians,
3046 *fe, this->finite_element_output.shape_function_to_row_table,
3047 hes, quadrature_points_fastest);
3052 template <
int dim,
int spacedim>
3053 template <
class InputVector>
3055 const InputVector &fe_function,
3056 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3058 bool quadrature_points_fastest)
const 3060 typedef typename InputVector::value_type Number;
3063 Assert (indices.size() % dofs_per_cell == 0,
3065 if (indices.size() <= 100)
3067 Number dof_values[100];
3068 for (
unsigned int i=0; i<indices.size(); ++i)
3069 dof_values[i] = get_vector_element (fe_function, indices[i]);
3070 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_hessians,
3071 *fe, this->finite_element_output.shape_function_to_row_table,
3072 hessians, quadrature_points_fastest,
3073 indices.size()/dofs_per_cell);
3078 for (
unsigned int i=0; i<indices.size(); ++i)
3079 dof_values[i] = get_vector_element (fe_function, indices[i]);
3080 internal::do_function_derivatives(dof_values.begin(),this->finite_element_output.shape_hessians,
3081 *fe, this->finite_element_output.shape_function_to_row_table,
3082 hessians, quadrature_points_fastest,
3083 indices.size()/dofs_per_cell);
3089 template <
int dim,
int spacedim>
3090 template <
class InputVector>
3092 const InputVector &fe_function,
3093 std::vector<typename InputVector::value_type> &laplacians)
const 3095 typedef typename InputVector::value_type Number;
3099 Assert (present_cell.get() != 0,
3100 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
3101 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3105 present_cell->get_interpolated_dof_values(fe_function, dof_values);
3106 internal::do_function_laplacians(dof_values.begin(), this->finite_element_output.shape_hessians,
3112 template <
int dim,
int spacedim>
3113 template <
class InputVector>
3115 const InputVector &fe_function,
3116 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3117 std::vector<typename InputVector::value_type> &laplacians)
const 3119 typedef typename InputVector::value_type Number;
3124 if (dofs_per_cell <= 100)
3126 Number dof_values[100];
3127 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3128 dof_values[i] = get_vector_element (fe_function, indices[i]);
3129 internal::do_function_laplacians(&dof_values[0], this->finite_element_output.shape_hessians,
3135 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3136 dof_values[i] = get_vector_element (fe_function, indices[i]);
3137 internal::do_function_laplacians(dof_values.
begin(), this->finite_element_output.shape_hessians,
3144 template <
int dim,
int spacedim>
3145 template <
class InputVector>
3147 const InputVector &fe_function,
3150 typedef typename InputVector::value_type Number;
3151 Assert (present_cell.get() != 0,
3152 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
3155 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3159 present_cell->get_interpolated_dof_values(fe_function, dof_values);
3160 internal::do_function_laplacians(dof_values.begin(), this->finite_element_output.shape_hessians,
3161 *fe, this->finite_element_output.shape_function_to_row_table,
3167 template <
int dim,
int spacedim>
3168 template <
class InputVector>
3170 const InputVector &fe_function,
3171 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3174 typedef typename InputVector::value_type Number;
3177 Assert (indices.size() % dofs_per_cell == 0,
3181 if (indices.size() <= 100)
3183 Number dof_values[100];
3184 for (
unsigned int i=0; i<indices.size(); ++i)
3185 dof_values[i] = get_vector_element (fe_function, indices[i]);
3186 internal::do_function_laplacians(&dof_values[0], this->finite_element_output.shape_hessians,
3187 *fe, this->finite_element_output.shape_function_to_row_table,
3189 indices.size()/dofs_per_cell);
3194 for (
unsigned int i=0; i<indices.size(); ++i)
3195 dof_values[i] = get_vector_element (fe_function, indices[i]);
3196 internal::do_function_laplacians(dof_values.begin(),this->finite_element_output.shape_hessians,
3197 *fe, this->finite_element_output.shape_function_to_row_table,
3199 indices.size()/dofs_per_cell);
3205 template <
int dim,
int spacedim>
3206 template <
class InputVector>
3208 const InputVector &fe_function,
3209 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3210 std::vector<std::vector<typename InputVector::value_type> > &laplacians,
3211 bool quadrature_points_fastest)
const 3213 typedef typename InputVector::value_type Number;
3214 Assert (indices.size() % dofs_per_cell == 0,
3218 if (indices.size() <= 100)
3220 Number dof_values[100];
3221 for (
unsigned int i=0; i<indices.size(); ++i)
3222 dof_values[i] = get_vector_element (fe_function, indices[i]);
3223 internal::do_function_laplacians(&dof_values[0], this->finite_element_output.shape_hessians,
3224 *fe, this->finite_element_output.shape_function_to_row_table,
3225 laplacians, quadrature_points_fastest,
3226 indices.size()/dofs_per_cell);
3231 for (
unsigned int i=0; i<indices.size(); ++i)
3232 dof_values[i] = get_vector_element (fe_function, indices[i]);
3233 internal::do_function_laplacians(dof_values.begin(),this->finite_element_output.shape_hessians,
3234 *fe, this->finite_element_output.shape_function_to_row_table,
3235 laplacians, quadrature_points_fastest,
3236 indices.size()/dofs_per_cell);
3242 template <
int dim,
int spacedim>
3243 template <
class InputVector>
3249 typedef typename InputVector::value_type Number;
3253 Assert (present_cell.get() != 0,
3254 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
3255 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3259 present_cell->get_interpolated_dof_values(fe_function, dof_values);
3260 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_3rd_derivatives,
3266 template <
int dim,
int spacedim>
3267 template <
class InputVector>
3269 const InputVector &fe_function,
3270 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3273 typedef typename InputVector::value_type Number;
3276 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3278 if (dofs_per_cell <= 100)
3280 Number dof_values[100];
3281 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3282 dof_values[i] = get_vector_element (fe_function, indices[i]);
3283 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_3rd_derivatives,
3289 for (
unsigned int i=0; i<dofs_per_cell; ++i)
3290 dof_values[i] = get_vector_element (fe_function, indices[i]);
3291 internal::do_function_derivatives(dof_values.
begin(), this->finite_element_output.shape_3rd_derivatives,
3299 template <
int dim,
int spacedim>
3300 template <
class InputVector>
3305 bool quadrature_points_fastest)
const 3307 typedef typename InputVector::value_type Number;
3310 Assert (present_cell.get() != 0,
3311 ExcMessage (
"FEValues object is not reinit'ed to any cell"));
3312 AssertDimension (fe_function.size(), present_cell->n_dofs_for_dof_handler());
3316 present_cell->get_interpolated_dof_values(fe_function, dof_values);
3318 internal::do_function_derivatives(dof_values.begin(), this->finite_element_output.shape_3rd_derivatives,
3319 *fe, this->finite_element_output.shape_function_to_row_table,
3320 third, quadrature_points_fastest);
3325 template <
int dim,
int spacedim>
3326 template <
class InputVector>
3328 const InputVector &fe_function,
3329 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
3331 bool quadrature_points_fastest)
const 3333 typedef typename InputVector::value_type Number;
3336 Assert (indices.size() % dofs_per_cell == 0,
3338 if (indices.size() <= 100)
3340 Number dof_values[100];
3341 for (
unsigned int i=0; i<indices.size(); ++i)
3342 dof_values[i] = get_vector_element (fe_function, indices[i]);
3343 internal::do_function_derivatives(&dof_values[0], this->finite_element_output.shape_3rd_derivatives,
3344 *fe, this->finite_element_output.shape_function_to_row_table,
3345 third_derivatives, quadrature_points_fastest,
3346 indices.size()/dofs_per_cell);
3351 for (
unsigned int i=0; i<indices.size(); ++i)
3352 dof_values[i] = get_vector_element (fe_function, indices[i]);
3353 internal::do_function_derivatives(dof_values.begin(),this->finite_element_output.shape_3rd_derivatives,
3354 *fe, this->finite_element_output.shape_function_to_row_table,
3355 third_derivatives, quadrature_points_fastest,
3356 indices.size()/dofs_per_cell);
3362 template <
int dim,
int spacedim>
3366 return *present_cell;
3371 template <
int dim,
int spacedim>
3372 const std::vector<Tensor<1,spacedim> > &
3377 return this->mapping_output.normal_vectors;
3382 template <
int dim,
int spacedim>
3383 std::vector<Point<spacedim> >
3390 std::vector<Point<spacedim> > tmp (this->mapping_output.normal_vectors.size());
3391 for (
unsigned int q=0; q<this->mapping_output.normal_vectors.size(); ++q)
3399 template <
int dim,
int spacedim>
3413 template <
int dim,
int spacedim>
3417 return (
sizeof(this->update_flags) +
3419 sizeof (cell_similarity) +
3433 template <
int dim,
int spacedim>
3445 flags |= mapping->requires_update_flags (flags);
3451 template <
int dim,
int spacedim>
3462 tria_listener_refinement.disconnect ();
3463 tria_listener_mesh_transform.disconnect ();
3464 present_cell.reset ();
3468 template <
int dim,
int spacedim>
3473 if (present_cell.get() != 0)
3475 if (&cell->get_triangulation() !=
3477 ->get_triangulation())
3484 invalidate_present_cell();
3485 tria_listener_refinement =
3486 cell->get_triangulation().signals.any_change.connect
3489 tria_listener_mesh_transform =
3490 cell->get_triangulation().signals.mesh_movement.connect
3500 tria_listener_refinement =
3501 cell->get_triangulation().signals.post_refinement.connect
3504 tria_listener_mesh_transform =
3505 cell->get_triangulation().signals.mesh_movement.connect
3512 template <
int dim,
int spacedim>
3540 if (this->present_cell.get() == 0)
3548 cell_similarity = (cell->is_translation_of
3558 (*this->present_cell)->direction_flag()
3559 != cell->direction_flag() )
3568 template <
int dim,
int spacedim>
3572 return cell_similarity;
3576 template <
int dim,
int spacedim>
3580 template <
int dim,
int spacedim>
3585 template <
int dim,
int spacedim>
3591 template <
int dim,
int spacedim>
3609 template <
int dim,
int spacedim>
3626 template <
int dim,
int spacedim>
3633 if (dim != spacedim-1)
3637 const UpdateFlags flags = this->compute_update_flags (update_flags);
3641 this->mapping_output.initialize(this->n_quadrature_points, flags);
3642 this->finite_element_output.initialize(this->n_quadrature_points, *this->fe, flags);
3652 this->finite_element_output);
3661 this->update_flags = flags;
3666 this->mapping_data.reset (mapping_get_data.
return_value());
3678 template <
typename Type,
typename Po
inter,
typename Iterator>
3680 reset_pointer_in_place_if_possible
3681 (std_cxx11::unique_ptr<Pointer> &present_cell,
3682 const Iterator &new_cell)
3687 if (present_cell.get()
3689 (
typeid(*present_cell.get()) ==
typeid(Type)))
3692 static_cast<const Type *
>(present_cell.get())->~Type();
3695 new(
const_cast<void *
>(
static_cast<const void *
>(present_cell.get()))) Type(new_cell);
3699 present_cell.reset (
new Type(new_cell));
3704 template <
int dim,
int spacedim>
3709 this->maybe_invalidate_previous_present_cell (cell);
3710 this->check_cell_similarity(cell);
3712 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::TriaCellIterator>
3713 (this->present_cell, cell);
3725 template <
int dim,
int spacedim>
3726 template <
template <
int,
int>
class DoFHandlerType,
bool lda>
3739 this->maybe_invalidate_previous_present_cell (cell);
3740 this->check_cell_similarity(cell);
3742 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template
3745 (this->present_cell, cell);
3757 template <
int dim,
int spacedim>
3766 this->cell_similarity
3767 = this->get_mapping().fill_fe_values(*this->present_cell,
3768 this->cell_similarity,
3770 *this->mapping_data,
3771 this->mapping_output);
3778 this->get_fe().fill_fe_values(*this->present_cell,
3779 this->cell_similarity,
3781 this->get_mapping(),
3782 *this->mapping_data,
3783 this->mapping_output,
3785 this->finite_element_output);
3790 template <
int dim,
int spacedim>
3802 template <
int dim,
int spacedim>
3804 const unsigned int dofs_per_cell,
3815 present_face_index (
numbers::invalid_unsigned_int),
3816 quadrature(quadrature)
3821 template <
int dim,
int spacedim>
3822 const std::vector<Tensor<1,spacedim> > &
3827 return this->mapping_output.boundary_forms;
3832 template <
int dim,
int spacedim>
3843 template <
int dim,
int spacedim>
3846 template <
int dim,
int spacedim>
3850 template <
int dim,
int spacedim>
3867 template <
int dim,
int spacedim>
3883 template <
int dim,
int spacedim>
3887 const UpdateFlags flags = this->compute_update_flags (update_flags);
3891 this->mapping_output.initialize(this->n_quadrature_points, flags);
3892 this->finite_element_output.initialize(this->n_quadrature_points, *this->fe, flags);
3902 this->finite_element_output);
3911 this->update_flags = flags;
3916 this->mapping_data.reset (mapping_get_data.
return_value());
3923 template <
int dim,
int spacedim>
3924 template <
template <
int,
int>
class DoFHandlerType,
bool lda>
3928 const unsigned int face_no)
3936 cell->get_dof_handler().get_fe()[cell->active_fe_index ()]),
3942 this->maybe_invalidate_previous_present_cell (cell);
3943 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template
3946 (this->present_cell, cell);
3953 do_reinit (face_no);
3958 template <
int dim,
int spacedim>
3960 const unsigned int face_no)
3965 this->maybe_invalidate_previous_present_cell (cell);
3966 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::TriaCellIterator>
3967 (this->present_cell, cell);
3974 do_reinit (face_no);
3979 template <
int dim,
int spacedim>
3984 this->present_face_index=cell->face_index(face_no);
3988 this->get_mapping().fill_fe_face_values(*this->present_cell,
3991 *this->mapping_data,
3992 this->mapping_output);
3995 this->get_fe().fill_fe_face_values(*this->present_cell,
3998 this->get_mapping(),
3999 *this->mapping_data,
4000 this->mapping_output,
4002 this->finite_element_output);
4009 template <
int dim,
int spacedim>
4012 template <
int dim,
int spacedim>
4017 template <
int dim,
int spacedim>
4034 template <
int dim,
int spacedim>
4050 template <
int dim,
int spacedim>
4054 const UpdateFlags flags = this->compute_update_flags (update_flags);
4058 this->mapping_output.initialize(this->n_quadrature_points, flags);
4059 this->finite_element_output.initialize(this->n_quadrature_points, *this->fe, flags);
4070 this->finite_element_output);
4079 this->update_flags = flags;
4084 this->mapping_data.reset (mapping_get_data.
return_value());
4090 template <
int dim,
int spacedim>
4091 template <
template <
int,
int>
class DoFHandlerType,
bool lda>
4094 const unsigned int face_no,
4095 const unsigned int subface_no)
4103 cell->get_dof_handler().get_fe()[cell->active_fe_index ()]),
4115 Assert (cell->face(face_no)->has_children() ||
4116 subface_no < GeometryInfo<dim>::max_children_per_face,
4118 Assert (!cell->face(face_no)->has_children() ||
4119 subface_no < cell->face(face_no)->number_of_children(),
4120 ExcIndexRange (subface_no, 0, cell->face(face_no)->number_of_children()));
4121 Assert (cell->has_children() ==
false,
4122 ExcMessage (
"You can't use subface data for cells that are " 4123 "already refined. Iterate over their children " 4124 "instead in these cases."));
4126 this->maybe_invalidate_previous_present_cell (cell);
4127 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template
4130 (this->present_cell, cell);
4137 do_reinit (face_no, subface_no);
4141 template <
int dim,
int spacedim>
4143 const unsigned int face_no,
4144 const unsigned int subface_no)
4148 Assert (subface_no < cell->face(face_no)->n_children(),
4149 ExcIndexRange (subface_no, 0, cell->face(face_no)->n_children()));
4151 this->maybe_invalidate_previous_present_cell (cell);
4152 reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::TriaCellIterator>
4153 (this->present_cell, cell);
4160 do_reinit (face_no, subface_no);
4165 template <
int dim,
int spacedim>
4167 const unsigned int subface_no)
4173 if (!cell->face(face_no)->has_children())
4177 this->present_face_index=cell->face_index(face_no);
4179 this->present_face_index=cell->face(face_no)->child_index(subface_no);
4187 switch (cell->subface_case(face_no))
4192 subface_index=cell->face(face_no)->child_index(subface_no);
4196 subface_index=cell->face(face_no)->child(subface_no/2)->child_index(subface_no%2);
4204 subface_index=cell->face(face_no)->child(0)->child_index(subface_no);
4207 subface_index=cell->face(face_no)->child_index(1);
4218 subface_index=cell->face(face_no)->child_index(0);
4222 subface_index=cell->face(face_no)->child(1)->child_index(subface_no-1);
4234 this->present_face_index=subface_index;
4240 this->get_mapping().fill_fe_subface_values(*this->present_cell,
4244 *this->mapping_data,
4245 this->mapping_output);
4248 this->get_fe().fill_fe_subface_values(*this->present_cell,
4252 this->get_mapping(),
4253 *this->mapping_data,
4254 this->mapping_output,
4256 this->finite_element_output);
4261 #define SPLIT_INSTANTIATIONS_COUNT 2 4262 #ifndef SPLIT_INSTANTIATIONS_INDEX 4263 #define SPLIT_INSTANTIATIONS_INDEX 0 4265 #include "fe_values.inst" 4267 DEAL_II_NAMESPACE_CLOSE
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
Transformed quadrature weights.
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
static const unsigned int invalid_unsigned_int
Cache(const FEValuesBase< dim, spacedim > &fe_values)
#define AssertDimension(dim1, dim2)
unsigned int n_nonzero_components(const unsigned int i) const
static ::ExceptionBase & ExcAccessToUninitializedField()
static const unsigned int n_independent_components
const unsigned int dofs_per_cell
const unsigned int component
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
FEValuesBase(const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe)
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
static TableIndices< rank > unrolled_to_component_indices(const unsigned int i)
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
Outer normal vector, not normalized.
const FiniteElement< dim, spacedim > & get_fe() const
Scalar & operator=(const Scalar< dim, spacedim > &)
void get_function_divergences(const InputVector &fe_function, std::vector< typename ProductType< divergence_type, typename InputVector::value_type >::type > &divergences) const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
TriaCellIterator(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
Transformed quadrature points.
void do_reinit(const unsigned int face_no)
virtual void get_interpolated_dof_values(const IndexSet &in, Vector< IndexSet::value_type > &out) const
std::size_t memory_consumption() const
void transform(std::vector< Tensor< 1, spacedim > > &transformed, const std::vector< Tensor< 1, dim > > &original, MappingType mapping) const 1
bool is_primitive() const
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Tensor & operator=(const Tensor< rank_, dim, OtherNumber > &rhs)
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
static unsigned int component_to_unrolled_index(const TableIndices< rank > &indices)
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
static ::ExceptionBase & ExcMessage(std::string arg1)
CellSimilarity::Similarity get_cell_similarity() const
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
unsigned int global_dof_index
Third derivatives of shape functions.
void get_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &third_derivatives) const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell, const unsigned int face_no)
Abstract base class for mapping classes.
const ComponentMask & get_nonzero_components(const unsigned int i) const
std::vector< ShapeFunctionData > shape_function_data
const unsigned int first_vector_component
const Triangulation< dim, spacedim >::cell_iterator cell
void invalidate_present_cell()
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Vector & operator=(const Vector< dim, spacedim > &)
static const char *const message_string
void get_function_symmetric_gradients(const InputVector &fe_function, std::vector< typename ProductType< symmetric_gradient_type, typename InputVector::value_type >::type > &symmetric_gradients) const
Second derivatives of shape functions.
Gradient of volume element.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
std::size_t memory_consumption() const
const unsigned int dofs_per_cell
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
const unsigned int n_quadrature_points
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void get_function_curls(const InputVector &fe_function, std::vector< typename ProductType< curl_type, typename InputVector::value_type >::type > &curls) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
ArrayView< ElementType > make_array_view(std::vector< ElementType > &vector)
unsigned int n_components() const
FEFaceValuesBase(const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim-1 > &quadrature)
void initialize(const UpdateFlags update_flags)
Shape function gradients.
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
void initialize(const UpdateFlags update_flags)
static ::ExceptionBase & ExcNotMultiple(int arg1, int arg2)
static ::ExceptionBase & ExcNotImplemented()
static unsigned int n_threads()
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
unsigned int size(const unsigned int i) const
bool is_element(const size_type index) const
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
void initialize(const UpdateFlags update_flags)
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
std::vector< Point< spacedim > > get_normal_vectors() const 1
Point< 3 > point(const gp_Pnt &p)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell)
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim-1 > &face_quadrature, const UpdateFlags update_flags)
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim-1 > &quadrature, const UpdateFlags update_flags)
SymmetricTensor & operator=(const Number d)
std::size_t memory_consumption() const
virtual types::global_dof_index n_dofs_for_dof_handler() const
std::vector< ShapeFunctionData > shape_function_data
Task< RT > new_task(const std_cxx11::function< RT()> &function)
static ::ExceptionBase & ExcInternalError()