Reference documentation for deal.II version 8.5.1
fe_evaluation.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2017 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii__matrix_free_fe_evaluation_h
18 #define dealii__matrix_free_fe_evaluation_h
19 
20 
21 #include <deal.II/base/array_view.h>
22 #include <deal.II/base/config.h>
23 #include <deal.II/base/exceptions.h>
24 #include <deal.II/base/smartpointer.h>
25 #include <deal.II/base/symmetric_tensor.h>
26 #include <deal.II/base/template_constraints.h>
27 #include <deal.II/base/vectorization.h>
28 #include <deal.II/matrix_free/mapping_data_on_the_fly.h>
29 #include <deal.II/matrix_free/matrix_free.h>
30 #include <deal.II/matrix_free/shape_info.h>
31 #include <deal.II/matrix_free/evaluation_kernels.h>
32 #include <deal.II/matrix_free/tensor_product_kernels.h>
33 
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 
38 
39 // forward declarations
40 namespace LinearAlgebra
41 {
42  namespace distributed
43  {
44  template <typename> class Vector;
45  }
46 }
47 namespace internal
48 {
50 }
51 
52 template <int dim, int fe_degree, int n_q_points_1d = fe_degree+1,
53  int n_components_ = 1, typename Number = double > class FEEvaluation;
54 
55 
79 template <int dim, int n_components_, typename Number>
81 {
82 public:
83  typedef Number number_type;
86  static const unsigned int dimension = dim;
87  static const unsigned int n_components = n_components_;
88 
97 
106  void reinit (const unsigned int cell);
107 
120  template <typename DoFHandlerType, bool level_dof_access>
122 
133  void reinit (const typename Triangulation<dim>::cell_iterator &cell);
134 
143  unsigned int get_cell_data_number() const;
144 
151  internal::MatrixFreeFunctions::CellType get_cell_type() const;
152 
157  get_shape_info() const;
158 
162  void
164 
166 
203  template <typename VectorType>
204  void read_dof_values (const VectorType &src,
205  const unsigned int first_index = 0);
206 
235  template <typename VectorType>
236  void read_dof_values_plain (const VectorType &src,
237  const unsigned int first_index = 0);
238 
262  template<typename VectorType>
263  void distribute_local_to_global (VectorType &dst,
264  const unsigned int first_index = 0) const;
265 
292  template<typename VectorType>
293  void set_dof_values (VectorType &dst,
294  const unsigned int first_index = 0) const;
295 
297 
315  value_type get_dof_value (const unsigned int dof) const;
316 
327  void submit_dof_value (const value_type val_in,
328  const unsigned int dof);
329 
341  value_type get_value (const unsigned int q_point) const;
342 
354  void submit_value (const value_type val_in,
355  const unsigned int q_point);
356 
366  gradient_type get_gradient (const unsigned int q_point) const;
367 
380  void submit_gradient(const gradient_type grad_in,
381  const unsigned int q_point);
382 
394  get_hessian (const unsigned int q_point) const;
395 
404  gradient_type get_hessian_diagonal (const unsigned int q_point) const;
405 
416  value_type get_laplacian (const unsigned int q_point) const;
417 
418 #ifdef DOXYGEN
419  // doxygen does not anyhow mention functions coming from partial template
420  // specialization of the base class, in this case FEEvaluationAccess<dim,dim>.
421  // For now, hack-in those functions manually only to fix documentation:
422 
426  VectorizedArray<Number> get_divergence (const unsigned int q_point) const;
427 
431  SymmetricTensor<2, dim, VectorizedArray<Number> > get_symmetric_gradient (const unsigned int q_point) const;
432 
436  Tensor<1,(dim==2?1:dim), VectorizedArray<Number> > get_curl (const unsigned int q_point) const;
437 
441  void submit_divergence (const VectorizedArray<Number> div_in, const unsigned int q_point);
442 
446  void submit_symmetric_gradient (const SymmetricTensor<2, dim, VectorizedArray<Number> > grad_in, const unsigned int q_point);
447 
451  void submit_curl (const Tensor<1, dim==2?1:dim, VectorizedArray<Number> > curl_in, const unsigned int q_point);
452 
453 #endif
454 
463  value_type integrate_value () const;
464 
469  VectorizedArray<Number> JxW(const unsigned int q_point) const;
470 
472 
486 
496 
507  const VectorizedArray<Number> *begin_values () const;
508 
520 
533 
546 
559  const VectorizedArray<Number> *begin_hessians () const;
560 
574 
580  const std::vector<unsigned int> &
582 
590  get_scratch_data() const;
591 
593 
594 protected:
595 
603  FEEvaluationBase (const MatrixFree<dim,Number> &matrix_free,
604  const unsigned int fe_no,
605  const unsigned int quad_no,
606  const unsigned int fe_degree,
607  const unsigned int n_q_points);
608 
643  template <int n_components_other>
644  FEEvaluationBase (const Mapping<dim> &mapping,
645  const FiniteElement<dim> &fe,
646  const Quadrature<1> &quadrature,
647  const UpdateFlags update_flags,
648  const unsigned int first_selected_component,
650 
657  FEEvaluationBase (const FEEvaluationBase &other);
658 
666 
673  template<typename VectorType, typename VectorOperation>
674  void read_write_operation (const VectorOperation &operation,
675  VectorType *vectors[]) const;
676 
684  template<typename VectorType>
685  void read_dof_values_plain (const VectorType *src_data[]);
686 
691 
698 
712 
725 
740 
752  VectorizedArray<Number> *hessians_quad[n_components][(dim*(dim+1))/2];
753 
757  const unsigned int quad_no;
758 
763  const unsigned int n_fe_components;
764 
769  const unsigned int active_fe_index;
770 
775  const unsigned int active_quad_index;
776 
781 
787  const internal::MatrixFreeFunctions::DoFInfo *dof_info;
788 
796 
804 
810 
816 
824 
829 
834 
841 
848 
853  unsigned int cell;
854 
861  internal::MatrixFreeFunctions::CellType cell_type;
862 
866  unsigned int cell_data_number;
867 
874 
881 
888 
895 
902 
909 
914  std_cxx1x::shared_ptr<internal::MatrixFreeFunctions::MappingDataOnTheFly<dim,Number> > mapped_geometry;
915 
920  std::vector<types::global_dof_index> old_style_dof_indices;
921 
926  const unsigned int first_selected_component;
927 
932  mutable std::vector<types::global_dof_index> local_dof_indices;
933 
934 private:
939  void set_data_pointers();
940 
944  template <int, int, typename> friend class FEEvaluationBase;
945  template <int, int, int, int, typename> friend class FEEvaluation;
946 };
947 
948 
949 
957 template <int dim, int n_components_, typename Number>
958 class FEEvaluationAccess : public FEEvaluationBase<dim,n_components_,Number>
959 {
960 public:
961  typedef Number number_type;
964  static const unsigned int dimension = dim;
965  static const unsigned int n_components = n_components_;
967 
968 protected:
976  FEEvaluationAccess (const MatrixFree<dim,Number> &matrix_free,
977  const unsigned int fe_no,
978  const unsigned int quad_no,
979  const unsigned int fe_degree,
980  const unsigned int n_q_points);
981 
986  template <int n_components_other>
987  FEEvaluationAccess (const Mapping<dim> &mapping,
988  const FiniteElement<dim> &fe,
989  const Quadrature<1> &quadrature,
990  const UpdateFlags update_flags,
991  const unsigned int first_selected_component,
993 
997  FEEvaluationAccess (const FEEvaluationAccess &other);
998 
1003 };
1004 
1005 
1006 
1007 
1016 template <int dim, typename Number>
1017 class FEEvaluationAccess<dim,1,Number> : public FEEvaluationBase<dim,1,Number>
1018 {
1019 public:
1020  typedef Number number_type;
1023  static const unsigned int dimension = dim;
1025 
1035  value_type get_dof_value (const unsigned int dof) const;
1036 
1041  void submit_dof_value (const value_type val_in,
1042  const unsigned int dof);
1043 
1051  value_type get_value (const unsigned int q_point) const;
1052 
1060  void submit_value (const value_type val_in,
1061  const unsigned int q_point);
1062 
1068  gradient_type get_gradient (const unsigned int q_point) const;
1069 
1078  void submit_gradient(const gradient_type grad_in,
1079  const unsigned int q_point);
1080 
1088  get_hessian (unsigned int q_point) const;
1089 
1094  gradient_type get_hessian_diagonal (const unsigned int q_point) const;
1095 
1100  value_type get_laplacian (const unsigned int q_point) const;
1101 
1110  value_type integrate_value () const;
1111 
1112 protected:
1120  FEEvaluationAccess (const MatrixFree<dim,Number> &matrix_free,
1121  const unsigned int fe_no,
1122  const unsigned int quad_no,
1123  const unsigned int fe_degree,
1124  const unsigned int n_q_points);
1125 
1130  template <int n_components_other>
1131  FEEvaluationAccess (const Mapping<dim> &mapping,
1132  const FiniteElement<dim> &fe,
1133  const Quadrature<1> &quadrature,
1134  const UpdateFlags update_flags,
1135  const unsigned int first_selected_component,
1137 
1141  FEEvaluationAccess (const FEEvaluationAccess &other);
1142 
1147 };
1148 
1149 
1150 
1160 template <int dim, typename Number>
1161 class FEEvaluationAccess<dim,dim,Number> : public FEEvaluationBase<dim,dim,Number>
1162 {
1163 public:
1164  typedef Number number_type;
1167  static const unsigned int dimension = dim;
1168  static const unsigned int n_components = dim;
1170 
1175  gradient_type get_gradient (const unsigned int q_point) const;
1176 
1181  VectorizedArray<Number> get_divergence (const unsigned int q_point) const;
1182 
1190  get_symmetric_gradient (const unsigned int q_point) const;
1191 
1197  get_curl (const unsigned int q_point) const;
1198 
1206  get_hessian (const unsigned int q_point) const;
1207 
1212  gradient_type get_hessian_diagonal (const unsigned int q_point) const;
1213 
1222  void submit_gradient(const gradient_type grad_in,
1223  const unsigned int q_point);
1224 
1233  void submit_gradient(const Tensor<1,dim,Tensor<1,dim,VectorizedArray<Number> > > grad_in,
1234  const unsigned int q_point);
1235 
1244  void submit_divergence (const VectorizedArray<Number> div_in,
1245  const unsigned int q_point);
1246 
1256  const unsigned int q_point);
1257 
1262  void submit_curl (const Tensor<1,dim==2?1:dim,VectorizedArray<Number> > curl_in,
1263  const unsigned int q_point);
1264 
1265 protected:
1273  FEEvaluationAccess (const MatrixFree<dim,Number> &matrix_free,
1274  const unsigned int fe_no,
1275  const unsigned int quad_no,
1276  const unsigned int dofs_per_cell,
1277  const unsigned int n_q_points);
1278 
1283  template <int n_components_other>
1284  FEEvaluationAccess (const Mapping<dim> &mapping,
1285  const FiniteElement<dim> &fe,
1286  const Quadrature<1> &quadrature,
1287  const UpdateFlags update_flags,
1288  const unsigned int first_selected_component,
1290 
1294  FEEvaluationAccess (const FEEvaluationAccess &other);
1295 
1300 };
1301 
1302 
1312 template <typename Number>
1313 class FEEvaluationAccess<1,1,Number> : public FEEvaluationBase<1,1,Number>
1314 {
1315 public:
1316  typedef Number number_type;
1319  static const unsigned int dimension = 1;
1321 
1331  value_type get_dof_value (const unsigned int dof) const;
1332 
1337  void submit_dof_value (const value_type val_in,
1338  const unsigned int dof);
1339 
1347  value_type get_value (const unsigned int q_point) const;
1348 
1356  void submit_value (const value_type val_in,
1357  const unsigned int q_point);
1358 
1364  gradient_type get_gradient (const unsigned int q_point) const;
1365 
1374  void submit_gradient(const gradient_type grad_in,
1375  const unsigned int q_point);
1376 
1384  get_hessian (unsigned int q_point) const;
1385 
1390  gradient_type get_hessian_diagonal (const unsigned int q_point) const;
1391 
1396  value_type get_laplacian (const unsigned int q_point) const;
1397 
1406  value_type integrate_value () const;
1407 
1408 protected:
1416  FEEvaluationAccess (const MatrixFree<1,Number> &matrix_free,
1417  const unsigned int fe_no,
1418  const unsigned int quad_no,
1419  const unsigned int fe_degree,
1420  const unsigned int n_q_points);
1421 
1426  template <int n_components_other>
1427  FEEvaluationAccess (const Mapping<1> &mapping,
1428  const FiniteElement<1> &fe,
1429  const Quadrature<1> &quadrature,
1430  const UpdateFlags update_flags,
1431  const unsigned int first_selected_component,
1433 
1437  FEEvaluationAccess (const FEEvaluationAccess &other);
1438 
1443 };
1444 
1445 
1446 
1944 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
1945  typename Number >
1946 class FEEvaluation : public FEEvaluationAccess<dim,n_components_,Number>
1947 {
1948 public:
1950  typedef Number number_type;
1951  typedef typename BaseClass::value_type value_type;
1952  typedef typename BaseClass::gradient_type gradient_type;
1953  static const unsigned int dimension = dim;
1954  static const unsigned int n_components = n_components_;
1955  static const unsigned int static_n_q_points = Utilities::fixed_int_power<n_q_points_1d,dim>::value;
1956  static const unsigned int tensor_dofs_per_cell = Utilities::fixed_int_power<fe_degree+1,dim>::value;
1957 
1964  FEEvaluation (const MatrixFree<dim,Number> &matrix_free,
1965  const unsigned int fe_no = 0,
1966  const unsigned int quad_no = 0);
1967 
1994  FEEvaluation (const Mapping<dim> &mapping,
1995  const FiniteElement<dim> &fe,
1996  const Quadrature<1> &quadrature,
1997  const UpdateFlags update_flags,
1998  const unsigned int first_selected_component = 0);
1999 
2005  FEEvaluation (const FiniteElement<dim> &fe,
2006  const Quadrature<1> &quadrature,
2007  const UpdateFlags update_flags,
2008  const unsigned int first_selected_component = 0);
2009 
2020  template <int n_components_other>
2021  FEEvaluation (const FiniteElement<dim> &fe,
2023  const unsigned int first_selected_component = 0);
2024 
2031  FEEvaluation (const FEEvaluation &other);
2032 
2039  FEEvaluation &operator= (const FEEvaluation &other);
2040 
2049  void evaluate (const bool evaluate_values,
2050  const bool evaluate_gradients,
2051  const bool evaluate_hessians = false);
2052 
2060  void integrate (const bool integrate_values,
2061  const bool integrate_gradients);
2062 
2067  quadrature_point (const unsigned int q_point) const;
2068 
2074  const unsigned int dofs_per_cell;
2075 
2081  const unsigned int n_q_points;
2082 
2083 private:
2088  void check_template_arguments(const unsigned int fe_no,
2089  const unsigned int first_selected_component);
2090 
2095  VectorizedArray<Number> *values_dofs_actual[],
2098  VectorizedArray<Number> *hessians_quad[][(dim*(dim+1))/2],
2100  const bool evaluate_val,
2101  const bool evaluate_grad,
2102  const bool evaluate_lapl);
2103 
2108  VectorizedArray<Number> *values_dofs_actual[],
2112  const bool evaluate_val,
2113  const bool evaluate_grad);
2114 };
2115 
2116 
2117 
2118 namespace internal
2119 {
2120  namespace MatrixFreeFunctions
2121  {
2122  // a helper function to compute the number of DoFs of a DGP element at compile
2123  // time, depending on the degree
2124  template <int dim, int degree>
2125  struct DGP_dofs_per_cell
2126  {
2127  // this division is always without remainder
2128  static const unsigned int value =
2129  (DGP_dofs_per_cell<dim-1,degree>::value * (degree+dim)) / dim;
2130  };
2131 
2132  // base specialization: 1d elements have 'degree+1' degrees of freedom
2133  template <int degree>
2134  struct DGP_dofs_per_cell<1,degree>
2135  {
2136  static const unsigned int value = degree+1;
2137  };
2138  }
2139 }
2140 
2141 
2142 /*----------------------- Inline functions ----------------------------------*/
2143 
2144 #ifndef DOXYGEN
2145 
2146 
2147 
2148 /*----------------------- FEEvaluationBase ----------------------------------*/
2149 
2150 template <int dim, int n_components_, typename Number>
2151 inline
2154  const unsigned int fe_no_in,
2155  const unsigned int quad_no_in,
2156  const unsigned int fe_degree,
2157  const unsigned int n_q_points)
2158  :
2159  scratch_data_array (data_in.acquire_scratch_data()),
2160  quad_no (quad_no_in),
2161  n_fe_components (data_in.get_dof_info(fe_no_in).n_components),
2162  active_fe_index (fe_degree != numbers::invalid_unsigned_int ?
2163  data_in.get_dof_info(fe_no_in).fe_index_from_degree(fe_degree)
2164  :
2165  0),
2166  active_quad_index (fe_degree != numbers::invalid_unsigned_int ?
2167  data_in.get_mapping_info().
2168  mapping_data_gen[quad_no_in].
2169  quad_index_from_n_q_points(n_q_points)
2170  :
2171  0),
2172  matrix_info (&data_in),
2173  dof_info (&data_in.get_dof_info(fe_no_in)),
2174  mapping_info (&data_in.get_mapping_info()),
2175  data (&data_in.get_shape_info
2176  (fe_no_in, quad_no_in, active_fe_index,
2177  active_quad_index)),
2178  cartesian_data (0),
2179  jacobian (0),
2180  J_value (0),
2181  quadrature_weights (mapping_info->mapping_data_gen[quad_no].
2182  quadrature_weights[active_quad_index].begin()),
2183  quadrature_points (0),
2184  jacobian_grad (0),
2185  jacobian_grad_upper(0),
2186  cell (numbers::invalid_unsigned_int),
2187  cell_type (internal::MatrixFreeFunctions::undefined),
2188  cell_data_number (numbers::invalid_unsigned_int),
2189  dof_values_initialized (false),
2190  values_quad_initialized (false),
2191  gradients_quad_initialized(false),
2192  hessians_quad_initialized (false),
2193  values_quad_submitted (false),
2194  gradients_quad_submitted (false),
2195  first_selected_component (0)
2196 {
2197  set_data_pointers();
2198  Assert (matrix_info->mapping_initialized() == true,
2199  ExcNotInitialized());
2200  AssertDimension (matrix_info->get_size_info().vectorization_length,
2203  dof_info->dofs_per_cell[active_fe_index]/n_fe_components);
2204  AssertDimension (data->n_q_points,
2205  mapping_info->mapping_data_gen[quad_no].n_q_points[active_quad_index]);
2206  Assert (n_fe_components == 1 ||
2207  n_components == 1 ||
2208  n_components == n_fe_components,
2209  ExcMessage ("The underlying FE is vector-valued. In this case, the "
2210  "template argument n_components must be a the same "
2211  "as the number of underlying vector components."));
2212 
2213 
2214  // do not check for correct dimensions of data fields here, should be done
2215  // in derived classes
2216 }
2217 
2218 
2219 
2220 template <int dim, int n_components_, typename Number>
2221 template <int n_components_other>
2222 inline
2224 ::FEEvaluationBase (const Mapping<dim> &mapping,
2225  const FiniteElement<dim> &fe,
2226  const Quadrature<1> &quadrature,
2227  const UpdateFlags update_flags,
2228  const unsigned int first_selected_component,
2230  :
2231  scratch_data_array (new AlignedVector<VectorizedArray<Number> >()),
2232  quad_no (numbers::invalid_unsigned_int),
2233  n_fe_components (n_components_),
2234  active_fe_index (numbers::invalid_unsigned_int),
2235  active_quad_index (numbers::invalid_unsigned_int),
2236  matrix_info (0),
2237  dof_info (0),
2238  mapping_info (0),
2239  // select the correct base element from the given FE component
2240  data (new internal::MatrixFreeFunctions::ShapeInfo<Number>(quadrature, fe, fe.component_to_base_index(first_selected_component).first)),
2241  cartesian_data (0),
2242  jacobian (0),
2243  J_value (0),
2244  quadrature_weights (0),
2245  quadrature_points (0),
2246  jacobian_grad (0),
2247  jacobian_grad_upper(0),
2248  cell (0),
2249  cell_type (internal::MatrixFreeFunctions::general),
2250  cell_data_number (numbers::invalid_unsigned_int),
2251  dof_values_initialized (false),
2252  values_quad_initialized (false),
2253  gradients_quad_initialized(false),
2254  hessians_quad_initialized (false),
2255  values_quad_submitted (false),
2256  gradients_quad_submitted (false),
2257  // keep the number of the selected component within the current base element
2258  // for reading dof values
2259  first_selected_component (fe.component_to_base_index(first_selected_component).second)
2260 {
2261  const unsigned int base_element_number =
2262  fe.component_to_base_index(first_selected_component).first;
2263  set_data_pointers();
2264 
2265  Assert(other == 0 || other->mapped_geometry.get() != 0, ExcInternalError());
2266  if (other != 0 &&
2267  other->mapped_geometry->get_quadrature() == quadrature)
2268  mapped_geometry = other->mapped_geometry;
2269  else
2270  mapped_geometry.reset(new internal::MatrixFreeFunctions::
2271  MappingDataOnTheFly<dim,Number>(mapping, quadrature,
2272  update_flags));
2273  jacobian = mapped_geometry->get_inverse_jacobians().begin();
2274  J_value = mapped_geometry->get_JxW_values().begin();
2275  quadrature_points = mapped_geometry->get_quadrature_points().begin();
2276 
2277  Assert(fe.element_multiplicity(base_element_number) == 1 ||
2278  fe.element_multiplicity(base_element_number)-first_selected_component >= n_components_,
2279  ExcMessage("The underlying element must at least contain as many "
2280  "components as requested by this class"));
2281  (void) base_element_number;
2282 }
2283 
2284 
2285 
2286 template <int dim, int n_components_, typename Number>
2287 inline
2290  :
2291  scratch_data_array (other.matrix_info == 0 ?
2292  new AlignedVector<VectorizedArray<Number> >() :
2293  other.matrix_info->acquire_scratch_data()),
2294  quad_no (other.quad_no),
2295  n_fe_components (other.n_fe_components),
2296  active_fe_index (other.active_fe_index),
2297  active_quad_index (other.active_quad_index),
2298  matrix_info (other.matrix_info),
2299  dof_info (other.dof_info),
2300  mapping_info (other.mapping_info),
2301  data (other.matrix_info == 0 ?
2302  new internal::MatrixFreeFunctions::ShapeInfo<Number>(*other.data) :
2303  other.data),
2304  cartesian_data (0),
2305  jacobian (0),
2306  J_value (0),
2307  quadrature_weights (mapping_info != 0 ?
2308  mapping_info->mapping_data_gen[quad_no].
2309  quadrature_weights[active_quad_index].begin()
2310  :
2311  0),
2312  quadrature_points (0),
2313  jacobian_grad (0),
2314  jacobian_grad_upper(0),
2315  cell (numbers::invalid_unsigned_int),
2316  cell_type (internal::MatrixFreeFunctions::general),
2317  cell_data_number (numbers::invalid_unsigned_int),
2318  dof_values_initialized (false),
2319  values_quad_initialized (false),
2320  gradients_quad_initialized(false),
2321  hessians_quad_initialized (false),
2322  values_quad_submitted (false),
2323  gradients_quad_submitted (false),
2324  first_selected_component (other.first_selected_component)
2325 {
2326  set_data_pointers();
2327 
2328  // Create deep copy of mapped geometry for use in parallel...
2329  if (other.mapped_geometry.get() != 0)
2330  {
2331  mapped_geometry.reset
2333  MappingDataOnTheFly<dim,Number>(other.mapped_geometry->get_fe_values().get_mapping(),
2334  other.mapped_geometry->get_quadrature(),
2335  other.mapped_geometry->get_fe_values().get_update_flags()));
2336  jacobian = mapped_geometry->get_inverse_jacobians().begin();
2337  J_value = mapped_geometry->get_JxW_values().begin();
2338  quadrature_points = mapped_geometry->get_quadrature_points().begin();
2339  cell = 0;
2340  }
2341 }
2342 
2343 
2344 
2345 template <int dim, int n_components_, typename Number>
2346 inline
2350 {
2351  AssertDimension(quad_no, other.quad_no);
2352  AssertDimension(n_fe_components, other.n_fe_components);
2353  AssertDimension(active_fe_index, other.active_fe_index);
2354  AssertDimension(active_quad_index, other.active_quad_index);
2355  AssertDimension(first_selected_component, other.first_selected_component);
2356 
2357  // release old memory
2358  if (matrix_info == 0)
2359  {
2360  delete data;
2361  delete scratch_data_array;
2362  }
2363  else
2364  {
2365  matrix_info->release_scratch_data(scratch_data_array);
2366  }
2367 
2368  matrix_info = other.matrix_info;
2369  dof_info = other.dof_info;
2370  mapping_info = other.mapping_info;
2371  if (other.matrix_info == 0)
2372  {
2374  scratch_data_array = new AlignedVector<VectorizedArray<Number> >();
2375  }
2376  else
2377  {
2378  data = other.data;
2379  scratch_data_array = matrix_info->acquire_scratch_data();
2380  }
2381  set_data_pointers();
2382 
2383  cartesian_data = 0;
2384  jacobian = 0;
2385  J_value = 0;
2386  quadrature_weights = mapping_info != 0 ?
2387  mapping_info->mapping_data_gen[quad_no].
2388  quadrature_weights[active_quad_index].begin()
2389  :
2390  0;
2391  quadrature_points = 0;
2392  jacobian_grad = 0;
2393  jacobian_grad_upper = 0;
2395  cell_type = internal::MatrixFreeFunctions::general;
2396  cell_data_number = numbers::invalid_unsigned_int;
2397 
2398  // Create deep copy of mapped geometry for use in parallel...
2399  if (other.mapped_geometry.get() != 0)
2400  {
2401  mapped_geometry.reset
2403  MappingDataOnTheFly<dim,Number>(other.mapped_geometry->get_fe_values().get_mapping(),
2404  other.mapped_geometry->get_quadrature(),
2405  other.mapped_geometry->get_fe_values().get_update_flags()));
2406  jacobian = mapped_geometry->get_inverse_jacobians().begin();
2407  J_value = mapped_geometry->get_JxW_values().begin();
2408  quadrature_points = mapped_geometry->get_quadrature_points().begin();
2409  cell = 0;
2410  }
2411 
2412  return *this;
2413 }
2414 
2415 
2416 
2417 template <int dim, int n_components_, typename Number>
2418 inline
2420 {
2421  if (matrix_info != 0)
2422  {
2423  try
2424  {
2425  matrix_info->release_scratch_data(scratch_data_array);
2426  }
2427  catch (...)
2428  {}
2429  }
2430  else
2431  {
2432  delete scratch_data_array;
2433  delete data;
2434  }
2435 }
2436 
2437 
2438 
2439 template <int dim, int n_components_, typename Number>
2440 inline
2441 void
2444 {
2445  Assert(scratch_data_array != NULL, ExcInternalError());
2446 
2447  const unsigned int tensor_dofs_per_cell =
2448  Utilities::fixed_power<dim>(this->data->fe_degree+1);
2449  const unsigned int dofs_per_cell = this->data->dofs_per_cell;
2450  const unsigned int n_quadrature_points = this->data->n_q_points;
2451 
2452  const unsigned int shift = std::max(tensor_dofs_per_cell+1, dofs_per_cell)*
2453  (n_components_+2) + 2 * n_quadrature_points;
2454  const unsigned int allocated_size = shift + n_components_ * dofs_per_cell
2455  + (n_components_*(dim*dim+2*dim+1)*n_quadrature_points);
2456  scratch_data_array->resize_fast(allocated_size);
2457 
2458  // set the pointers to the correct position in the data array
2459  for (unsigned int c=0; c<n_components_; ++c)
2460  {
2461  this->values_dofs[c] = scratch_data_array->begin() + c*dofs_per_cell;
2462  this->values_quad[c] = scratch_data_array->begin() +
2463  n_components*dofs_per_cell+c*n_quadrature_points;
2464  for (unsigned int d=0; d<dim; ++d)
2465  this->gradients_quad[c][d] = scratch_data_array->begin() +
2466  n_components*(dofs_per_cell+n_quadrature_points) +
2467  (c*dim+d)*n_quadrature_points;
2468  for (unsigned int d=0; d<(dim*dim+dim)/2; ++d)
2469  this->hessians_quad[c][d] = scratch_data_array->begin() +
2470  n_components*((dim+1)*n_quadrature_points + dofs_per_cell) +
2471  (c*(dim*dim+dim)+d)*n_quadrature_points;
2472  }
2473  scratch_data = scratch_data_array->begin() + n_components_ * dofs_per_cell
2474  + (n_components_*(dim*dim+2*dim+1)*n_quadrature_points);
2475 }
2476 
2477 
2478 
2479 template <int dim, int n_components_, typename Number>
2480 inline
2481 void
2482 FEEvaluationBase<dim,n_components_,Number>::reinit (const unsigned int cell_in)
2483 {
2484  Assert (mapped_geometry == 0,
2485  ExcMessage("FEEvaluation was initialized without a matrix-free object."
2486  " Integer indexing is not possible"));
2487  if (mapped_geometry != 0)
2488  return;
2489  Assert (dof_info != 0, ExcNotInitialized());
2490  Assert (mapping_info != 0, ExcNotInitialized());
2491  AssertIndexRange (cell_in, dof_info->row_starts.size()-1);
2492  AssertDimension (((dof_info->cell_active_fe_index.size() > 0) ?
2493  dof_info->cell_active_fe_index[cell_in] : 0),
2494  active_fe_index);
2495  cell = cell_in;
2496  cell_type = mapping_info->get_cell_type(cell);
2497  cell_data_number = mapping_info->get_cell_data_index(cell);
2498 
2499  if (mapping_info->quadrature_points_initialized == true)
2500  {
2501  AssertIndexRange (cell_data_number, mapping_info->
2502  mapping_data_gen[quad_no].rowstart_q_points.size());
2503  const unsigned int index = mapping_info->mapping_data_gen[quad_no].
2504  rowstart_q_points[cell];
2505  AssertIndexRange (index, mapping_info->mapping_data_gen[quad_no].
2506  quadrature_points.size());
2507  quadrature_points =
2508  &mapping_info->mapping_data_gen[quad_no].quadrature_points[index];
2509  }
2510 
2511  if (cell_type == internal::MatrixFreeFunctions::cartesian)
2512  {
2513  cartesian_data = &mapping_info->cartesian_data[cell_data_number].first;
2514  J_value = &mapping_info->cartesian_data[cell_data_number].second;
2515  }
2516  else if (cell_type == internal::MatrixFreeFunctions::affine)
2517  {
2518  jacobian = &mapping_info->affine_data[cell_data_number].first;
2519  J_value = &mapping_info->affine_data[cell_data_number].second;
2520  }
2521  else
2522  {
2523  const unsigned int rowstart = mapping_info->
2524  mapping_data_gen[quad_no].rowstart_jacobians[cell_data_number];
2525  AssertIndexRange (rowstart, mapping_info->
2526  mapping_data_gen[quad_no].jacobians.size());
2527  jacobian =
2528  &mapping_info->mapping_data_gen[quad_no].jacobians[rowstart];
2529  if (mapping_info->JxW_values_initialized == true)
2530  {
2531  AssertIndexRange (rowstart, mapping_info->
2532  mapping_data_gen[quad_no].JxW_values.size());
2533  J_value = &(mapping_info->mapping_data_gen[quad_no].
2534  JxW_values[rowstart]);
2535  }
2536  if (mapping_info->second_derivatives_initialized == true)
2537  {
2538  AssertIndexRange(rowstart, mapping_info->
2539  mapping_data_gen[quad_no].jacobians_grad_diag.size());
2540  jacobian_grad = &mapping_info->mapping_data_gen[quad_no].
2541  jacobians_grad_diag[rowstart];
2542  AssertIndexRange(rowstart, mapping_info->
2543  mapping_data_gen[quad_no].jacobians_grad_upper.size());
2544  jacobian_grad_upper = &mapping_info->mapping_data_gen[quad_no].
2545  jacobians_grad_upper[rowstart];
2546  }
2547  }
2548 
2549 #ifdef DEBUG
2550  dof_values_initialized = false;
2551  values_quad_initialized = false;
2552  gradients_quad_initialized = false;
2553  hessians_quad_initialized = false;
2554 #endif
2555 }
2556 
2557 
2558 
2559 template <int dim, int n_components_, typename Number>
2560 template <typename DoFHandlerType, bool level_dof_access>
2561 inline
2562 void
2565 {
2566  Assert(matrix_info == 0,
2567  ExcMessage("Cannot use initialization from cell iterator if "
2568  "initialized from MatrixFree object. Use variant for "
2569  "on the fly computation with arguments as for FEValues "
2570  "instead"));
2571  Assert(mapped_geometry.get() != 0, ExcNotInitialized());
2572  mapped_geometry->reinit(static_cast<typename Triangulation<dim>::cell_iterator>(cell));
2573  local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2574  if (level_dof_access)
2575  cell->get_mg_dof_indices(local_dof_indices);
2576  else
2577  cell->get_dof_indices(local_dof_indices);
2578 }
2579 
2580 
2581 
2582 template <int dim, int n_components_, typename Number>
2583 inline
2584 void
2586 ::reinit (const typename Triangulation<dim>::cell_iterator &cell)
2587 {
2588  Assert(matrix_info == 0,
2589  ExcMessage("Cannot use initialization from cell iterator if "
2590  "initialized from MatrixFree object. Use variant for "
2591  "on the fly computation with arguments as for FEValues "
2592  "instead"));
2593  Assert(mapped_geometry.get() != 0, ExcNotInitialized());
2594  mapped_geometry->reinit(cell);
2595 }
2596 
2597 
2598 
2599 template <int dim, int n_components_, typename Number>
2600 inline
2601 unsigned int
2603 ::get_cell_data_number () const
2604 {
2606  return cell_data_number;
2607 }
2608 
2609 
2610 
2611 template <int dim, int n_components_, typename Number>
2612 inline
2613 internal::MatrixFreeFunctions::CellType
2615 {
2617  return cell_type;
2618 }
2619 
2620 
2621 
2622 template <int dim, int n_components_, typename Number>
2623 inline
2626 {
2627  Assert(data != 0, ExcInternalError());
2628  return *data;
2629 }
2630 
2631 
2632 
2633 template <int dim, int n_components_, typename Number>
2634 inline
2635 void
2638 {
2639  AssertDimension(JxW_values.size(), data->n_q_points);
2640  Assert (this->J_value != 0, ExcNotInitialized());
2641  if (this->cell_type == internal::MatrixFreeFunctions::cartesian ||
2642  this->cell_type == internal::MatrixFreeFunctions::affine)
2643  {
2644  Assert (this->mapping_info != 0, ExcNotImplemented());
2645  VectorizedArray<Number> J = this->J_value[0];
2646  for (unsigned int q=0; q<this->data->n_q_points; ++q)
2647  JxW_values[q] = J * this->quadrature_weights[q];
2648  }
2649  else
2650  for (unsigned int q=0; q<data->n_q_points; ++q)
2651  JxW_values[q] = this->J_value[q];
2652 }
2653 
2654 
2655 
2656 template <int dim, int n_components_, typename Number>
2657 inline
2659 FEEvaluationBase<dim,n_components_,Number>::JxW(const unsigned int q_point) const
2660 {
2661  Assert (this->J_value != 0, ExcNotInitialized());
2662  if (this->cell_type == internal::MatrixFreeFunctions::cartesian ||
2663  this->cell_type == internal::MatrixFreeFunctions::affine)
2664  {
2665  Assert (this->mapping_info != 0, ExcInternalError());
2666  return this->J_value[0] * this->quadrature_weights[q_point];
2667  }
2668  else
2669  return this->J_value[q_point];
2670 }
2671 
2672 
2673 
2674 namespace internal
2675 {
2676  // write access to generic vectors that have operator ().
2677  template <typename VectorType>
2678  inline
2679  typename VectorType::value_type &
2680  vector_access (VectorType &vec,
2681  const unsigned int entry)
2682  {
2683  return vec(entry);
2684  }
2685 
2686 
2687 
2688  // read access to generic vectors that have operator ().
2689  template <typename VectorType>
2690  inline
2691  typename VectorType::value_type
2692  vector_access (const VectorType &vec,
2693  const unsigned int entry)
2694  {
2695  return vec(entry);
2696  }
2697 
2698 
2699 
2700  // write access to distributed MPI vectors that have a local_element(uint)
2701  // method to access data in local index space, which is what we use in
2702  // DoFInfo and hence in read_dof_values etc.
2703  template <typename Number>
2704  inline
2705  Number &
2706  vector_access (LinearAlgebra::distributed::Vector<Number> &vec,
2707  const unsigned int entry)
2708  {
2709  return vec.local_element(entry);
2710  }
2711 
2712 
2713 
2714  // read access to distributed MPI vectors that have a local_element(uint)
2715  // method to access data in local index space, which is what we use in
2716  // DoFInfo and hence in read_dof_values etc.
2717  template <typename Number>
2718  inline
2719  Number
2720  vector_access (const LinearAlgebra::distributed::Vector<Number> &vec,
2721  const unsigned int entry)
2722  {
2723  return vec.local_element(entry);
2724  }
2725 
2726 
2727 
2728  // this is to make sure that the parallel partitioning in the
2729  // LinearAlgebra::distributed::Vector is really the same as stored in
2730  // MatrixFree
2731  template <typename VectorType>
2732  inline
2733  void check_vector_compatibility (const VectorType &vec,
2734  const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2735  {
2736  (void) vec;
2737  (void) dof_info;
2738 
2739  AssertDimension (vec.size(),
2740  dof_info.vector_partitioner->size());
2741  }
2742 
2743  template <typename Number>
2744  inline
2745  void check_vector_compatibility (const LinearAlgebra::distributed::Vector<Number> &vec,
2746  const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2747  {
2748  (void) vec;
2749  (void) dof_info;
2750  Assert (vec.partitioners_are_compatible(*dof_info.vector_partitioner),
2751  ExcMessage("The parallel layout of the given vector is not "
2752  "compatible with the parallel partitioning in MatrixFree. "
2753  "Use MatrixFree::initialize_dof_vector to get a "
2754  "compatible vector."));
2755  }
2756 
2757  // A class to use the same code to read from and write to vector
2758  template <typename Number>
2759  struct VectorReader
2760  {
2761  template <typename VectorType>
2762  void process_dof (const unsigned int index,
2763  VectorType &vec,
2764  Number &res) const
2765  {
2766  res = vector_access (const_cast<const VectorType &>(vec), index);
2767  }
2768 
2769  template <typename VectorType>
2770  void process_dof_global (const types::global_dof_index index,
2771  VectorType &vec,
2772  Number &res) const
2773  {
2774  res = const_cast<const VectorType &>(vec)(index);
2775  }
2776 
2777  void pre_constraints (const Number &,
2778  Number &res) const
2779  {
2780  res = Number();
2781  }
2782 
2783  template <typename VectorType>
2784  void process_constraint (const unsigned int index,
2785  const Number weight,
2786  VectorType &vec,
2787  Number &res) const
2788  {
2789  res += weight * vector_access (const_cast<const VectorType &>(vec), index);
2790  }
2791 
2792  void post_constraints (const Number &sum,
2793  Number &write_pos) const
2794  {
2795  write_pos = sum;
2796  }
2797 
2798  void process_empty (Number &res) const
2799  {
2800  res = Number();
2801  }
2802  };
2803 
2804  // A class to use the same code to read from and write to vector
2805  template <typename Number>
2806  struct VectorDistributorLocalToGlobal
2807  {
2808  template <typename VectorType>
2809  void process_dof (const unsigned int index,
2810  VectorType &vec,
2811  Number &res) const
2812  {
2813  vector_access (vec, index) += res;
2814  }
2815 
2816  template <typename VectorType>
2817  void process_dof_global (const types::global_dof_index index,
2818  VectorType &vec,
2819  Number &res) const
2820  {
2821  vec(index) += res;
2822  }
2823 
2824  void pre_constraints (const Number &input,
2825  Number &res) const
2826  {
2827  res = input;
2828  }
2829 
2830  template <typename VectorType>
2831  void process_constraint (const unsigned int index,
2832  const Number weight,
2833  VectorType &vec,
2834  Number &res) const
2835  {
2836  vector_access (vec, index) += weight * res;
2837  }
2838 
2839  void post_constraints (const Number &,
2840  Number &) const
2841  {
2842  }
2843 
2844  void process_empty (Number &) const
2845  {
2846  }
2847  };
2848 
2849 
2850  // A class to use the same code to read from and write to vector
2851  template <typename Number>
2852  struct VectorSetter
2853  {
2854  template <typename VectorType>
2855  void process_dof (const unsigned int index,
2856  VectorType &vec,
2857  Number &res) const
2858  {
2859  vector_access (vec, index) = res;
2860  }
2861 
2862  template <typename VectorType>
2863  void process_dof_global (const types::global_dof_index index,
2864  VectorType &vec,
2865  Number &res) const
2866  {
2867  vec(index) = res;
2868  }
2869 
2870  void pre_constraints (const Number &,
2871  Number &) const
2872  {
2873  }
2874 
2875  template <typename VectorType>
2876  void process_constraint (const unsigned int,
2877  const Number,
2878  VectorType &,
2879  Number &) const
2880  {
2881  }
2882 
2883  void post_constraints (const Number &,
2884  Number &) const
2885  {
2886  }
2887 
2888  void process_empty (Number &) const
2889  {
2890  }
2891  };
2892 
2893  // allows to select between block vectors and non-block vectors, which
2894  // allows to use a unified interface for extracting blocks on block vectors
2895  // and doing nothing on usual vectors
2896  template <typename VectorType, bool>
2897  struct BlockVectorSelector {};
2898 
2899  template <typename VectorType>
2900  struct BlockVectorSelector<VectorType,true>
2901  {
2902  typedef typename VectorType::BlockType BaseVectorType;
2903 
2904  static BaseVectorType *get_vector_component (VectorType &vec,
2905  const unsigned int component)
2906  {
2907  AssertIndexRange (component, vec.n_blocks());
2908  return &vec.block(component);
2909  }
2910  };
2911 
2912  template <typename VectorType>
2913  struct BlockVectorSelector<VectorType,false>
2914  {
2915  typedef VectorType BaseVectorType;
2916 
2917  static BaseVectorType *get_vector_component (VectorType &vec,
2918  const unsigned int)
2919  {
2920  return &vec;
2921  }
2922  };
2923 
2924  template <typename VectorType>
2925  struct BlockVectorSelector<std::vector<VectorType>,false>
2926  {
2927  typedef VectorType BaseVectorType;
2928 
2929  static BaseVectorType *get_vector_component (std::vector<VectorType> &vec,
2930  const unsigned int component)
2931  {
2932  AssertIndexRange (component, vec.size());
2933  return &vec[component];
2934  }
2935  };
2936 
2937  template <typename VectorType>
2938  struct BlockVectorSelector<std::vector<VectorType *>,false>
2939  {
2940  typedef VectorType BaseVectorType;
2941 
2942  static BaseVectorType *get_vector_component (std::vector<VectorType *> &vec,
2943  const unsigned int component)
2944  {
2945  AssertIndexRange (component, vec.size());
2946  return vec[component];
2947  }
2948  };
2949 }
2950 
2951 
2952 
2953 template <int dim, int n_components_, typename Number>
2954 template<typename VectorType, typename VectorOperation>
2955 inline
2956 void
2958 ::read_write_operation (const VectorOperation &operation,
2959  VectorType *src[]) const
2960 {
2961  // This functions processes all the functions read_dof_values,
2962  // distribute_local_to_global, and set_dof_values with the same code. The
2963  // distinction between these three cases is made by the input
2964  // VectorOperation that either reads values from a vector and puts the data
2965  // into the local data field or write local data into the vector. Certain
2966  // operations are no-ops for the given use case.
2967 
2968  // Case 1: No MatrixFree object given, simple case because we do not need to
2969  // process constraints and need not care about vectorization
2970  if (matrix_info == 0)
2971  {
2972  Assert (!local_dof_indices.empty(), ExcNotInitialized());
2973 
2974  unsigned int index = first_selected_component * this->data->dofs_per_cell;
2975  for (unsigned int comp = 0; comp<n_components; ++comp)
2976  {
2977  for (unsigned int i=0; i<this->data->dofs_per_cell; ++i, ++index)
2978  {
2979  operation.process_dof_global(local_dof_indices[this->data->lexicographic_numbering[index]],
2980  *src[0], values_dofs[comp][i][0]);
2981  for (unsigned int v=1; v<VectorizedArray<Number>::n_array_elements; ++v)
2982  operation.process_empty(values_dofs[comp][i][v]);
2983  }
2984  }
2985  return;
2986  }
2987 
2988  Assert (dof_info != 0, ExcNotInitialized());
2989  Assert (matrix_info->indices_initialized() == true,
2990  ExcNotInitialized());
2992 
2993  // loop over all local dofs. ind_local holds local number on cell, index
2994  // iterates over the elements of index_local_to_global and dof_indices
2995  // points to the global indices stored in index_local_to_global
2996  const unsigned int *dof_indices = dof_info->begin_indices(cell);
2997  const std::pair<unsigned short,unsigned short> *indicators =
2998  dof_info->begin_indicators(cell);
2999  const std::pair<unsigned short,unsigned short> *indicators_end =
3000  dof_info->end_indicators(cell);
3001  unsigned int ind_local = 0;
3002  const unsigned int dofs_per_cell = this->data->dofs_per_cell;
3003 
3004  const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
3005  const bool at_irregular_cell = n_irreg_components_filled > 0;
3006 
3007  // scalar case (or case when all components have the same degrees of freedom
3008  // and sit on a different vector each)
3009  if (n_fe_components == 1)
3010  {
3011  const unsigned int n_local_dofs =
3013  for (unsigned int comp=0; comp<n_components; ++comp)
3014  internal::check_vector_compatibility (*src[comp], *dof_info);
3015  Number *local_data [n_components];
3016  for (unsigned int comp=0; comp<n_components; ++comp)
3017  local_data[comp] =
3018  const_cast<Number *>(&values_dofs[comp][0][0]);
3019 
3020  // standard case where there are sufficiently many cells to fill all
3021  // vectors
3022  if (at_irregular_cell == false)
3023  {
3024  // check whether there is any constraint on the current cell
3025  if (indicators != indicators_end)
3026  {
3027  for ( ; indicators != indicators_end; ++indicators)
3028  {
3029  // run through values up to next constraint
3030  for (unsigned int j=0; j<indicators->first; ++j)
3031  for (unsigned int comp=0; comp<n_components; ++comp)
3032  operation.process_dof (dof_indices[j], *src[comp],
3033  local_data[comp][ind_local+j]);
3034 
3035  ind_local += indicators->first;
3036  dof_indices += indicators->first;
3037 
3038  // constrained case: build the local value as a linear
3039  // combination of the global value according to constraints
3040  Number value [n_components];
3041  for (unsigned int comp=0; comp<n_components; ++comp)
3042  operation.pre_constraints (local_data[comp][ind_local],
3043  value[comp]);
3044 
3045  const Number *data_val =
3046  matrix_info->constraint_pool_begin(indicators->second);
3047  const Number *end_pool =
3048  matrix_info->constraint_pool_end(indicators->second);
3049  for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3050  for (unsigned int comp=0; comp<n_components; ++comp)
3051  operation.process_constraint (*dof_indices, *data_val,
3052  *src[comp], value[comp]);
3053 
3054  for (unsigned int comp=0; comp<n_components; ++comp)
3055  operation.post_constraints (value[comp],
3056  local_data[comp][ind_local]);
3057 
3058  ind_local++;
3059  }
3060 
3061  // get the dof values past the last constraint
3062  for (; ind_local < n_local_dofs; ++dof_indices, ++ind_local)
3063  {
3064  for (unsigned int comp=0; comp<n_components; ++comp)
3065  operation.process_dof (*dof_indices, *src[comp],
3066  local_data[comp][ind_local]);
3067  }
3068  }
3069  else
3070  {
3071  // no constraint at all: compiler can unroll at least the
3072  // vectorization loop
3073  AssertDimension (dof_info->end_indices(cell)-dof_indices,
3074  static_cast<int>(n_local_dofs));
3075  for (unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
3076  for (unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
3077  for (unsigned int comp=0; comp<n_components; ++comp)
3078  operation.process_dof (dof_indices[j+v], *src[comp],
3079  local_data[comp][j+v]);
3080  }
3081  }
3082 
3083  // non-standard case: need to fill in zeros for those components that
3084  // are not present (a bit more expensive), but there is not more than
3085  // one such cell
3086  else
3087  {
3088  Assert (n_irreg_components_filled > 0, ExcInternalError());
3089  for ( ; indicators != indicators_end; ++indicators)
3090  {
3091  for (unsigned int j=0; j<indicators->first; ++j)
3092  {
3093  // non-constrained case: copy the data from the global
3094  // vector, src, to the local one, local_src.
3095  for (unsigned int comp=0; comp<n_components; ++comp)
3096  operation.process_dof (dof_indices[j], *src[comp],
3097  local_data[comp][ind_local]);
3098 
3099  // here we jump over all the components that are artificial
3100  ++ind_local;
3102  >= n_irreg_components_filled)
3103  {
3104  for (unsigned int comp=0; comp<n_components; ++comp)
3105  operation.process_empty (local_data[comp][ind_local]);
3106  ++ind_local;
3107  }
3108  }
3109  dof_indices += indicators->first;
3110 
3111  // constrained case: build the local value as a linear
3112  // combination of the global value according to constraint
3113  Number value [n_components];
3114  for (unsigned int comp=0; comp<n_components; ++comp)
3115  operation.pre_constraints (local_data[comp][ind_local],
3116  value[comp]);
3117 
3118  const Number *data_val =
3119  matrix_info->constraint_pool_begin(indicators->second);
3120  const Number *end_pool =
3121  matrix_info->constraint_pool_end(indicators->second);
3122 
3123  for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3124  for (unsigned int comp=0; comp<n_components; ++comp)
3125  operation.process_constraint (*dof_indices, *data_val,
3126  *src[comp], value[comp]);
3127 
3128  for (unsigned int comp=0; comp<n_components; ++comp)
3129  operation.post_constraints (value[comp],
3130  local_data[comp][ind_local]);
3131  ind_local++;
3133  >= n_irreg_components_filled)
3134  {
3135  for (unsigned int comp=0; comp<n_components; ++comp)
3136  operation.process_empty (local_data[comp][ind_local]);
3137  ++ind_local;
3138  }
3139  }
3140  for (; ind_local<n_local_dofs; ++dof_indices)
3141  {
3142  Assert (dof_indices != dof_info->end_indices(cell),
3143  ExcInternalError());
3144 
3145  // non-constrained case: copy the data from the global vector,
3146  // src, to the local one, local_dst.
3147  for (unsigned int comp=0; comp<n_components; ++comp)
3148  operation.process_dof (*dof_indices, *src[comp],
3149  local_data[comp][ind_local]);
3150  ++ind_local;
3152  >= n_irreg_components_filled)
3153  {
3154  for (unsigned int comp=0; comp<n_components; ++comp)
3155  operation.process_empty(local_data[comp][ind_local]);
3156  ++ind_local;
3157  }
3158  }
3159  }
3160  }
3161  else
3162  // case with vector-valued finite elements where all components are
3163  // included in one single vector. Assumption: first come all entries to
3164  // the first component, then all entries to the second one, and so
3165  // on. This is ensured by the way MatrixFree reads out the indices.
3166  {
3167  internal::check_vector_compatibility (*src[0], *dof_info);
3168  Assert (n_fe_components == n_components_, ExcNotImplemented());
3169  const unsigned int n_local_dofs =
3171  Number *local_data =
3172  const_cast<Number *>(&values_dofs[0][0][0]);
3173  if (at_irregular_cell == false)
3174  {
3175  // check whether there is any constraint on the current cell
3176  if (indicators != indicators_end)
3177  {
3178  for ( ; indicators != indicators_end; ++indicators)
3179  {
3180  // run through values up to next constraint
3181  for (unsigned int j=0; j<indicators->first; ++j)
3182  operation.process_dof (dof_indices[j], *src[0],
3183  local_data[ind_local+j]);
3184  ind_local += indicators->first;
3185  dof_indices += indicators->first;
3186 
3187  // constrained case: build the local value as a linear
3188  // combination of the global value according to constraints
3189  Number value;
3190  operation.pre_constraints (local_data[ind_local], value);
3191 
3192  const Number *data_val =
3193  matrix_info->constraint_pool_begin(indicators->second);
3194  const Number *end_pool =
3195  matrix_info->constraint_pool_end(indicators->second);
3196 
3197  for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3198  operation.process_constraint (*dof_indices, *data_val,
3199  *src[0], value);
3200 
3201  operation.post_constraints (value, local_data[ind_local]);
3202  ind_local++;
3203  }
3204 
3205  // get the dof values past the last constraint
3206  for (; ind_local<n_local_dofs; ++dof_indices, ++ind_local)
3207  operation.process_dof (*dof_indices, *src[0],
3208  local_data[ind_local]);
3209  Assert (dof_indices == dof_info->end_indices(cell),
3210  ExcInternalError());
3211  }
3212  else
3213  {
3214  // no constraint at all: compiler can unroll at least the
3215  // vectorization loop
3216  AssertDimension (dof_info->end_indices(cell)-dof_indices,
3217  static_cast<int>(n_local_dofs));
3218  for (unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
3219  for (unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
3220  operation.process_dof (dof_indices[j+v], *src[0],
3221  local_data[j+v]);
3222  }
3223  }
3224 
3225  // non-standard case: need to fill in zeros for those components that
3226  // are not present (a bit more expensive), but there is not more than
3227  // one such cell
3228  else
3229  {
3230  Assert (n_irreg_components_filled > 0, ExcInternalError());
3231  for ( ; indicators != indicators_end; ++indicators)
3232  {
3233  for (unsigned int j=0; j<indicators->first; ++j)
3234  {
3235  // non-constrained case: copy the data from the global
3236  // vector, src, to the local one, local_src.
3237  operation.process_dof (dof_indices[j], *src[0],
3238  local_data[ind_local]);
3239 
3240  // here we jump over all the components that are artificial
3241  ++ind_local;
3243  >= n_irreg_components_filled)
3244  {
3245  operation.process_empty (local_data[ind_local]);
3246  ++ind_local;
3247  }
3248  }
3249  dof_indices += indicators->first;
3250 
3251  // constrained case: build the local value as a linear
3252  // combination of the global value according to constraint
3253  Number value;
3254  operation.pre_constraints (local_data[ind_local], value);
3255 
3256  const Number *data_val =
3257  matrix_info->constraint_pool_begin(indicators->second);
3258  const Number *end_pool =
3259  matrix_info->constraint_pool_end(indicators->second);
3260 
3261  for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3262  operation.process_constraint (*dof_indices, *data_val,
3263  *src[0], value);
3264 
3265  operation.post_constraints (value, local_data[ind_local]);
3266  ind_local++;
3268  >= n_irreg_components_filled)
3269  {
3270  operation.process_empty (local_data[ind_local]);
3271  ++ind_local;
3272  }
3273  }
3274  for (; ind_local<n_local_dofs; ++dof_indices)
3275  {
3276  Assert (dof_indices != dof_info->end_indices(cell),
3277  ExcInternalError());
3278 
3279  // non-constrained case: copy the data from the global vector,
3280  // src, to the local one, local_dst.
3281  operation.process_dof (*dof_indices, *src[0],
3282  local_data[ind_local]);
3283  ++ind_local;
3285  >= n_irreg_components_filled)
3286  {
3287  operation.process_empty (local_data[ind_local]);
3288  ++ind_local;
3289  }
3290  }
3291  }
3292  }
3293 }
3294 
3295 
3296 
3297 template <int dim, int n_components_, typename Number>
3298 template<typename VectorType>
3299 inline
3300 void
3302 ::read_dof_values (const VectorType &src,
3303  const unsigned int first_index)
3304 {
3305  // select between block vectors and non-block vectors. Note that the number
3306  // of components is checked in the internal data
3307  typename internal::BlockVectorSelector<VectorType,
3308  IsBlockVector<VectorType>::value>::BaseVectorType *src_data[n_components];
3309  for (unsigned int d=0; d<n_components; ++d)
3310  src_data[d] = internal::BlockVectorSelector<VectorType, IsBlockVector<VectorType>::value>::get_vector_component(const_cast<VectorType &>(src), d+first_index);
3311 
3312  internal::VectorReader<Number> reader;
3313  read_write_operation (reader, src_data);
3314 
3315 #ifdef DEBUG
3316  dof_values_initialized = true;
3317 #endif
3318 }
3319 
3320 
3321 
3322 template <int dim, int n_components_, typename Number>
3323 template<typename VectorType>
3324 inline
3325 void
3327 ::read_dof_values_plain (const VectorType &src,
3328  const unsigned int first_index)
3329 {
3330  // select between block vectors and non-block vectors. Note that the number
3331  // of components is checked in the internal data
3332  const typename internal::BlockVectorSelector<VectorType,
3333  IsBlockVector<VectorType>::value>::BaseVectorType *src_data[n_components];
3334  for (unsigned int d=0; d<n_components; ++d)
3335  src_data[d] = internal::BlockVectorSelector<VectorType, IsBlockVector<VectorType>::value>::get_vector_component(const_cast<VectorType &>(src), d+first_index);
3336 
3337  read_dof_values_plain (src_data);
3338 }
3339 
3340 
3341 
3342 template <int dim, int n_components_, typename Number>
3343 template<typename VectorType>
3344 inline
3345 void
3347 ::distribute_local_to_global (VectorType &dst,
3348  const unsigned int first_index) const
3349 {
3350  Assert (dof_values_initialized==true,
3352 
3353  // select between block vectors and non-block vectors. Note that the number
3354  // of components is checked in the internal data
3355  typename internal::BlockVectorSelector<VectorType,
3356  IsBlockVector<VectorType>::value>::BaseVectorType *dst_data[n_components];
3357  for (unsigned int d=0; d<n_components; ++d)
3358  dst_data[d] = internal::BlockVectorSelector<VectorType, IsBlockVector<VectorType>::value>::get_vector_component(dst, d+first_index);
3359 
3360  internal::VectorDistributorLocalToGlobal<Number> distributor;
3361  read_write_operation (distributor, dst_data);
3362 }
3363 
3364 
3365 
3366 template <int dim, int n_components_, typename Number>
3367 template<typename VectorType>
3368 inline
3369 void
3371 ::set_dof_values (VectorType &dst,
3372  const unsigned int first_index) const
3373 {
3374  Assert (dof_values_initialized==true,
3376 
3377  // select between block vectors and non-block vectors. Note that the number
3378  // of components is checked in the internal data
3379  typename internal::BlockVectorSelector<VectorType,
3380  IsBlockVector<VectorType>::value>::BaseVectorType *dst_data[n_components];
3381  for (unsigned int d=0; d<n_components; ++d)
3382  dst_data[d] = internal::BlockVectorSelector<VectorType, IsBlockVector<VectorType>::value>::get_vector_component(dst, d+first_index);
3383 
3384  internal::VectorSetter<Number> setter;
3385  read_write_operation (setter, dst_data);
3386 }
3387 
3388 
3389 
3390 template <int dim, int n_components_, typename Number>
3391 template<typename VectorType>
3392 inline
3393 void
3395 ::read_dof_values_plain (const VectorType *src[])
3396 {
3397  // Case without MatrixFree initialization object
3398  if (matrix_info == 0)
3399  {
3400  internal::VectorReader<Number> reader;
3401  read_write_operation (reader, src);
3402  return;
3403  }
3404 
3405  // this is different from the other three operations because we do not use
3406  // constraints here, so this is a separate function.
3407  Assert (dof_info != 0, ExcNotInitialized());
3408  Assert (matrix_info->indices_initialized() == true,
3409  ExcNotInitialized());
3411  Assert (dof_info->store_plain_indices == true, ExcNotInitialized());
3412 
3413  // loop over all local dofs. ind_local holds local number on cell, index
3414  // iterates over the elements of index_local_to_global and dof_indices
3415  // points to the global indices stored in index_local_to_global
3416  const unsigned int *dof_indices = dof_info->begin_indices_plain(cell);
3417  const unsigned int dofs_per_cell = this->data->dofs_per_cell;
3418 
3419  const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
3420  const bool at_irregular_cell = n_irreg_components_filled > 0;
3421 
3422  // scalar case (or case when all components have the same degrees of freedom
3423  // and sit on a different vector each)
3424  if (n_fe_components == 1)
3425  {
3426  const unsigned int n_local_dofs =
3428  for (unsigned int comp=0; comp<n_components; ++comp)
3429  internal::check_vector_compatibility (*src[comp], *dof_info);
3430  Number *local_src_number [n_components];
3431  for (unsigned int comp=0; comp<n_components; ++comp)
3432  local_src_number[comp] = &values_dofs[comp][0][0];
3433 
3434  // standard case where there are sufficiently many cells to fill all
3435  // vectors
3436  if (at_irregular_cell == false)
3437  {
3438  for (unsigned int j=0; j<n_local_dofs; ++j)
3439  for (unsigned int comp=0; comp<n_components; ++comp)
3440  local_src_number[comp][j] =
3441  internal::vector_access (*src[comp], dof_indices[j]);
3442  }
3443 
3444  // non-standard case: need to fill in zeros for those components that
3445  // are not present (a bit more expensive), but there is not more than
3446  // one such cell
3447  else
3448  {
3449  Assert (n_irreg_components_filled > 0, ExcInternalError());
3450  for (unsigned int ind_local=0; ind_local<n_local_dofs;
3451  ++dof_indices)
3452  {
3453  // non-constrained case: copy the data from the global vector,
3454  // src, to the local one, local_dst.
3455  for (unsigned int comp=0; comp<n_components; ++comp)
3456  local_src_number[comp][ind_local] =
3457  internal::vector_access (*src[comp], *dof_indices);
3458  ++ind_local;
3459  while (ind_local % VectorizedArray<Number>::n_array_elements >= n_irreg_components_filled)
3460  {
3461  for (unsigned int comp=0; comp<n_components; ++comp)
3462  local_src_number[comp][ind_local] = 0.;
3463  ++ind_local;
3464  }
3465  }
3466  }
3467  }
3468  else
3469  // case with vector-valued finite elements where all components are
3470  // included in one single vector. Assumption: first come all entries to
3471  // the first component, then all entries to the second one, and so
3472  // on. This is ensured by the way MatrixFree reads out the indices.
3473  {
3474  internal::check_vector_compatibility (*src[0], *dof_info);
3475  Assert (n_fe_components == n_components_, ExcNotImplemented());
3476  const unsigned int n_local_dofs =
3478  Number *local_src_number = &values_dofs[0][0][0];
3479  if (at_irregular_cell == false)
3480  {
3481  for (unsigned int j=0; j<n_local_dofs; ++j)
3482  local_src_number[j] =
3483  internal::vector_access (*src[0], dof_indices[j]);
3484  }
3485 
3486  // non-standard case: need to fill in zeros for those components that
3487  // are not present (a bit more expensive), but there is not more than
3488  // one such cell
3489  else
3490  {
3491  Assert (n_irreg_components_filled > 0, ExcInternalError());
3492  for (unsigned int ind_local=0; ind_local<n_local_dofs; ++dof_indices)
3493  {
3494  // non-constrained case: copy the data from the global vector,
3495  // src, to the local one, local_dst.
3496  local_src_number[ind_local] =
3497  internal::vector_access (*src[0], *dof_indices);
3498  ++ind_local;
3499  while (ind_local % VectorizedArray<Number>::n_array_elements >= n_irreg_components_filled)
3500  {
3501  local_src_number[ind_local] = 0.;
3502  ++ind_local;
3503  }
3504  }
3505  }
3506  }
3507 
3508 #ifdef DEBUG
3509  dof_values_initialized = true;
3510 #endif
3511 }
3512 
3513 
3514 
3515 
3516 /*------------------------------ access to data fields ----------------------*/
3517 
3518 template <int dim, int n_components, typename Number>
3519 inline
3520 const std::vector<unsigned int> &
3523 {
3524  return data->lexicographic_numbering;
3525 }
3526 
3527 
3528 
3529 template <int dim, int n_components, typename Number>
3530 inline
3533 get_scratch_data() const
3534 {
3537  scratch_data);
3538 }
3539 
3540 
3541 
3542 template <int dim, int n_components, typename Number>
3543 inline
3546 begin_dof_values () const
3547 {
3548  return &values_dofs[0][0];
3549 }
3550 
3551 
3552 
3553 template <int dim, int n_components, typename Number>
3554 inline
3558 {
3559 #ifdef DEBUG
3560  dof_values_initialized = true;
3561 #endif
3562  return &values_dofs[0][0];
3563 }
3564 
3565 
3566 
3567 template <int dim, int n_components, typename Number>
3568 inline
3571 begin_values () const
3572 {
3574  ExcNotInitialized());
3575  return &values_quad[0][0];
3576 }
3577 
3578 
3579 
3580 template <int dim, int n_components, typename Number>
3581 inline
3584 begin_values ()
3585 {
3586 #ifdef DEBUG
3587  values_quad_submitted = true;
3588 #endif
3589  return &values_quad[0][0];
3590 }
3591 
3592 
3593 
3594 template <int dim, int n_components, typename Number>
3595 inline
3598 begin_gradients () const
3599 {
3601  ExcNotInitialized());
3602  return &gradients_quad[0][0][0];
3603 }
3604 
3605 
3606 
3607 template <int dim, int n_components, typename Number>
3608 inline
3611 begin_gradients ()
3612 {
3613 #ifdef DEBUG
3614  gradients_quad_submitted = true;
3615 #endif
3616  return &gradients_quad[0][0][0];
3617 }
3618 
3619 
3620 
3621 template <int dim, int n_components, typename Number>
3622 inline
3625 begin_hessians () const
3626 {
3628  return &hessians_quad[0][0][0];
3629 }
3630 
3631 
3632 
3633 template <int dim, int n_components, typename Number>
3634 inline
3637 begin_hessians ()
3638 {
3639  return &hessians_quad[0][0][0];
3640 }
3641 
3642 
3643 
3644 template <int dim, int n_components_, typename Number>
3645 inline
3648 ::get_dof_value (const unsigned int dof) const
3649 {
3650  AssertIndexRange (dof, this->data->dofs_per_cell);
3652  for (unsigned int comp=0; comp<n_components; comp++)
3653  return_value[comp] = this->values_dofs[comp][dof];
3654  return return_value;
3655 }
3656 
3657 
3658 
3659 template <int dim, int n_components_, typename Number>
3660 inline
3663 ::get_value (const unsigned int q_point) const
3664 {
3665  Assert (this->values_quad_initialized==true,
3667  AssertIndexRange (q_point, this->data->n_q_points);
3669  for (unsigned int comp=0; comp<n_components; comp++)
3670  return_value[comp] = this->values_quad[comp][q_point];
3671  return return_value;
3672 }
3673 
3674 
3675 
3676 template <int dim, int n_components_, typename Number>
3677 inline
3680 ::get_gradient (const unsigned int q_point) const
3681 {
3682  Assert (this->gradients_quad_initialized==true,
3684  AssertIndexRange (q_point, this->data->n_q_points);
3685 
3687 
3688  // Cartesian cell
3689  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3690  {
3691  for (unsigned int comp=0; comp<n_components; comp++)
3692  for (unsigned int d=0; d<dim; ++d)
3693  grad_out[comp][d] = (this->gradients_quad[comp][d][q_point] *
3694  cartesian_data[0][d]);
3695  }
3696  // cell with general/affine Jacobian
3697  else
3698  {
3700  this->cell_type == internal::MatrixFreeFunctions::general ?
3701  jacobian[q_point] : jacobian[0];
3702  for (unsigned int comp=0; comp<n_components; comp++)
3703  {
3704  for (unsigned int d=0; d<dim; ++d)
3705  {
3706  grad_out[comp][d] = (jac[d][0] *
3707  this->gradients_quad[comp][0][q_point]);
3708  for (unsigned int e=1; e<dim; ++e)
3709  grad_out[comp][d] += (jac[d][e] *
3710  this->gradients_quad[comp][e][q_point]);
3711  }
3712  }
3713  }
3714  return grad_out;
3715 }
3716 
3717 
3718 
3719 namespace internal
3720 {
3721  // compute tmp = hess_unit(u) * J^T. do this manually because we do not
3722  // store the lower diagonal because of symmetry
3723  template <typename Number>
3724  inline
3725  void
3726  hessian_unit_times_jac (const Tensor<2,1,VectorizedArray<Number> > &jac,
3727  const VectorizedArray<Number> *const hessians_quad[1],
3728  const unsigned int q_point,
3729  VectorizedArray<Number> (&tmp)[1][1])
3730  {
3731  tmp[0][0] = jac[0][0] * hessians_quad[0][q_point];
3732  }
3733 
3734  template <typename Number>
3735  inline
3736  void
3737  hessian_unit_times_jac (const Tensor<2,2,VectorizedArray<Number> > &jac,
3738  const VectorizedArray<Number> *const hessians_quad[3],
3739  const unsigned int q_point,
3740  VectorizedArray<Number> (&tmp)[2][2])
3741  {
3742  for (unsigned int d=0; d<2; ++d)
3743  {
3744  tmp[0][d] = (jac[d][0] * hessians_quad[0][q_point] +
3745  jac[d][1] * hessians_quad[2][q_point]);
3746  tmp[1][d] = (jac[d][0] * hessians_quad[2][q_point] +
3747  jac[d][1] * hessians_quad[1][q_point]);
3748  }
3749  }
3750 
3751  template <typename Number>
3752  inline
3753  void
3754  hessian_unit_times_jac (const Tensor<2,3,VectorizedArray<Number> > &jac,
3755  const VectorizedArray<Number> *const hessians_quad[6],
3756  const unsigned int q_point,
3757  VectorizedArray<Number> (&tmp)[3][3])
3758  {
3759  for (unsigned int d=0; d<3; ++d)
3760  {
3761  tmp[0][d] = (jac[d][0] * hessians_quad[0][q_point] +
3762  jac[d][1] * hessians_quad[3][q_point] +
3763  jac[d][2] * hessians_quad[4][q_point]);
3764  tmp[1][d] = (jac[d][0] * hessians_quad[3][q_point] +
3765  jac[d][1] * hessians_quad[1][q_point] +
3766  jac[d][2] * hessians_quad[5][q_point]);
3767  tmp[2][d] = (jac[d][0] * hessians_quad[4][q_point] +
3768  jac[d][1] * hessians_quad[5][q_point] +
3769  jac[d][2] * hessians_quad[2][q_point]);
3770  }
3771  }
3772 }
3773 
3774 
3775 
3776 template <int dim, int n_components_, typename Number>
3777 inline
3780 ::get_hessian (const unsigned int q_point) const
3781 {
3782  Assert (this->hessians_quad_initialized==true,
3784  AssertIndexRange (q_point, this->data->n_q_points);
3785 
3787 
3788  // Cartesian cell
3789  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3790  {
3791  const Tensor<1,dim,VectorizedArray<Number> > &jac = cartesian_data[0];
3792  for (unsigned int comp=0; comp<n_components; comp++)
3793  for (unsigned int d=0; d<dim; ++d)
3794  {
3795  hessian_out[comp][d][d] = (this->hessians_quad[comp][d][q_point] *
3796  jac[d] * jac[d]);
3797  switch (dim)
3798  {
3799  case 1:
3800  break;
3801  case 2:
3802  hessian_out[comp][0][1] = (this->hessians_quad[comp][2][q_point] *
3803  jac[0] * jac[1]);
3804  break;
3805  case 3:
3806  hessian_out[comp][0][1] = (this->hessians_quad[comp][3][q_point] *
3807  jac[0] * jac[1]);
3808  hessian_out[comp][0][2] = (this->hessians_quad[comp][4][q_point] *
3809  jac[0] * jac[2]);
3810  hessian_out[comp][1][2] = (this->hessians_quad[comp][5][q_point] *
3811  jac[1] * jac[2]);
3812  break;
3813  default:
3814  Assert (false, ExcNotImplemented());
3815  }
3816  for (unsigned int e=d+1; e<dim; ++e)
3817  hessian_out[comp][e][d] = hessian_out[comp][d][e];
3818  }
3819  }
3820  // cell with general Jacobian
3821  else if (this->cell_type == internal::MatrixFreeFunctions::general)
3822  {
3823  Assert (this->mapping_info->second_derivatives_initialized == true,
3824  ExcNotInitialized());
3825  const Tensor<2,dim,VectorizedArray<Number> > &jac = jacobian[q_point];
3826  const Tensor<2,dim,VectorizedArray<Number> > &jac_grad = jacobian_grad[q_point];
3827  const Tensor<1,(dim>1?dim*(dim-1)/2:1),
3829  & jac_grad_UT = jacobian_grad_upper[q_point];
3830  for (unsigned int comp=0; comp<n_components; comp++)
3831  {
3832  // compute laplacian before the gradient because it needs to access
3833  // unscaled gradient data
3834  VectorizedArray<Number> tmp[dim][dim];
3835  internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3836  q_point, tmp);
3837 
3838  // compute first part of hessian, J * tmp = J * hess_unit(u) * J^T
3839  for (unsigned int d=0; d<dim; ++d)
3840  for (unsigned int e=d; e<dim; ++e)
3841  {
3842  hessian_out[comp][d][e] = jac[d][0] * tmp[0][e];
3843  for (unsigned int f=1; f<dim; ++f)
3844  hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3845  }
3846 
3847  // add diagonal part of J' * grad(u)
3848  for (unsigned int d=0; d<dim; ++d)
3849  for (unsigned int e=0; e<dim; ++e)
3850  hessian_out[comp][d][d] += (jac_grad[d][e] *
3851  this->gradients_quad[comp][e][q_point]);
3852 
3853  // add off-diagonal part of J' * grad(u)
3854  for (unsigned int d=0, count=0; d<dim; ++d)
3855  for (unsigned int e=d+1; e<dim; ++e, ++count)
3856  for (unsigned int f=0; f<dim; ++f)
3857  hessian_out[comp][d][e] += (jac_grad_UT[count][f] *
3858  this->gradients_quad[comp][f][q_point]);
3859 
3860  // take symmetric part
3861  for (unsigned int d=0; d<dim; ++d)
3862  for (unsigned int e=d+1; e<dim; ++e)
3863  hessian_out[comp][e][d] = hessian_out[comp][d][e];
3864  }
3865  }
3866  // cell with general Jacobian, but constant within the cell
3867  else // if (this->cell_type == internal::MatrixFreeFunctions::affine)
3868  {
3869  const Tensor<2,dim,VectorizedArray<Number> > &jac = jacobian[0];
3870  for (unsigned int comp=0; comp<n_components; comp++)
3871  {
3872  // compute laplacian before the gradient because it needs to access
3873  // unscaled gradient data
3874  VectorizedArray<Number> tmp[dim][dim];
3875  internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3876  q_point, tmp);
3877 
3878  // compute first part of hessian, J * tmp = J * hess_unit(u) * J^T
3879  for (unsigned int d=0; d<dim; ++d)
3880  for (unsigned int e=d; e<dim; ++e)
3881  {
3882  hessian_out[comp][d][e] = jac[d][0] * tmp[0][e];
3883  for (unsigned int f=1; f<dim; ++f)
3884  hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3885  }
3886 
3887  // no J' * grad(u) part here because the Jacobian is constant
3888  // throughout the cell and hence, its derivative is zero
3889 
3890  // take symmetric part
3891  for (unsigned int d=0; d<dim; ++d)
3892  for (unsigned int e=d+1; e<dim; ++e)
3893  hessian_out[comp][e][d] = hessian_out[comp][d][e];
3894  }
3895  }
3897 }
3898 
3899 
3900 
3901 template <int dim, int n_components_, typename Number>
3902 inline
3905 ::get_hessian_diagonal (const unsigned int q_point) const
3906 {
3907  Assert (this->hessians_quad_initialized==true,
3909  AssertIndexRange (q_point, this->data->n_q_points);
3910 
3912 
3913  // Cartesian cell
3914  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3915  {
3916  const Tensor<1,dim,VectorizedArray<Number> > &jac = cartesian_data[0];
3917  for (unsigned int comp=0; comp<n_components; comp++)
3918  for (unsigned int d=0; d<dim; ++d)
3919  hessian_out[comp][d] = (this->hessians_quad[comp][d][q_point] *
3920  jac[d] * jac[d]);
3921  }
3922  // cell with general Jacobian
3923  else if (this->cell_type == internal::MatrixFreeFunctions::general)
3924  {
3925  Assert (this->mapping_info->second_derivatives_initialized == true,
3926  ExcNotInitialized());
3927  const Tensor<2,dim,VectorizedArray<Number> > &jac = jacobian[q_point];
3928  const Tensor<2,dim,VectorizedArray<Number> > &jac_grad = jacobian_grad[q_point];
3929  for (unsigned int comp=0; comp<n_components; comp++)
3930  {
3931  // compute laplacian before the gradient because it needs to access
3932  // unscaled gradient data
3933  VectorizedArray<Number> tmp[dim][dim];
3934  internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3935  q_point, tmp);
3936 
3937  // compute only the trace part of hessian, J * tmp = J *
3938  // hess_unit(u) * J^T
3939  for (unsigned int d=0; d<dim; ++d)
3940  {
3941  hessian_out[comp][d] = jac[d][0] * tmp[0][d];
3942  for (unsigned int f=1; f<dim; ++f)
3943  hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3944  }
3945 
3946  for (unsigned int d=0; d<dim; ++d)
3947  for (unsigned int e=0; e<dim; ++e)
3948  hessian_out[comp][d] += (jac_grad[d][e] *
3949  this->gradients_quad[comp][e][q_point]);
3950  }
3951  }
3952  // cell with general Jacobian, but constant within the cell
3953  else // if (this->cell_type == internal::MatrixFreeFunctions::affine)
3954  {
3955  const Tensor<2,dim,VectorizedArray<Number> > &jac = jacobian[0];
3956  for (unsigned int comp=0; comp<n_components; comp++)
3957  {
3958  // compute laplacian before the gradient because it needs to access
3959  // unscaled gradient data
3960  VectorizedArray<Number> tmp[dim][dim];
3961  internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3962  q_point, tmp);
3963 
3964  // compute only the trace part of hessian, J * tmp = J *
3965  // hess_unit(u) * J^T
3966  for (unsigned int d=0; d<dim; ++d)
3967  {
3968  hessian_out[comp][d] = jac[d][0] * tmp[0][d];
3969  for (unsigned int f=1; f<dim; ++f)
3970  hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3971  }
3972  }
3973  }
3974  return hessian_out;
3975 }
3976 
3977 
3978 
3979 template <int dim, int n_components_, typename Number>
3980 inline
3983 ::get_laplacian (const unsigned int q_point) const
3984 {
3985  Assert (this->hessians_quad_initialized==true,
3987  AssertIndexRange (q_point, this->data->n_q_points);
3990  = get_hessian_diagonal(q_point);
3991  for (unsigned int comp=0; comp<n_components; ++comp)
3992  {
3993  laplacian_out[comp] = hess_diag[comp][0];
3994  for (unsigned int d=1; d<dim; ++d)
3995  laplacian_out[comp] += hess_diag[comp][d];
3996  }
3997  return laplacian_out;
3998 }
3999 
4000 
4001 
4002 template <int dim, int n_components_, typename Number>
4003 inline
4004 void
4006 ::submit_dof_value (const Tensor<1,n_components_,VectorizedArray<Number> > val_in,
4007  const unsigned int dof)
4008 {
4009 #ifdef DEBUG
4010  this->dof_values_initialized = true;
4011 #endif
4012  AssertIndexRange (dof, this->data->dofs_per_cell);
4013  for (unsigned int comp=0; comp<n_components; comp++)
4014  this->values_dofs[comp][dof] = val_in[comp];
4015 }
4016 
4017 
4018 
4019 template <int dim, int n_components_, typename Number>
4020 inline
4021 void
4023 ::submit_value (const Tensor<1,n_components_,VectorizedArray<Number> > val_in,
4024  const unsigned int q_point)
4025 {
4026 #ifdef DEBUG
4028  AssertIndexRange (q_point, this->data->n_q_points);
4029  this->values_quad_submitted = true;
4030 #endif
4031  if (this->cell_type == internal::MatrixFreeFunctions::general)
4032  {
4033  const VectorizedArray<Number> JxW = J_value[q_point];
4034  for (unsigned int comp=0; comp<n_components; ++comp)
4035  this->values_quad[comp][q_point] = val_in[comp] * JxW;
4036  }
4037  else //if (this->cell_type < internal::MatrixFreeFunctions::general)
4038  {
4039  const VectorizedArray<Number> JxW = J_value[0] * quadrature_weights[q_point];
4040  for (unsigned int comp=0; comp<n_components; ++comp)
4041  this->values_quad[comp][q_point] = val_in[comp] * JxW;
4042  }
4043 }
4044 
4045 
4046 
4047 template <int dim, int n_components_, typename Number>
4048 inline
4049 void
4051 ::submit_gradient (const Tensor<1,n_components_,
4052  Tensor<1,dim,VectorizedArray<Number> > >grad_in,
4053  const unsigned int q_point)
4054 {
4055 #ifdef DEBUG
4057  AssertIndexRange (q_point, this->data->n_q_points);
4058  this->gradients_quad_submitted = true;
4059 #endif
4060  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4061  {
4062  const VectorizedArray<Number> JxW = J_value[0] * quadrature_weights[q_point];
4063  for (unsigned int comp=0; comp<n_components; comp++)
4064  for (unsigned int d=0; d<dim; ++d)
4065  this->gradients_quad[comp][d][q_point] = (grad_in[comp][d] *
4066  cartesian_data[0][d] * JxW);
4067  }
4068  else
4069  {
4071  this->cell_type == internal::MatrixFreeFunctions::general ?
4072  jacobian[q_point] : jacobian[0];
4073  const VectorizedArray<Number> JxW =
4074  this->cell_type == internal::MatrixFreeFunctions::general ?
4075  J_value[q_point] : J_value[0] * quadrature_weights[q_point];
4076  for (unsigned int comp=0; comp<n_components; ++comp)
4077  for (unsigned int d=0; d<dim; ++d)
4078  {
4079  VectorizedArray<Number> new_val = jac[0][d] * grad_in[comp][0];
4080  for (unsigned int e=1; e<dim; ++e)
4081  new_val += (jac[e][d] * grad_in[comp][e]);
4082  this->gradients_quad[comp][d][q_point] = new_val * JxW;
4083  }
4084  }
4085 }
4086 
4087 
4088 
4089 template <int dim, int n_components_, typename Number>
4090 inline
4093 ::integrate_value () const
4094 {
4095 #ifdef DEBUG
4097  Assert (this->values_quad_submitted == true,
4099 #endif
4101  for (unsigned int comp=0; comp<n_components; ++comp)
4102  return_value[comp] = this->values_quad[comp][0];
4103  const unsigned int n_q_points = this->data->n_q_points;
4104  for (unsigned int q=1; q<n_q_points; ++q)
4105  for (unsigned int comp=0; comp<n_components; ++comp)
4106  return_value[comp] += this->values_quad[comp][q];
4107  return (return_value);
4108 }
4109 
4110 
4111 
4112 /*----------------------- FEEvaluationAccess --------------------------------*/
4113 
4114 
4115 template <int dim, int n_components_, typename Number>
4116 inline
4119  const unsigned int fe_no,
4120  const unsigned int quad_no_in,
4121  const unsigned int fe_degree,
4122  const unsigned int n_q_points)
4123  :
4124  FEEvaluationBase <dim,n_components_,Number>
4125  (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4126 {}
4127 
4128 
4129 
4130 template <int dim, int n_components_, typename Number>
4131 template <int n_components_other>
4132 inline
4134 ::FEEvaluationAccess (const Mapping<dim> &mapping,
4135  const FiniteElement<dim> &fe,
4136  const Quadrature<1> &quadrature,
4137  const UpdateFlags update_flags,
4138  const unsigned int first_selected_component,
4140  :
4141  FEEvaluationBase <dim,n_components_,Number>(mapping, fe, quadrature, update_flags,
4142  first_selected_component, other)
4143 {}
4144 
4145 
4146 
4147 template <int dim, int n_components_, typename Number>
4148 inline
4151  :
4152  FEEvaluationBase <dim,n_components_,Number>(other)
4153 {}
4154 
4155 
4156 
4157 template <int dim, int n_components_, typename Number>
4158 inline
4162 {
4164  return *this;
4165 }
4166 
4167 
4168 
4169 /*-------------------- FEEvaluationAccess scalar ----------------------------*/
4170 
4171 
4172 template <int dim, typename Number>
4173 inline
4176  const unsigned int fe_no,
4177  const unsigned int quad_no_in,
4178  const unsigned int fe_degree,
4179  const unsigned int n_q_points)
4180  :
4181  FEEvaluationBase <dim,1,Number>
4182  (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4183 {}
4184 
4185 
4186 
4187 template <int dim, typename Number>
4188 template <int n_components_other>
4189 inline
4191 ::FEEvaluationAccess (const Mapping<dim> &mapping,
4192  const FiniteElement<dim> &fe,
4193  const Quadrature<1> &quadrature,
4194  const UpdateFlags update_flags,
4195  const unsigned int first_selected_component,
4197  :
4198  FEEvaluationBase <dim,1,Number> (mapping, fe, quadrature, update_flags,
4199  first_selected_component, other)
4200 {}
4201 
4202 
4203 
4204 template <int dim, typename Number>
4205 inline
4208  :
4209  FEEvaluationBase <dim,1,Number>(other)
4210 {}
4211 
4212 
4213 
4214 template <int dim, typename Number>
4215 inline
4219 {
4221  return *this;
4222 }
4223 
4224 
4225 
4226 template <int dim, typename Number>
4227 inline
4230 ::get_dof_value (const unsigned int dof) const
4231 {
4232  AssertIndexRange (dof, this->data->dofs_per_cell);
4233  return this->values_dofs[0][dof];
4234 }
4235 
4236 
4237 
4238 template <int dim, typename Number>
4239 inline
4242 ::get_value (const unsigned int q_point) const
4243 {
4244  Assert (this->values_quad_initialized==true,
4246  AssertIndexRange (q_point, this->data->n_q_points);
4247  return this->values_quad[0][q_point];
4248 }
4249 
4250 
4251 
4252 template <int dim, typename Number>
4253 inline
4256 ::get_gradient (const unsigned int q_point) const
4257 {
4258  // could use the base class gradient, but that involves too many inefficient
4259  // initialization operations on tensors
4260 
4261  Assert (this->gradients_quad_initialized==true,
4263  AssertIndexRange (q_point, this->data->n_q_points);
4264 
4266 
4267  // Cartesian cell
4268  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4269  {
4270  for (unsigned int d=0; d<dim; ++d)
4271  grad_out[d] = (this->gradients_quad[0][d][q_point] *
4272  this->cartesian_data[0][d]);
4273  }
4274  // cell with general/constant Jacobian
4275  else
4276  {
4278  this->cell_type == internal::MatrixFreeFunctions::general ?
4279  this->jacobian[q_point] : this->jacobian[0];
4280  for (unsigned int d=0; d<dim; ++d)
4281  {
4282  grad_out[d] = (jac[d][0] * this->gradients_quad[0][0][q_point]);
4283  for (unsigned int e=1; e<dim; ++e)
4284  grad_out[d] += (jac[d][e] * this->gradients_quad[0][e][q_point]);
4285  }
4286  }
4287  return grad_out;
4288 }
4289 
4290 
4291 
4292 template <int dim, typename Number>
4293 inline
4296 ::get_hessian (const unsigned int q_point) const
4297 {
4298  return BaseClass::get_hessian(q_point)[0];
4299 }
4300 
4301 
4302 
4303 template <int dim, typename Number>
4304 inline
4307 ::get_hessian_diagonal (const unsigned int q_point) const
4308 {
4309  return BaseClass::get_hessian_diagonal(q_point)[0];
4310 }
4311 
4312 
4313 
4314 template <int dim, typename Number>
4315 inline
4318 ::get_laplacian (const unsigned int q_point) const
4319 {
4320  return BaseClass::get_laplacian(q_point)[0];
4321 }
4322 
4323 
4324 
4325 template <int dim, typename Number>
4326 inline
4327 void
4330  const unsigned int dof)
4331 {
4332 #ifdef DEBUG
4333  this->dof_values_initialized = true;
4334  AssertIndexRange (dof, this->data->dofs_per_cell);
4335 #endif
4336  this->values_dofs[0][dof] = val_in;
4337 }
4338 
4339 
4340 
4341 template <int dim, typename Number>
4342 inline
4343 void
4346  const unsigned int q_point)
4347 {
4348 #ifdef DEBUG
4350  AssertIndexRange (q_point, this->data->n_q_points);
4351  this->values_quad_submitted = true;
4352 #endif
4353  if (this->cell_type == internal::MatrixFreeFunctions::general)
4354  {
4355  const VectorizedArray<Number> JxW = this->J_value[q_point];
4356  this->values_quad[0][q_point] = val_in * JxW;
4357  }
4358  else //if (this->cell_type < internal::MatrixFreeFunctions::general)
4359  {
4360  const VectorizedArray<Number> JxW = this->J_value[0] * this->quadrature_weights[q_point];
4361  this->values_quad[0][q_point] = val_in * JxW;
4362  }
4363 }
4364 
4365 
4366 
4367 template <int dim, typename Number>
4368 inline
4369 void
4371 ::submit_gradient (const Tensor<1,dim,VectorizedArray<Number> > grad_in,
4372  const unsigned int q_point)
4373 {
4374 #ifdef DEBUG
4376  AssertIndexRange (q_point, this->data->n_q_points);
4377  this->gradients_quad_submitted = true;
4378 #endif
4379  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4380  {
4381  const VectorizedArray<Number> JxW = this->J_value[0] * this->quadrature_weights[q_point];
4382  for (unsigned int d=0; d<dim; ++d)
4383  this->gradients_quad[0][d][q_point] = (grad_in[d] *
4384  this->cartesian_data[0][d] *
4385  JxW);
4386  }
4387  // general/affine cell type
4388  else
4389  {
4391  this->cell_type == internal::MatrixFreeFunctions::general ?
4392  this->jacobian[q_point] : this->jacobian[0];
4393  const VectorizedArray<Number> JxW =
4394  this->cell_type == internal::MatrixFreeFunctions::general ?
4395  this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4396  for (unsigned int d=0; d<dim; ++d)
4397  {
4398  VectorizedArray<Number> new_val = jac[0][d] * grad_in[0];
4399  for (unsigned int e=1; e<dim; ++e)
4400  new_val += jac[e][d] * grad_in[e];
4401  this->gradients_quad[0][d][q_point] = new_val * JxW;
4402  }
4403  }
4404 }
4405 
4406 
4407 
4408 template <int dim, typename Number>
4409 inline
4412 ::integrate_value () const
4413 {
4414  return BaseClass::integrate_value()[0];
4415 }
4416 
4417 
4418 
4419 
4420 /*----------------- FEEvaluationAccess vector-valued ------------------------*/
4421 
4422 
4423 template <int dim, typename Number>
4424 inline
4427  const unsigned int fe_no,
4428  const unsigned int quad_no_in,
4429  const unsigned int fe_degree,
4430  const unsigned int n_q_points)
4431  :
4432  FEEvaluationBase <dim,dim,Number>
4433  (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4434 {}
4435 
4436 
4437 
4438 template <int dim, typename Number>
4439 template <int n_components_other>
4440 inline
4442 ::FEEvaluationAccess (const Mapping<dim> &mapping,
4443  const FiniteElement<dim> &fe,
4444  const Quadrature<1> &quadrature,
4445  const UpdateFlags update_flags,
4446  const unsigned int first_selected_component,
4448  :
4449  FEEvaluationBase <dim,dim,Number> (mapping, fe, quadrature, update_flags,
4450  first_selected_component, other)
4451 {}
4452 
4453 
4454 
4455 template <int dim, typename Number>
4456 inline
4459  :
4460  FEEvaluationBase <dim,dim,Number>(other)
4461 {}
4462 
4463 
4464 
4465 template <int dim, typename Number>
4466 inline
4470 {
4472  return *this;
4473 }
4474 
4475 
4476 
4477 template <int dim, typename Number>
4478 inline
4481 ::get_gradient (const unsigned int q_point) const
4482 {
4483  return BaseClass::get_gradient (q_point);
4484 }
4485 
4486 
4487 
4488 template <int dim, typename Number>
4489 inline
4492 ::get_divergence (const unsigned int q_point) const
4493 {
4494  Assert (this->gradients_quad_initialized==true,
4496  AssertIndexRange (q_point, this->data->n_q_points);
4497 
4498  VectorizedArray<Number> divergence;
4499 
4500  // Cartesian cell
4501  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4502  {
4503  divergence = (this->gradients_quad[0][0][q_point] *
4504  this->cartesian_data[0][0]);
4505  for (unsigned int d=1; d<dim; ++d)
4506  divergence += (this->gradients_quad[d][d][q_point] *
4507  this->cartesian_data[0][d]);
4508  }
4509  // cell with general/constant Jacobian
4510  else
4511  {
4513  this->cell_type == internal::MatrixFreeFunctions::general ?
4514  this->jacobian[q_point] : this->jacobian[0];
4515  divergence = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4516  for (unsigned int e=1; e<dim; ++e)
4517  divergence += (jac[0][e] * this->gradients_quad[0][e][q_point]);
4518  for (unsigned int d=1; d<dim; ++d)
4519  for (unsigned int e=0; e<dim; ++e)
4520  divergence += (jac[d][e] * this->gradients_quad[d][e][q_point]);
4521  }
4522  return divergence;
4523 }
4524 
4525 
4526 
4527 template <int dim, typename Number>
4528 inline
4531 ::get_symmetric_gradient (const unsigned int q_point) const
4532 {
4533  // copy from generic function into dim-specialization function
4534  const Tensor<2,dim,VectorizedArray<Number> > grad = get_gradient(q_point);
4535  VectorizedArray<Number> symmetrized [(dim*dim+dim)/2];
4537  for (unsigned int d=0; d<dim; ++d)
4538  symmetrized[d] = grad[d][d];
4539  switch (dim)
4540  {
4541  case 1:
4542  break;
4543  case 2:
4544  symmetrized[2] = grad[0][1] + grad[1][0];
4545  symmetrized[2] *= half;
4546  break;
4547  case 3:
4548  symmetrized[3] = grad[0][1] + grad[1][0];
4549  symmetrized[3] *= half;
4550  symmetrized[4] = grad[0][2] + grad[2][0];
4551  symmetrized[4] *= half;
4552  symmetrized[5] = grad[1][2] + grad[2][1];
4553  symmetrized[5] *= half;
4554  break;
4555  default:
4556  Assert (false, ExcNotImplemented());
4557  }
4558  return SymmetricTensor<2,dim,VectorizedArray<Number> > (symmetrized);
4559 }
4560 
4561 
4562 
4563 template <int dim, typename Number>
4564 inline
4567 ::get_curl (const unsigned int q_point) const
4568 {
4569  // copy from generic function into dim-specialization function
4570  const Tensor<2,dim,VectorizedArray<Number> > grad = get_gradient(q_point);
4572  switch (dim)
4573  {
4574  case 1:
4575  Assert (false,
4576  ExcMessage("Computing the curl in 1d is not a useful operation"));
4577  break;
4578  case 2:
4579  curl[0] = grad[1][0] - grad[0][1];
4580  break;
4581  case 3:
4582  curl[0] = grad[2][1] - grad[1][2];
4583  curl[1] = grad[0][2] - grad[2][0];
4584  curl[2] = grad[1][0] - grad[0][1];
4585  break;
4586  default:
4587  Assert (false, ExcNotImplemented());
4588  }
4589  return curl;
4590 }
4591 
4592 
4593 
4594 template <int dim, typename Number>
4595 inline
4598 ::get_hessian_diagonal (const unsigned int q_point) const
4599 {
4600  Assert (this->hessians_quad_initialized==true,
4602  AssertIndexRange (q_point, this->data->n_q_points);
4603 
4604  return BaseClass::get_hessian_diagonal (q_point);
4605 }
4606 
4607 
4608 
4609 template <int dim, typename Number>
4610 inline
4613 ::get_hessian (const unsigned int q_point) const
4614 {
4615  Assert (this->hessians_quad_initialized==true,
4617  AssertIndexRange (q_point, this->data->n_q_points);
4618  return BaseClass::get_hessian(q_point);
4619 }
4620 
4621 
4622 
4623 template <int dim, typename Number>
4624 inline
4625 void
4627 ::submit_gradient (const Tensor<2,dim,VectorizedArray<Number> > grad_in,
4628  const unsigned int q_point)
4629 {
4630  BaseClass::submit_gradient (grad_in, q_point);
4631 }
4632 
4633 
4634 
4635 template <int dim, typename Number>
4636 inline
4637 void
4640  grad_in,
4641  const unsigned int q_point)
4642 {
4643  BaseClass::submit_gradient(grad_in, q_point);
4644 }
4645 
4646 
4647 
4648 template <int dim, typename Number>
4649 inline
4650 void
4653  const unsigned int q_point)
4654 {
4655 #ifdef DEBUG
4657  AssertIndexRange (q_point, this->data->n_q_points);
4658  this->gradients_quad_submitted = true;
4659 #endif
4660  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4661  {
4662  const VectorizedArray<Number> fac = this->J_value[0] *
4663  this->quadrature_weights[q_point] * div_in;
4664  for (unsigned int d=0; d<dim; ++d)
4665  {
4666  this->gradients_quad[d][d][q_point] = (fac *
4667  this->cartesian_data[0][d]);
4668  for (unsigned int e=d+1; e<dim; ++e)
4669  {
4670  this->gradients_quad[d][e][q_point] = VectorizedArray<Number>();
4671  this->gradients_quad[e][d][q_point] = VectorizedArray<Number>();
4672  }
4673  }
4674  }
4675  else
4676  {
4678  this->cell_type == internal::MatrixFreeFunctions::general ?
4679  this->jacobian[q_point] : this->jacobian[0];
4680  const VectorizedArray<Number> fac =
4681  (this->cell_type == internal::MatrixFreeFunctions::general ?
4682  this->J_value[q_point] : this->J_value[0] *
4683  this->quadrature_weights[q_point]) * div_in;
4684  for (unsigned int d=0; d<dim; ++d)
4685  {
4686  for (unsigned int e=0; e<dim; ++e)
4687  this->gradients_quad[d][e][q_point] = jac[d][e] * fac;
4688  }
4689  }
4690 }
4691 
4692 
4693 
4694 template <int dim, typename Number>
4695 inline
4696 void
4699  sym_grad,
4700  const unsigned int q_point)
4701 {
4702  // could have used base class operator, but that involves some overhead
4703  // which is inefficient. it is nice to have the symmetric tensor because
4704  // that saves some operations
4705 #ifdef DEBUG
4707  AssertIndexRange (q_point, this->data->n_q_points);
4708  this->gradients_quad_submitted = true;
4709 #endif
4710  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4711  {
4712  const VectorizedArray<Number> JxW = this->J_value[0] * this->quadrature_weights[q_point];
4713  for (unsigned int d=0; d<dim; ++d)
4714  this->gradients_quad[d][d][q_point] = (sym_grad.access_raw_entry(d) *
4715  JxW *
4716  this->cartesian_data[0][d]);
4717  for (unsigned int e=0, counter=dim; e<dim; ++e)
4718  for (unsigned int d=e+1; d<dim; ++d, ++counter)
4719  {
4720  const VectorizedArray<Number> value = sym_grad.access_raw_entry(counter) * JxW;
4721  this->gradients_quad[e][d][q_point] = (value *
4722  this->cartesian_data[0][d]);
4723  this->gradients_quad[d][e][q_point] = (value *
4724  this->cartesian_data[0][e]);
4725  }
4726  }
4727  // general/affine cell type
4728  else
4729  {
4730  const VectorizedArray<Number> JxW =
4731  this->cell_type == internal::MatrixFreeFunctions::general ?
4732  this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4734  this->cell_type == internal::MatrixFreeFunctions::general ?
4735  this->jacobian[q_point] : this->jacobian[0];
4736  VectorizedArray<Number> weighted [dim][dim];
4737  for (unsigned int i=0; i<dim; ++i)
4738  weighted[i][i] = sym_grad.access_raw_entry(i) * JxW;
4739  for (unsigned int i=0, counter=dim; i<dim; ++i)
4740  for (unsigned int j=i+1; j<dim; ++j, ++counter)
4741  {
4742  const VectorizedArray<Number> value = sym_grad.access_raw_entry(counter) * JxW;
4743  weighted[i][j] = value;
4744  weighted[j][i] = value;
4745  }
4746  for (unsigned int comp=0; comp<dim; ++comp)
4747  for (unsigned int d=0; d<dim; ++d)
4748  {
4749  VectorizedArray<Number> new_val = jac[0][d] * weighted[comp][0];
4750  for (unsigned int e=1; e<dim; ++e)
4751  new_val += jac[e][d] * weighted[comp][e];
4752  this->gradients_quad[comp][d][q_point] = new_val;
4753  }
4754  }
4755 }
4756 
4757 
4758 
4759 template <int dim, typename Number>
4760 inline
4761 void
4763 ::submit_curl (const Tensor<1,dim==2?1:dim,VectorizedArray<Number> > curl,
4764  const unsigned int q_point)
4765 {
4767  switch (dim)
4768  {
4769  case 1:
4770  Assert (false,
4771  ExcMessage("Testing by the curl in 1d is not a useful operation"));
4772  break;
4773  case 2:
4774  grad[1][0] = curl[0];
4775  grad[0][1] = -curl[0];
4776  break;
4777  case 3:
4778  grad[2][1] = curl[0];
4779  grad[1][2] = -curl[0];
4780  grad[0][2] = curl[1];
4781  grad[2][0] = -curl[1];
4782  grad[1][0] = curl[2];
4783  grad[0][1] = -curl[2];
4784  break;
4785  default:
4786  Assert (false, ExcNotImplemented());
4787  }
4788  submit_gradient (grad, q_point);
4789 }
4790 
4791 
4792 /*-------------------- FEEvaluationAccess scalar for 1d ----------------------------*/
4793 
4794 
4795 template <typename Number>
4796 inline
4799  const unsigned int fe_no,
4800  const unsigned int quad_no_in,
4801  const unsigned int fe_degree,
4802  const unsigned int n_q_points)
4803  :
4804  FEEvaluationBase <1,1,Number>
4805  (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4806 {}
4807 
4808 
4809 
4810 template <typename Number>
4811 template <int n_components_other>
4812 inline
4814 ::FEEvaluationAccess (const Mapping<1> &mapping,
4815  const FiniteElement<1> &fe,
4816  const Quadrature<1> &quadrature,
4817  const UpdateFlags update_flags,
4818  const unsigned int first_selected_component,
4820  :
4821  FEEvaluationBase <1,1,Number> (mapping, fe, quadrature, update_flags,
4822  first_selected_component, other)
4823 {}
4824 
4825 
4826 
4827 template <typename Number>
4828 inline
4831  :
4832  FEEvaluationBase <1,1,Number>(other)
4833 {}
4834 
4835 
4836 
4837 template <typename Number>
4838 inline
4842 {
4844  return *this;
4845 }
4846 
4847 
4848 
4849 template <typename Number>
4850 inline
4853 ::get_dof_value (const unsigned int dof) const
4854 {
4855  AssertIndexRange (dof, this->data->dofs_per_cell);
4856  return this->values_dofs[0][dof];
4857 }
4858 
4859 
4860 
4861 template <typename Number>
4862 inline
4865 ::get_value (const unsigned int q_point) const
4866 {
4867  Assert (this->values_quad_initialized==true,
4869  AssertIndexRange (q_point, this->data->n_q_points);
4870  return this->values_quad[0][q_point];
4871 }
4872 
4873 
4874 
4875 template <typename Number>
4876 inline
4879 ::get_gradient (const unsigned int q_point) const
4880 {
4881  // could use the base class gradient, but that involves too many inefficient
4882  // initialization operations on tensors
4883 
4884  Assert (this->gradients_quad_initialized==true,
4886  AssertIndexRange (q_point, this->data->n_q_points);
4887 
4889 
4890  // Cartesian cell
4891  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4892  {
4893  grad_out[0] = (this->gradients_quad[0][0][q_point] *
4894  this->cartesian_data[0][0]);
4895  }
4896  // cell with general/constant Jacobian
4897  else
4898  {
4899  const Tensor<2,1,VectorizedArray<Number> > &jac =
4900  this->cell_type == internal::MatrixFreeFunctions::general ?
4901  this->jacobian[q_point] : this->jacobian[0];
4902 
4903  grad_out[0] = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4904  }
4905  return grad_out;
4906 }
4907 
4908 
4909 
4910 template <typename Number>
4911 inline
4914 ::get_hessian (const unsigned int q_point) const
4915 {
4916  return BaseClass::get_hessian(q_point)[0];
4917 }
4918 
4919 
4920 
4921 template <typename Number>
4922 inline
4925 ::get_hessian_diagonal (const unsigned int q_point) const
4926 {
4927  return BaseClass::get_hessian_diagonal(q_point)[0];
4928 }
4929 
4930 
4931 
4932 template <typename Number>
4933 inline
4936 ::get_laplacian (const unsigned int q_point) const
4937 {
4938  return BaseClass::get_laplacian(q_point)[0];
4939 }
4940 
4941 
4942 
4943 template <typename Number>
4944 inline
4945 void
4948  const unsigned int dof)
4949 {
4950 #ifdef DEBUG
4951  this->dof_values_initialized = true;
4952  AssertIndexRange (dof, this->data->dofs_per_cell);
4953 #endif
4954  this->values_dofs[0][dof] = val_in;
4955 }
4956 
4957 
4958 
4959 template <typename Number>
4960 inline
4961 void
4964  const unsigned int q_point)
4965 {
4966 #ifdef DEBUG
4968  AssertIndexRange (q_point, this->data->n_q_points);
4969  this->values_quad_submitted = true;
4970 #endif
4971  if (this->cell_type == internal::MatrixFreeFunctions::general)
4972  {
4973  const VectorizedArray<Number> JxW = this->J_value[q_point];
4974  this->values_quad[0][q_point] = val_in * JxW;
4975  }
4976  else //if (this->cell_type < internal::MatrixFreeFunctions::general)
4977  {
4978  const VectorizedArray<Number> JxW = this->J_value[0] * this->quadrature_weights[q_point];
4979  this->values_quad[0][q_point] = val_in * JxW;
4980  }
4981 }
4982 
4983 
4984 
4985 template <typename Number>
4986 inline
4987 void
4989 ::submit_gradient (const Tensor<1,1,VectorizedArray<Number> > grad_in,
4990  const unsigned int q_point)
4991 {
4992 #ifdef DEBUG
4994  AssertIndexRange (q_point, this->data->n_q_points);
4995  this->gradients_quad_submitted = true;
4996 #endif
4997  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4998  {
4999  const VectorizedArray<Number> JxW = this->J_value[0] * this->quadrature_weights[q_point];
5000  this->gradients_quad[0][0][q_point] = (grad_in[0] *
5001  this->cartesian_data[0][0] *
5002  JxW);
5003  }
5004  // general/affine cell type
5005  else
5006  {
5007  const Tensor<2,1,VectorizedArray<Number> > &jac =
5008  this->cell_type == internal::MatrixFreeFunctions::general ?
5009  this->jacobian[q_point] : this->jacobian[0];
5010  const VectorizedArray<Number> JxW =
5011  this->cell_type == internal::MatrixFreeFunctions::general ?
5012  this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
5013 
5014  this->gradients_quad[0][0][q_point] = jac[0][0] * grad_in[0] * JxW;
5015  }
5016 }
5017 
5018 
5019 
5020 template <typename Number>
5021 inline
5024 ::integrate_value () const
5025 {
5026  return BaseClass::integrate_value()[0];
5027 }
5028 
5029 
5030 
5031 
5032 /*-------------------------- FEEvaluation -----------------------------------*/
5033 
5034 
5035 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5036  typename Number>
5037 inline
5039 ::FEEvaluation (const MatrixFree<dim,Number> &data_in,
5040  const unsigned int fe_no,
5041  const unsigned int quad_no)
5042  :
5043  BaseClass (data_in, fe_no, quad_no, fe_degree, static_n_q_points),
5044  dofs_per_cell (this->data->dofs_per_cell),
5045  n_q_points (this->data->n_q_points)
5046 {
5047  check_template_arguments(fe_no, 0);
5048 }
5049 
5050 
5051 
5052 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5053  typename Number>
5054 inline
5056 ::FEEvaluation (const Mapping<dim> &mapping,
5057  const FiniteElement<dim> &fe,
5058  const Quadrature<1> &quadrature,
5059  const UpdateFlags update_flags,
5060  const unsigned int first_selected_component)
5061  :
5062  BaseClass (mapping, fe, quadrature, update_flags,
5063  first_selected_component,
5064  static_cast<FEEvaluationBase<dim,1,Number>*>(0)),
5065  dofs_per_cell (this->data->dofs_per_cell),
5066  n_q_points (this->data->n_q_points)
5067 {
5068  check_template_arguments(numbers::invalid_unsigned_int, 0);
5069 }
5070 
5071 
5072 
5073 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5074  typename Number>
5075 inline
5078  const Quadrature<1> &quadrature,
5079  const UpdateFlags update_flags,
5080  const unsigned int first_selected_component)
5081  :
5082  BaseClass (StaticMappingQ1<dim>::mapping, fe, quadrature, update_flags,
5083  first_selected_component,
5084  static_cast<FEEvaluationBase<dim,1,Number>*>(0)),
5085  dofs_per_cell (this->data->dofs_per_cell),
5086  n_q_points (this->data->n_q_points)
5087 {
5088  check_template_arguments(numbers::invalid_unsigned_int, 0);
5089 }
5090 
5091 
5092 
5093 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5094  typename Number>
5095 template <int n_components_other>
5096 inline
5100  const unsigned int first_selected_component)
5101  :
5102  BaseClass (other.mapped_geometry->get_fe_values().get_mapping(),
5103  fe, other.mapped_geometry->get_quadrature(),
5104  other.mapped_geometry->get_fe_values().get_update_flags(),
5105  first_selected_component, &other),
5106  dofs_per_cell (this->data->dofs_per_cell),
5107  n_q_points (this->data->n_q_points)
5108 {
5109  check_template_arguments(numbers::invalid_unsigned_int, 0);
5110 }
5111 
5112 
5113 
5114 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5115  typename Number>
5116 inline
5118 ::FEEvaluation (const FEEvaluation &other)
5119  :
5120  BaseClass (other),
5121  dofs_per_cell (this->data->dofs_per_cell),
5122  n_q_points (this->data->n_q_points)
5123 {
5124  check_template_arguments(numbers::invalid_unsigned_int, 0);
5125 }
5126 
5127 
5128 
5129 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5130  typename Number>
5131 inline
5134 ::operator= (const FEEvaluation &other)
5135 {
5136  BaseClass::operator=(other);
5137  check_template_arguments(numbers::invalid_unsigned_int, 0);
5138  return *this;
5139 }
5140 
5141 
5142 
5143 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5144  typename Number>
5145 inline
5146 void
5148 ::check_template_arguments(const unsigned int fe_no,
5149  const unsigned int first_selected_component)
5150 {
5151  const unsigned int fe_degree_templ = fe_degree != -1 ? fe_degree : 0;
5152  const unsigned int n_q_points_1d_templ = n_q_points_1d > 0 ? n_q_points_1d : 1;
5153  if (fe_degree == -1)
5154  {
5155  if (this->data->element_type == internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
5156  {
5157  evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5158  dim, -1, 0, n_components_, Number>::evaluate;
5159  integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5160  dim, -1, 0, n_components_, Number>::integrate;
5161  }
5162  else if (this->data->element_type == internal::MatrixFreeFunctions::truncated_tensor)
5163  {
5164  evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5165  dim, -1, 0, n_components_, Number>::evaluate;
5166  integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5167  dim, -1, 0, n_components_, Number>::integrate;
5168  }
5169  else
5170  {
5171  evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5172  dim, -1, 0, n_components_, Number>::evaluate;
5173  integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5174  dim, -1, 0, n_components_, Number>::integrate;
5175  }
5176  }
5177  else
5178  switch (this->data->element_type)
5179  {
5180  case internal::MatrixFreeFunctions::tensor_symmetric:
5181  evaluate_funct =
5182  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
5183  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5184  Number>::evaluate;
5185  integrate_funct =
5186  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
5187  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5188  Number>::integrate;
5189  break;
5190 
5191  case internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0:
5192  evaluate_funct =
5193  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5194  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5195  Number>::evaluate;
5196  integrate_funct =
5197  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5198  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5199  Number>::integrate;
5200  break;
5201 
5202  case internal::MatrixFreeFunctions::tensor_general:
5203  evaluate_funct =
5204  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5205  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5206  Number>::evaluate;
5207  integrate_funct =
5208  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5209  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5210  Number>::integrate;
5211  break;
5212 
5213  case internal::MatrixFreeFunctions::tensor_gausslobatto:
5214  evaluate_funct =
5215  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
5216  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5217  Number>::evaluate;
5218  integrate_funct =
5219  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
5220  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5221  Number>::integrate;
5222  break;
5223 
5224  case internal::MatrixFreeFunctions::truncated_tensor:
5225  evaluate_funct =
5226  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5227  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5228  Number>::evaluate;
5229  integrate_funct =
5230  internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5231  dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5232  Number>::integrate;
5233  break;
5234 
5235  default:
5236  AssertThrow(false, ExcNotImplemented());
5237  }
5238 
5239  (void)fe_no;
5240  (void)first_selected_component;
5241 
5242 #ifdef DEBUG
5243  // print error message when the dimensions do not match. Propose a possible
5244  // fix
5245  if ((fe_degree != -1 && static_cast<unsigned int>(fe_degree) != this->data->fe_degree)
5246  ||
5247  (fe_degree != -1 && static_n_q_points != this->data->n_q_points))
5248  {
5249  std::string message =
5250  "-------------------------------------------------------\n";
5251  message += "Illegal arguments in constructor/wrong template arguments!\n";
5252  message += " Called --> FEEvaluation<dim,";
5253  message += Utilities::int_to_string(fe_degree) + ",";
5254  message += Utilities::int_to_string(n_q_points_1d);
5255  message += "," + Utilities::int_to_string(n_components);
5256  message += ",Number>(data";
5257  if (fe_no != numbers::invalid_unsigned_int)
5258  {
5259  message += ", " + Utilities::int_to_string(fe_no) + ", ";
5260  message += Utilities::int_to_string(this->quad_no);
5261  }
5262  message += ")\n";
5263 
5264  // check whether some other vector component has the correct number of
5265  // points
5266  unsigned int proposed_dof_comp = numbers::invalid_unsigned_int,
5267  proposed_quad_comp = numbers::invalid_unsigned_int;
5268  if (fe_no != numbers::invalid_unsigned_int)
5269  {
5270  if (static_cast<unsigned int>(fe_degree) == this->data->fe_degree)
5271  proposed_dof_comp = fe_no;
5272  else
5273  for (unsigned int no=0; no<this->matrix_info->n_components(); ++no)
5274  if (this->matrix_info->get_shape_info(no,0,this->active_fe_index,0).fe_degree
5275  == static_cast<unsigned int>(fe_degree))
5276  {
5277  proposed_dof_comp = no;
5278  break;
5279  }
5280  if (static_n_q_points ==
5281  this->mapping_info->mapping_data_gen[this->quad_no].n_q_points[this->active_quad_index])
5282  proposed_quad_comp = this->quad_no;
5283  else
5284  for (unsigned int no=0; no<this->mapping_info->mapping_data_gen.size(); ++no)
5285  if (this->mapping_info->mapping_data_gen[no].n_q_points[this->active_quad_index]
5286  == static_n_q_points)
5287  {
5288  proposed_quad_comp = no;
5289  break;
5290  }
5291  }
5292  if (proposed_dof_comp != numbers::invalid_unsigned_int &&
5293  proposed_quad_comp != numbers::invalid_unsigned_int)
5294  {
5295  if (proposed_dof_comp != fe_no)
5296  message += "Wrong vector component selection:\n";
5297  else
5298  message += "Wrong quadrature formula selection:\n";
5299  message += " Did you mean FEEvaluation<dim,";
5300  message += Utilities::int_to_string(fe_degree) + ",";
5301  message += Utilities::int_to_string(n_q_points_1d);
5302  message += "," + Utilities::int_to_string(n_components);
5303  message += ",Number>(data";
5304  if (fe_no != numbers::invalid_unsigned_int)
5305  {
5306  message += ", " + Utilities::int_to_string(proposed_dof_comp) + ", ";
5307  message += Utilities::int_to_string(proposed_quad_comp);
5308  }
5309  message += ")?\n";
5310  std::string correct_pos;
5311  if (proposed_dof_comp != fe_no)
5312  correct_pos = " ^ ";
5313  else
5314  correct_pos = " ";
5315  if (proposed_quad_comp != this->quad_no)
5316  correct_pos += " ^\n";
5317  else
5318  correct_pos += " \n";
5319  message += " " + correct_pos;
5320  }
5321  // ok, did not find the numbers specified by the template arguments in
5322  // the given list. Suggest correct template arguments
5323  const unsigned int proposed_n_q_points_1d = static_cast<unsigned int>(std::pow(1.001*this->data->n_q_points,1./dim));
5324  message += "Wrong template arguments:\n";
5325  message += " Did you mean FEEvaluation<dim,";
5326  message += Utilities::int_to_string(this->data->fe_degree) + ",";
5327  message += Utilities::int_to_string(proposed_n_q_points_1d);
5328  message += "," + Utilities::int_to_string(n_components);
5329  message += ",Number>(data";
5330  if (fe_no != numbers::invalid_unsigned_int)
5331  {
5332  message += ", " + Utilities::int_to_string(fe_no) + ", ";
5333  message += Utilities::int_to_string(this->quad_no);
5334  }
5335  message += ")?\n";
5336  std::string correct_pos;
5337  if (this->data->fe_degree != static_cast<unsigned int>(fe_degree))
5338  correct_pos = " ^";
5339  else
5340  correct_pos = " ";
5341  if (proposed_n_q_points_1d != n_q_points_1d)
5342  correct_pos += " ^\n";
5343  else
5344  correct_pos += " \n";
5345  message += " " + correct_pos;
5346 
5347  Assert (static_cast<unsigned int>(fe_degree) == this->data->fe_degree &&
5348  static_n_q_points == this->data->n_q_points,
5349  ExcMessage(message));
5350  }
5351  if (fe_no != numbers::invalid_unsigned_int)
5352  {
5353  AssertDimension (n_q_points,
5354  this->mapping_info->mapping_data_gen[this->quad_no].
5355  n_q_points[this->active_quad_index]);
5356  AssertDimension (this->data->dofs_per_cell * this->n_fe_components,
5357  this->dof_info->dofs_per_cell[this->active_fe_index]);
5358  }
5359 #endif
5360 }
5361 
5362 
5363 
5364 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5365  typename Number>
5366 inline
5369 ::quadrature_point (const unsigned int q) const
5370 {
5371  Assert (this->mapping_info->quadrature_points_initialized == true,
5372  ExcNotInitialized());
5373  Assert (this->quadrature_points != 0, ExcNotInitialized());
5374  AssertIndexRange (q, n_q_points);
5375 
5376  // Cartesian mesh: not all quadrature points are stored, only the
5377  // diagonal. Hence, need to find the tensor product index and retrieve the
5378  // value from that
5379  if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
5380  {
5382  switch (dim)
5383  {
5384  case 1:
5385  return this->quadrature_points[q];
5386  case 2:
5387  point[0] = this->quadrature_points[q%n_q_points_1d][0];
5388  point[1] = this->quadrature_points[q/n_q_points_1d][1];
5389  return point;
5390  case 3:
5391  point[0] = this->quadrature_points[q%n_q_points_1d][0];
5392  point[1] = this->quadrature_points[(q/n_q_points_1d)%n_q_points_1d][1];
5393  point[2] = this->quadrature_points[q/(n_q_points_1d*n_q_points_1d)][2];
5394  return point;
5395  default:
5396  Assert (false, ExcNotImplemented());
5397  return point;
5398  }
5399  }
5400  // all other cases: just return the respective data as it is fully stored
5401  else
5402  return this->quadrature_points[q];
5403 }
5404 
5405 
5406 
5407 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5408  typename Number>
5409 inline
5410 void
5412 ::evaluate (const bool evaluate_val,
5413  const bool evaluate_grad,
5414  const bool evaluate_lapl)
5415 {
5416  Assert (this->dof_values_initialized == true,
5418  Assert(this->matrix_info != 0 ||
5419  this->mapped_geometry->is_initialized(), ExcNotInitialized());
5420 
5421  // Select algorithm matching the element type at run time (the function
5422  // pointer is easy to predict, so negligible in cost)
5423  evaluate_funct (*this->data, &this->values_dofs[0], this->values_quad,
5424  this->gradients_quad, this->hessians_quad, this->scratch_data,
5425  evaluate_val, evaluate_grad, evaluate_lapl);
5426 
5427 #ifdef DEBUG
5428  if (evaluate_val == true)
5429  this->values_quad_initialized = true;
5430  if (evaluate_grad == true)
5431  this->gradients_quad_initialized = true;
5432  if (evaluate_lapl == true)
5433  this->hessians_quad_initialized = true;
5434 #endif
5435 }
5436 
5437 
5438 
5439 template <int dim, int fe_degree, int n_q_points_1d, int n_components_,
5440  typename Number>
5441 inline
5442 void
5444 ::integrate (bool integrate_val,bool integrate_grad)
5445 {
5446  if (integrate_val == true)
5447  Assert (this->values_quad_submitted == true,
5449  if (integrate_grad == true)
5450  Assert (this->gradients_quad_submitted == true,
5452  Assert(this->matrix_info != 0 ||
5453  this->mapped_geometry->is_initialized(), ExcNotInitialized());
5454 
5455  // Select algorithm matching the element type at run time (the function
5456  // pointer is easy to predict, so negligible in cost)
5457  integrate_funct (*this->data, this->values_dofs, this->values_quad,
5458  this->gradients_quad, this->scratch_data,
5459  integrate_val, integrate_grad);
5460 
5461 #ifdef DEBUG
5462  this->dof_values_initialized = true;
5463 #endif
5464 }
5465 
5466 
5467 
5468 #endif // ifndef DOXYGEN
5469 
5470 
5471 DEAL_II_NAMESPACE_CLOSE
5472 
5473 #endif
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info(const unsigned int fe_component=0, const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0, const unsigned int hp_active_quad_index=0) const
void set_data_pointers()
static const unsigned int invalid_unsigned_int
Definition: types.h:170
const Number * constraint_pool_begin(const unsigned int pool_index) const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
AlignedVector< std::pair< Tensor< 1, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > cartesian_data
Definition: mapping_info.h:146
unsigned int get_cell_data_number() const
Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArray< Number > > > get_hessian(const unsigned int q_point) const
void check_template_arguments(const unsigned int fe_no, const unsigned int first_selected_component)
const MatrixFree< dim, Number > * matrix_info
FEEvaluation & operator=(const FEEvaluation &other)
void reinit(const unsigned int cell)
const Tensor< 1,(dim >1?dim *(dim-1)/2:1), Tensor< 1, dim, VectorizedArray< Number > > > * jacobian_grad_upper
AlignedVector< std::pair< Tensor< 2, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > affine_data
Definition: mapping_info.h:161
unsigned int cell_data_number
void read_write_operation(const VectorOperation &operation, VectorType *vectors[]) const
static ::ExceptionBase & ExcAccessToUninitializedField()
const unsigned int first_selected_component
value_type get_laplacian(const unsigned int q_point) const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
const unsigned int dofs_per_cell
unsigned int n_components() const
VectorizedArray< Number > JxW(const unsigned int q_point) const
void evaluate(const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians=false)
const unsigned int n_fe_components
void read_dof_values(const VectorType &src, const unsigned int first_index=0)
void integrate(const bool integrate_values, const bool integrate_gradients)
void read_dof_values_plain(const VectorType &src, const unsigned int first_index=0)
bool gradients_quad_initialized
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
const VectorizedArray< Number > * begin_gradients() const
AlignedVector< VectorizedArray< Number > > * scratch_data_array
std::vector< types::global_dof_index > old_style_dof_indices
unsigned int get_cell_data_index(const unsigned int cell_chunk_no) const
Definition: mapping_info.h:369
const VectorizedArray< Number > * begin_hessians() const
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian
friend class FEEvaluationBase
STL namespace.
static ::ExceptionBase & ExcNotInitialized()
std::vector< MappingInfoDependent > mapping_data_gen
Definition: mapping_info.h:290
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
FEEvaluationAccess(const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no, const unsigned int quad_no, const unsigned int fe_degree, const unsigned int n_q_points)
const internal::MatrixFreeFunctions::DoFInfo * dof_info
void submit_gradient(const gradient_type grad_in, const unsigned int q_point)
SymmetricTensor< 2, dim, VectorizedArray< Number > > get_symmetric_gradient(const unsigned int q_point) const
void submit_symmetric_gradient(const SymmetricTensor< 2, dim, VectorizedArray< Number > > grad_in, const unsigned int q_point)
const VectorizedArray< Number > * begin_dof_values() const
Definition: point.h:89
VectorizedArray< Number > * scratch_data
const Tensor< 1, dim, VectorizedArray< Number > > * cartesian_data
void submit_dof_value(const value_type val_in, const unsigned int dof)
const internal::MatrixFreeFunctions::ShapeInfo< Number > * data
VectorizedArray< Number > * values_quad[n_components]
AlignedVector< VectorizedArray< Number > > * acquire_scratch_data() const
FEEvaluationBase & operator=(const FEEvaluationBase &other)
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int n_q_points
void fill_JxW_values(AlignedVector< VectorizedArray< Number > > &JxW_values) const
value_type get_dof_value(const unsigned int dof) const
void submit_divergence(const VectorizedArray< Number > div_in, const unsigned int q_point)
const unsigned int active_fe_index
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
Definition: fe.h:2959
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int global_dof_index
Definition: types.h:88
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian_grad
#define Assert(cond, exc)
Definition: exceptions.h:313
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:2863
UpdateFlags
FEEvaluation(const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no=0, const unsigned int quad_no=0)
const unsigned int active_quad_index
#define DeclException0(Exception0)
Definition: exceptions.h:541
const Number * constraint_pool_end(const unsigned int pool_index) const
const internal::MatrixFreeFunctions::MappingInfo< dim, Number > * mapping_info
std::vector< types::global_dof_index > local_dof_indices
const unsigned int quad_no
const Point< dim, VectorizedArray< Number > > * quadrature_points
void submit_curl(const Tensor< 1, dim==2?1:dim, VectorizedArray< Number > > curl_in, const unsigned int q_point)
bool mapping_initialized() const
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:85
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0) const
const VectorizedArray< Number > * quadrature_weights
gradient_type get_hessian_diagonal(const unsigned int q_point) const
bool hessians_quad_initialized
Tensor< 1,(dim==2?1:dim), VectorizedArray< Number > > get_curl(const unsigned int q_point) const
const internal::MatrixFreeFunctions::SizeInfo & get_size_info() const
FEEvaluationAccess & operator=(const FEEvaluationAccess &other)
void set_dof_values(VectorType &dst, const unsigned int first_index=0) const
void release_scratch_data(const AlignedVector< VectorizedArray< Number > > *memory) const
unsigned int cell
value_type integrate_value() const
void submit_value(const value_type val_in, const unsigned int q_point)
internal::MatrixFreeFunctions::CellType cell_type
void(* evaluate_funct)(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], VectorizedArray< Number > *hessians_quad[][(dim *(dim+1))/2], VectorizedArray< Number > *scratch_data, const bool evaluate_val, const bool evaluate_grad, const bool evaluate_lapl)
VectorizedArray< Number > get_divergence(const unsigned int q_point) const
gradient_type get_gradient(const unsigned int q_point) const
VectorizedArray< Number > * values_dofs[n_components]
Definition: mpi.h:41
ArrayView< VectorizedArray< Number > > get_scratch_data() const
iterator end()
const VectorizedArray< Number > * begin_values() const
VectorizedArray< Number > * hessians_quad[n_components][(dim *(dim+1))/2]
bool partitioners_are_compatible(const Utilities::MPI::Partitioner &part) const
Number local_element(const size_type local_index) const
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info() const
static ::ExceptionBase & ExcNotImplemented()
bool indices_initialized() const
std::vector< unsigned int > lexicographic_numbering
Definition: shape_info.h:190
bool values_quad_initialized
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
const std::vector< unsigned int > & get_internal_dof_numbering() const
value_type get_value(const unsigned int q_point) const
const VectorizedArray< Number > * J_value
CellType get_cell_type(const unsigned int cell_chunk_no) const
Definition: mapping_info.h:356
Point< 3 > point(const gp_Pnt &p)
Definition: utilities.cc:176
VectorizedArray< Number > * gradients_quad[n_components][dim]
internal::MatrixFreeFunctions::CellType get_cell_type() const
Point< dim, VectorizedArray< Number > > quadrature_point(const unsigned int q_point) const
void(* integrate_funct)(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], VectorizedArray< Number > *scratch_data, const bool evaluate_val, const bool evaluate_grad)
bool gradients_quad_submitted
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > make_vectorized_array(const Number &u)
std_cxx1x::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
static ::ExceptionBase & ExcInternalError()