17 #ifndef dealii__matrix_free_fe_evaluation_h 18 #define dealii__matrix_free_fe_evaluation_h 21 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/template_constraints.h> 24 #include <deal.II/base/symmetric_tensor.h> 25 #include <deal.II/base/vectorization.h> 26 #include <deal.II/base/smartpointer.h> 27 #include <deal.II/matrix_free/matrix_free.h> 28 #include <deal.II/matrix_free/shape_info.h> 29 #include <deal.II/matrix_free/mapping_data_on_the_fly.h> 32 DEAL_II_NAMESPACE_OPEN
41 template <
typename>
class Vector;
49 template <
int dim,
int fe_degree,
int n_q_points_1d = fe_degree+1,
50 int n_components_ = 1,
typename Number =
double >
class FEEvaluation;
76 template <
int dim,
int n_components_,
typename Number>
80 typedef Number number_type;
83 static const unsigned int dimension = dim;
84 static const unsigned int n_components = n_components_;
98 void reinit (
const unsigned int cell);
112 template <
typename DoFHandlerType,
bool level_dof_access>
135 unsigned int get_cell_data_number()
const;
143 internal::MatrixFreeFunctions::CellType get_cell_type()
const;
149 get_shape_info()
const;
188 template <
typename VectorType>
189 void read_dof_values (
const VectorType &src);
199 template <
typename VectorType>
200 void read_dof_values (
const std::vector<VectorType> &src,
201 const unsigned int first_index=0);
207 template <
typename VectorType>
208 void read_dof_values (
const std::vector<VectorType *> &src,
209 const unsigned int first_index=0);
232 template <
typename VectorType>
233 void read_dof_values_plain (
const VectorType &src);
242 template <
typename VectorType>
243 void read_dof_values_plain (
const std::vector<VectorType> &src,
244 const unsigned int first_index=0);
251 template <
typename VectorType>
252 void read_dof_values_plain (
const std::vector<VectorType *> &src,
253 const unsigned int first_index=0);
271 template<
typename VectorType>
272 void distribute_local_to_global (VectorType &dst)
const;
283 template<
typename VectorType>
284 void distribute_local_to_global (std::vector<VectorType> &dst,
285 const unsigned int first_index=0)
const;
291 template<
typename VectorType>
292 void distribute_local_to_global (std::vector<VectorType *> &dst,
293 const unsigned int first_index=0)
const;
311 template<
typename VectorType>
312 void set_dof_values (VectorType &dst)
const;
323 template<
typename VectorType>
324 void set_dof_values (std::vector<VectorType> &dst,
325 const unsigned int first_index=0)
const;
331 template<
typename VectorType>
332 void set_dof_values (std::vector<VectorType *> &dst,
333 const unsigned int first_index=0)
const;
354 value_type get_dof_value (
const unsigned int dof)
const;
366 void submit_dof_value (
const value_type val_in,
367 const unsigned int dof);
380 value_type get_value (
const unsigned int q_point)
const;
393 void submit_value (
const value_type val_in,
394 const unsigned int q_point);
405 gradient_type get_gradient (
const unsigned int q_point)
const;
419 void submit_gradient(
const gradient_type grad_in,
420 const unsigned int q_point);
433 get_hessian (
const unsigned int q_point)
const;
443 gradient_type get_hessian_diagonal (
const unsigned int q_point)
const;
455 value_type get_laplacian (
const unsigned int q_point)
const;
465 value_type integrate_value ()
const;
576 const std::vector<unsigned int> &
577 get_internal_dof_numbering()
const;
591 const unsigned int fe_no,
592 const unsigned int quad_no,
593 const unsigned int fe_degree,
594 const unsigned int n_q_points);
630 template <
int n_components_other>
635 const unsigned int first_selected_component,
644 FEEvaluationBase (
const FEEvaluationBase &other);
652 template<
typename VectorType,
typename VectorOperation>
654 VectorType *vectors[])
const;
663 template<
typename VectorType>
664 void read_dof_values_plain (
const VectorType *src_data[]);
737 const internal::MatrixFreeFunctions::DoFInfo *
dof_info;
870 std_cxx1x::shared_ptr<internal::MatrixFreeFunctions::MappingDataOnTheFly<dim,Number> >
mapped_geometry;
894 template <
int,
int,
int,
int,
typename>
friend class FEEvaluation;
906 template <
int dim,
int n_components_,
typename Number>
910 typedef Number number_type;
913 static const unsigned int dimension = dim;
914 static const unsigned int n_components = n_components_;
926 const unsigned int fe_no,
927 const unsigned int quad_no,
928 const unsigned int fe_degree,
929 const unsigned int n_q_points);
935 template <
int n_components_other>
940 const unsigned int first_selected_component,
946 FEEvaluationAccess (
const FEEvaluationAccess &other);
960 template <
int dim,
typename Number>
964 typedef Number number_type;
967 static const unsigned int dimension = dim;
979 value_type get_dof_value (
const unsigned int dof)
const;
985 void submit_dof_value (
const value_type val_in,
986 const unsigned int dof);
995 value_type get_value (
const unsigned int q_point)
const;
1004 void submit_value (
const value_type val_in,
1005 const unsigned int q_point);
1012 gradient_type get_gradient (
const unsigned int q_point)
const;
1022 void submit_gradient(
const gradient_type grad_in,
1023 const unsigned int q_point);
1032 get_hessian (
unsigned int q_point)
const;
1038 gradient_type get_hessian_diagonal (
const unsigned int q_point)
const;
1044 value_type get_laplacian (
const unsigned int q_point)
const;
1054 value_type integrate_value ()
const;
1065 const unsigned int fe_no,
1066 const unsigned int quad_no,
1067 const unsigned int fe_degree,
1068 const unsigned int n_q_points);
1074 template <
int n_components_other>
1079 const unsigned int first_selected_component,
1085 FEEvaluationAccess (
const FEEvaluationAccess &other);
1099 template <
int dim,
typename Number>
1103 typedef Number number_type;
1106 static const unsigned int dimension = dim;
1107 static const unsigned int n_components = dim;
1114 gradient_type get_gradient (
const unsigned int q_point)
const;
1129 get_symmetric_gradient (
const unsigned int q_point)
const;
1136 get_curl (
const unsigned int q_point)
const;
1145 get_hessian (
const unsigned int q_point)
const;
1151 gradient_type get_hessian_diagonal (
const unsigned int q_point)
const;
1161 void submit_gradient(
const gradient_type grad_in,
1162 const unsigned int q_point);
1173 const unsigned int q_point);
1184 const unsigned int q_point);
1195 const unsigned int q_point);
1202 const unsigned int q_point);
1213 const unsigned int fe_no,
1214 const unsigned int quad_no,
1215 const unsigned int dofs_per_cell,
1216 const unsigned int n_q_points);
1222 template <
int n_components_other>
1227 const unsigned int first_selected_component,
1233 FEEvaluationAccess (
const FEEvaluationAccess &other);
1246 template <
typename Number>
1250 typedef Number number_type;
1253 static const unsigned int dimension = 1;
1265 value_type get_dof_value (
const unsigned int dof)
const;
1271 void submit_dof_value (
const value_type val_in,
1272 const unsigned int dof);
1281 value_type get_value (
const unsigned int q_point)
const;
1290 void submit_value (
const value_type val_in,
1291 const unsigned int q_point);
1298 gradient_type get_gradient (
const unsigned int q_point)
const;
1308 void submit_gradient(
const gradient_type grad_in,
1309 const unsigned int q_point);
1318 get_hessian (
unsigned int q_point)
const;
1324 gradient_type get_hessian_diagonal (
const unsigned int q_point)
const;
1330 value_type get_laplacian (
const unsigned int q_point)
const;
1340 value_type integrate_value ()
const;
1351 const unsigned int fe_no,
1352 const unsigned int quad_no,
1353 const unsigned int fe_degree,
1354 const unsigned int n_q_points);
1360 template <
int n_components_other>
1365 const unsigned int first_selected_component,
1371 FEEvaluationAccess (
const FEEvaluationAccess &other);
1528 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
1534 typedef Number number_type;
1537 static const unsigned int dimension = dim;
1538 static const unsigned int n_components = n_components_;
1549 const unsigned int fe_no = 0,
1550 const unsigned int quad_no = 0);
1582 const unsigned int first_selected_component = 0);
1592 const unsigned int first_selected_component = 0);
1604 template <
int n_components_other>
1607 const unsigned int first_selected_component = 0);
1615 FEEvaluation (
const FEEvaluation &other);
1625 void evaluate (
const bool evaluate_val,
1626 const bool evaluate_grad,
1627 const bool evaluate_hess =
false);
1636 void integrate (
const bool integrate_val,
1637 const bool integrate_grad);
1643 quadrature_point (
const unsigned int q_point)
const;
1662 void check_template_arguments(
const unsigned int fe_no);
1667 void set_data_pointers();
1677 const bool evaluate_val,
1678 const bool evaluate_grad,
1679 const bool evaluate_lapl);
1688 const bool evaluate_val,
1689 const bool evaluate_grad);
1696 namespace MatrixFreeFunctions
1700 template <
int dim,
int degree>
1701 struct DGP_dofs_per_cell
1704 static const unsigned int value =
1705 (DGP_dofs_per_cell<dim-1,degree>::value * (degree+dim)) / dim;
1709 template <
int degree>
1710 struct DGP_dofs_per_cell<1,degree>
1712 static const unsigned int value = degree+1;
1726 template <
int dim,
int n_components_,
typename Number>
1730 const unsigned int fe_no_in,
1731 const unsigned int quad_no_in,
1732 const unsigned int fe_degree,
1733 const unsigned int n_q_points)
1735 quad_no (quad_no_in),
1736 n_fe_components (data_in.
get_dof_info(fe_no_in).n_components),
1737 active_fe_index (data_in.
get_dof_info(fe_no_in).fe_index_from_degree
1739 active_quad_index (data_in.get_mapping_info().
1740 mapping_data_gen[quad_no_in].
1741 quad_index_from_n_q_points(n_q_points)),
1742 matrix_info (&data_in),
1744 mapping_info (&data_in.get_mapping_info()),
1746 (fe_no_in, quad_no_in, active_fe_index,
1747 active_quad_index)),
1752 quadrature_weights[active_quad_index].begin()),
1753 quadrature_points (0),
1755 jacobian_grad_upper(0),
1757 cell_type (internal::MatrixFreeFunctions::undefined),
1758 cell_data_number (0),
1759 first_selected_component (0)
1761 for (
unsigned int c=0; c<n_components_; ++c)
1765 for (
unsigned int d=0; d<dim; ++d)
1766 gradients_quad[c][d] = 0;
1767 for (
unsigned int d=0; d<(dim*dim+dim)/2; ++d)
1768 hessians_quad[c][d] = 0;
1771 ExcNotInitialized());
1775 dof_info->dofs_per_cell[active_fe_index]/n_fe_components);
1778 Assert (n_fe_components == 1 ||
1779 n_components == 1 ||
1780 n_components == n_fe_components,
1781 ExcMessage (
"The underlying FE is vector-valued. In this case, the " 1782 "template argument n_components must be a the same " 1783 "as the number of underlying vector components."));
1792 template <
int dim,
int n_components_,
typename Number>
1793 template <
int n_components_other>
1800 const unsigned int first_selected_component,
1804 n_fe_components (n_components_),
1805 active_fe_index (-1),
1806 active_quad_index (-1),
1812 data (stored_shape_info.get()),
1816 quadrature_weights (0),
1817 quadrature_points (0),
1819 jacobian_grad_upper(0),
1821 cell_type (internal::MatrixFreeFunctions::general),
1822 cell_data_number (0),
1827 const unsigned int base_element_number =
1829 for (
unsigned int c=0; c<n_components_; ++c)
1833 for (
unsigned int d=0; d<dim; ++d)
1834 gradients_quad[c][d] = 0;
1835 for (
unsigned int d=0; d<(dim*dim+dim)/2; ++d)
1836 hessians_quad[c][d] = 0;
1845 MappingDataOnTheFly<dim,Number>(mapping, quadrature,
1847 jacobian = mapped_geometry->get_inverse_jacobians().begin();
1848 J_value = mapped_geometry->get_JxW_values().begin();
1849 quadrature_points = mapped_geometry->get_quadrature_points().begin();
1853 ExcMessage(
"The underlying element must at least contain as many " 1854 "components as requested by this class"));
1859 template <
int dim,
int n_components_,
typename Number>
1885 for (
unsigned int c=0; c<n_components_; ++c)
1889 for (
unsigned int d=0; d<dim; ++d)
1890 gradients_quad[c][d] = 0;
1891 for (
unsigned int d=0; d<(dim*dim+dim)/2; ++d)
1892 hessians_quad[c][d] = 0;
1898 mapped_geometry.reset
1900 MappingDataOnTheFly<dim,Number>(other.
mapped_geometry->get_fe_values().get_mapping(),
1903 jacobian = mapped_geometry->get_inverse_jacobians().begin();
1904 J_value = mapped_geometry->get_JxW_values().begin();
1905 quadrature_points = mapped_geometry->get_quadrature_points().begin();
1911 template <
int dim,
int n_components_,
typename Number>
1916 Assert (mapped_geometry == 0,
1917 ExcMessage(
"FEEvaluation was initialized without a matrix-free object." 1918 " Integer indexing is not possible"));
1919 if (mapped_geometry != 0)
1921 Assert (dof_info != 0, ExcNotInitialized());
1922 Assert (mapping_info != 0, ExcNotInitialized());
1925 dof_info->cell_active_fe_index[cell_in] : 0),
1934 mapping_data_gen[quad_no].rowstart_q_points.size());
1936 rowstart_q_points[cell];
1938 quadrature_points.size());
1943 if (cell_type == internal::MatrixFreeFunctions::cartesian)
1945 cartesian_data = &mapping_info->
cartesian_data[cell_data_number].first;
1946 J_value = &mapping_info->
cartesian_data[cell_data_number].second;
1948 else if (cell_type == internal::MatrixFreeFunctions::affine)
1950 jacobian = &mapping_info->
affine_data[cell_data_number].first;
1951 J_value = &mapping_info->
affine_data[cell_data_number].second;
1955 const unsigned int rowstart = mapping_info->
1956 mapping_data_gen[quad_no].rowstart_jacobians[cell_data_number];
1958 mapping_data_gen[quad_no].jacobians.size());
1964 mapping_data_gen[quad_no].JxW_values.size());
1966 JxW_values[rowstart]);
1971 mapping_data_gen[quad_no].jacobians_grad_diag.size());
1973 jacobians_grad_diag[rowstart];
1975 mapping_data_gen[quad_no].jacobians_grad_upper.size());
1977 jacobians_grad_upper[rowstart];
1982 dof_values_initialized =
false;
1983 values_quad_initialized =
false;
1984 gradients_quad_initialized =
false;
1985 hessians_quad_initialized =
false;
1991 template <
int dim,
int n_components_,
typename Number>
1992 template <
typename DoFHandlerType,
bool level_dof_access>
1999 ExcMessage(
"Cannot use initialization from cell iterator if " 2000 "initialized from MatrixFree object. Use variant for " 2001 "on the fly computation with arguments as for FEValues " 2003 Assert(mapped_geometry.get() != 0, ExcNotInitialized());
2005 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2006 if (level_dof_access)
2007 cell->get_mg_dof_indices(local_dof_indices);
2009 cell->get_dof_indices(local_dof_indices);
2014 template <
int dim,
int n_components_,
typename Number>
2021 ExcMessage(
"Cannot use initialization from cell iterator if " 2022 "initialized from MatrixFree object. Use variant for " 2023 "on the fly computation with arguments as for FEValues " 2025 Assert(mapped_geometry.get() != 0, ExcNotInitialized());
2026 mapped_geometry->reinit(cell);
2031 template <
int dim,
int n_components_,
typename Number>
2038 return cell_data_number;
2043 template <
int dim,
int n_components_,
typename Number>
2045 internal::MatrixFreeFunctions::CellType
2054 template <
int dim,
int n_components_,
typename Number>
2059 Assert(data != 0, ExcInternalError());
2065 template <
int dim,
int n_components_,
typename Number>
2072 Assert (this->J_value != 0, ExcNotImplemented());
2073 if (this->cell_type == internal::MatrixFreeFunctions::cartesian ||
2074 this->cell_type == internal::MatrixFreeFunctions::affine)
2076 Assert (this->mapping_info != 0, ExcNotImplemented());
2078 for (
unsigned int q=0; q<this->data->
n_q_points; ++q)
2079 JxW_values[q] = J * this->quadrature_weights[q];
2082 for (
unsigned int q=0; q<data->
n_q_points; ++q)
2083 JxW_values[q] = this->J_value[q];
2091 template <
typename VectorType>
2093 typename VectorType::value_type &
2094 vector_access (VectorType &vec,
2095 const unsigned int entry)
2103 template <
typename VectorType>
2105 typename VectorType::value_type
2106 vector_access (
const VectorType &vec,
2107 const unsigned int entry)
2117 template <
typename Number>
2121 const unsigned int entry)
2131 template <
typename Number>
2135 const unsigned int entry)
2144 template <
typename VectorType>
2146 void check_vector_compatibility (
const VectorType &vec,
2147 const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2150 dof_info.vector_partitioner->size());
2153 template <
typename Number>
2156 const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2159 ExcMessage(
"The parallel layout of the given vector is not " 2160 "compatible with the parallel partitioning in MatrixFree. " 2161 "Use MatrixFree::initialize_dof_vector to get a " 2162 "compatible vector."));
2166 template <
typename Number>
2169 template <
typename VectorType>
2170 void process_dof (
const unsigned int index,
2174 res = vector_access (const_cast<const VectorType &>(vec), index);
2177 template <
typename VectorType>
2182 res =
const_cast<const VectorType &
>(vec)(index);
2185 void pre_constraints (
const Number &,
2191 template <
typename VectorType>
2192 void process_constraint (
const unsigned int index,
2193 const Number weight,
2197 res += weight * vector_access (const_cast<const VectorType &>(vec), index);
2200 void post_constraints (
const Number &sum,
2201 Number &write_pos)
const 2206 void process_empty (Number &res)
const 2213 template <
typename Number>
2214 struct VectorDistributorLocalToGlobal
2216 template <
typename VectorType>
2217 void process_dof (
const unsigned int index,
2221 vector_access (vec, index) += res;
2224 template <
typename VectorType>
2232 void pre_constraints (
const Number &input,
2238 template <
typename VectorType>
2239 void process_constraint (
const unsigned int index,
2240 const Number weight,
2244 vector_access (vec, index) += weight * res;
2247 void post_constraints (
const Number &,
2252 void process_empty (Number &)
const 2259 template <
typename Number>
2262 template <
typename VectorType>
2263 void process_dof (
const unsigned int index,
2267 vector_access (vec, index) = res;
2270 template <
typename VectorType>
2278 void pre_constraints (
const Number &,
2283 template <
typename VectorType>
2284 void process_constraint (
const unsigned int,
2291 void post_constraints (
const Number &,
2296 void process_empty (Number &)
const 2304 template <
typename VectorType,
bool>
2305 struct BlockVectorSelector {};
2307 template <
typename VectorType>
2308 struct BlockVectorSelector<VectorType,true>
2310 typedef typename VectorType::BlockType BaseVectorType;
2312 static BaseVectorType *get_vector_component (VectorType &vec,
2313 const unsigned int component)
2316 return &vec.block(component);
2320 template <
typename VectorType>
2321 struct BlockVectorSelector<VectorType,false>
2323 typedef VectorType BaseVectorType;
2325 static BaseVectorType *get_vector_component (VectorType &vec,
2335 template <
int dim,
int n_components_,
typename Number>
2336 template<
typename VectorType,
typename VectorOperation>
2341 VectorType *src[])
const 2352 if (matrix_info == 0)
2354 Assert (!local_dof_indices.empty(), ExcNotInitialized());
2356 unsigned int index = first_selected_component * this->data->
dofs_per_cell;
2357 for (
unsigned int comp = 0; comp<n_components; ++comp)
2359 for (
unsigned int i=0; i<this->data->
dofs_per_cell; ++i, ++index)
2362 *src[0], values_dofs[comp][i][0]);
2363 for (
unsigned int v=1; v<VectorizedArray<Number>::n_array_elements; ++v)
2364 operation.process_empty(values_dofs[comp][i][v]);
2370 Assert (dof_info != 0, ExcNotInitialized());
2372 ExcNotInitialized());
2378 const unsigned int *dof_indices = dof_info->begin_indices(cell);
2379 const std::pair<unsigned short,unsigned short> *indicators =
2380 dof_info->begin_indicators(cell);
2381 const std::pair<unsigned short,unsigned short> *indicators_end =
2382 dof_info->end_indicators(cell);
2383 unsigned int ind_local = 0;
2384 const unsigned int dofs_per_cell = this->data->
dofs_per_cell;
2386 const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
2387 const bool at_irregular_cell = n_irreg_components_filled > 0;
2391 if (n_fe_components == 1)
2393 const unsigned int n_local_dofs =
2395 for (
unsigned int comp=0; comp<n_components; ++comp)
2396 internal::check_vector_compatibility (*src[comp], *dof_info);
2397 Number *local_data [n_components];
2398 for (
unsigned int comp=0; comp<n_components; ++comp)
2400 const_cast<Number *>(&values_dofs[comp][0][0]);
2404 if (at_irregular_cell ==
false)
2407 if (indicators != indicators_end)
2409 for ( ; indicators != indicators_end; ++indicators)
2412 for (
unsigned int j=0; j<indicators->first; ++j)
2413 for (
unsigned int comp=0; comp<n_components; ++comp)
2414 operation.process_dof (dof_indices[j], *src[comp],
2415 local_data[comp][ind_local+j]);
2417 ind_local += indicators->first;
2418 dof_indices += indicators->first;
2422 Number value [n_components];
2423 for (
unsigned int comp=0; comp<n_components; ++comp)
2424 operation.pre_constraints (local_data[comp][ind_local],
2427 const Number *data_val =
2429 const Number *end_pool =
2431 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
2432 for (
unsigned int comp=0; comp<n_components; ++comp)
2433 operation.process_constraint (*dof_indices, *data_val,
2434 *src[comp], value[comp]);
2436 for (
unsigned int comp=0; comp<n_components; ++comp)
2437 operation.post_constraints (value[comp],
2438 local_data[comp][ind_local]);
2444 for (; ind_local < n_local_dofs; ++dof_indices, ++ind_local)
2446 for (
unsigned int comp=0; comp<n_components; ++comp)
2447 operation.process_dof (*dof_indices, *src[comp],
2448 local_data[comp][ind_local]);
2456 static_cast<int>(n_local_dofs));
2457 for (
unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
2458 for (
unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
2459 for (
unsigned int comp=0; comp<n_components; ++comp)
2460 operation.process_dof (dof_indices[j+v], *src[comp],
2461 local_data[comp][j+v]);
2470 Assert (n_irreg_components_filled > 0, ExcInternalError());
2471 for ( ; indicators != indicators_end; ++indicators)
2473 for (
unsigned int j=0; j<indicators->first; ++j)
2477 for (
unsigned int comp=0; comp<n_components; ++comp)
2478 operation.process_dof (dof_indices[j], *src[comp],
2479 local_data[comp][ind_local]);
2484 >= n_irreg_components_filled)
2486 for (
unsigned int comp=0; comp<n_components; ++comp)
2487 operation.process_empty (local_data[comp][ind_local]);
2491 dof_indices += indicators->first;
2495 Number value [n_components];
2496 for (
unsigned int comp=0; comp<n_components; ++comp)
2497 operation.pre_constraints (local_data[comp][ind_local],
2500 const Number *data_val =
2502 const Number *end_pool =
2505 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
2506 for (
unsigned int comp=0; comp<n_components; ++comp)
2507 operation.process_constraint (*dof_indices, *data_val,
2508 *src[comp], value[comp]);
2510 for (
unsigned int comp=0; comp<n_components; ++comp)
2511 operation.post_constraints (value[comp],
2512 local_data[comp][ind_local]);
2515 >= n_irreg_components_filled)
2517 for (
unsigned int comp=0; comp<n_components; ++comp)
2518 operation.process_empty (local_data[comp][ind_local]);
2522 for (; ind_local<n_local_dofs; ++dof_indices)
2524 Assert (dof_indices != dof_info->end_indices(cell),
2525 ExcInternalError());
2529 for (
unsigned int comp=0; comp<n_components; ++comp)
2530 operation.process_dof (*dof_indices, *src[comp],
2531 local_data[comp][ind_local]);
2534 >= n_irreg_components_filled)
2536 for (
unsigned int comp=0; comp<n_components; ++comp)
2537 operation.process_empty(local_data[comp][ind_local]);
2549 internal::check_vector_compatibility (*src[0], *dof_info);
2550 Assert (n_fe_components == n_components_, ExcNotImplemented());
2551 const unsigned int n_local_dofs =
2553 Number *local_data =
2554 const_cast<Number *
>(&values_dofs[0][0][0]);
2555 if (at_irregular_cell ==
false)
2558 if (indicators != indicators_end)
2560 for ( ; indicators != indicators_end; ++indicators)
2563 for (
unsigned int j=0; j<indicators->first; ++j)
2564 operation.process_dof (dof_indices[j], *src[0],
2565 local_data[ind_local+j]);
2566 ind_local += indicators->first;
2567 dof_indices += indicators->first;
2572 operation.pre_constraints (local_data[ind_local], value);
2574 const Number *data_val =
2576 const Number *end_pool =
2579 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
2580 operation.process_constraint (*dof_indices, *data_val,
2583 operation.post_constraints (value, local_data[ind_local]);
2588 for (; ind_local<n_local_dofs; ++dof_indices, ++ind_local)
2589 operation.process_dof (*dof_indices, *src[0],
2590 local_data[ind_local]);
2591 Assert (dof_indices == dof_info->end_indices(cell),
2592 ExcInternalError());
2599 static_cast<int>(n_local_dofs));
2600 for (
unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
2601 for (
unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
2602 operation.process_dof (dof_indices[j+v], *src[0],
2612 Assert (n_irreg_components_filled > 0, ExcInternalError());
2613 for ( ; indicators != indicators_end; ++indicators)
2615 for (
unsigned int j=0; j<indicators->first; ++j)
2619 operation.process_dof (dof_indices[j], *src[0],
2620 local_data[ind_local]);
2625 >= n_irreg_components_filled)
2627 operation.process_empty (local_data[ind_local]);
2631 dof_indices += indicators->first;
2636 operation.pre_constraints (local_data[ind_local], value);
2638 const Number *data_val =
2640 const Number *end_pool =
2643 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
2644 operation.process_constraint (*dof_indices, *data_val,
2647 operation.post_constraints (value, local_data[ind_local]);
2650 >= n_irreg_components_filled)
2652 operation.process_empty (local_data[ind_local]);
2656 for (; ind_local<n_local_dofs; ++dof_indices)
2658 Assert (dof_indices != dof_info->end_indices(cell),
2659 ExcInternalError());
2663 operation.process_dof (*dof_indices, *src[0],
2664 local_data[ind_local]);
2667 >= n_irreg_components_filled)
2669 operation.process_empty (local_data[ind_local]);
2679 template <
int dim,
int n_components_,
typename Number>
2680 template<
typename VectorType>
2688 typename internal::BlockVectorSelector<VectorType,
2690 for (
unsigned int d=0; d<n_components; ++d)
2693 internal::VectorReader<Number> reader;
2694 read_write_operation (reader, src_data);
2697 dof_values_initialized =
true;
2703 template <
int dim,
int n_components_,
typename Number>
2704 template<
typename VectorType>
2709 const unsigned int first_index)
2712 Assert (n_fe_components == 1, ExcNotImplemented());
2713 Assert ((n_fe_components == 1 ?
2714 (first_index+n_components <= src.size()) :
true),
2715 ExcIndexRange (first_index + n_components_, 0, src.size()));
2717 VectorType *src_data [n_components];
2718 for (
unsigned int comp=0; comp<n_components; ++comp)
2719 src_data[comp] = const_cast<VectorType *>(&src[comp+first_index]);
2721 internal::VectorReader<Number> reader;
2722 read_write_operation (reader, src_data);
2725 dof_values_initialized =
true;
2731 template <
int dim,
int n_components_,
typename Number>
2732 template<
typename VectorType>
2737 const unsigned int first_index)
2740 Assert (n_fe_components == 1, ExcNotImplemented());
2741 Assert ((n_fe_components == 1 ?
2742 (first_index+n_components <= src.size()) :
true),
2743 ExcIndexRange (first_index + n_components_, 0, src.size()));
2745 const VectorType *src_data [n_components];
2746 for (
unsigned int comp=0; comp<n_components; ++comp)
2747 src_data[comp] = const_cast<VectorType *>(src[comp+first_index]);
2749 internal::VectorReader<Number> reader;
2750 read_write_operation (reader, src_data);
2753 dof_values_initialized =
true;
2759 template <
int dim,
int n_components_,
typename Number>
2760 template<
typename VectorType>
2768 const typename internal::BlockVectorSelector<VectorType,
2770 for (
unsigned int d=0; d<n_components; ++d)
2773 read_dof_values_plain (src_data);
2778 template <
int dim,
int n_components_,
typename Number>
2779 template<
typename VectorType>
2784 const unsigned int first_index)
2787 Assert (n_fe_components == 1, ExcNotImplemented());
2788 Assert ((n_fe_components == 1 ?
2789 (first_index+n_components <= src.size()) :
true),
2790 ExcIndexRange (first_index + n_components_, 0, src.size()));
2791 const VectorType *src_data [n_components];
2792 for (
unsigned int comp=0; comp<n_components; ++comp)
2793 src_data[comp] = &src[comp+first_index];
2794 read_dof_values_plain (src_data);
2799 template <
int dim,
int n_components_,
typename Number>
2800 template<
typename VectorType>
2805 const unsigned int first_index)
2808 Assert (n_fe_components == 1, ExcNotImplemented());
2809 Assert ((n_fe_components == 1 ?
2810 (first_index+n_components <= src.size()) :
true),
2811 ExcIndexRange (first_index + n_components_, 0, src.size()));
2812 const VectorType *src_data [n_components];
2813 for (
unsigned int comp=0; comp<n_components; ++comp)
2814 src_data[comp] = src[comp+first_index];
2815 read_dof_values_plain (src_data);
2820 template <
int dim,
int n_components_,
typename Number>
2821 template<
typename VectorType>
2827 Assert (dof_values_initialized==
true,
2828 internal::ExcAccessToUninitializedField());
2832 typename internal::BlockVectorSelector<VectorType,
2834 for (
unsigned int d=0; d<n_components; ++d)
2837 internal::VectorDistributorLocalToGlobal<Number> distributor;
2838 read_write_operation (distributor, dst_data);
2843 template <
int dim,
int n_components_,
typename Number>
2844 template<
typename VectorType>
2849 const unsigned int first_index)
const 2852 Assert (n_fe_components == 1, ExcNotImplemented());
2853 Assert ((n_fe_components == 1 ?
2854 (first_index+n_components <= dst.size()) :
true),
2855 ExcIndexRange (first_index + n_components_, 0, dst.size()));
2856 Assert (dof_values_initialized==
true,
2857 internal::ExcAccessToUninitializedField());
2859 VectorType *dst_data [n_components];
2860 for (
unsigned int comp=0; comp<n_components; ++comp)
2861 dst_data[comp] = &dst[comp+first_index];
2863 internal::VectorDistributorLocalToGlobal<Number> distributor;
2864 read_write_operation (distributor, dst_data);
2869 template <
int dim,
int n_components_,
typename Number>
2870 template<
typename VectorType>
2875 const unsigned int first_index)
const 2878 Assert (n_fe_components == 1, ExcNotImplemented());
2879 Assert ((n_fe_components == 1 ?
2880 (first_index+n_components <= dst.size()) :
true),
2881 ExcIndexRange (first_index + n_components_, 0, dst.size()));
2882 Assert (dof_values_initialized==
true,
2883 internal::ExcAccessToUninitializedField());
2885 VectorType *dst_data [n_components];
2886 for (
unsigned int comp=0; comp<n_components; ++comp)
2887 dst_data[comp] = dst[comp+first_index];
2889 internal::VectorDistributorLocalToGlobal<Number> distributor;
2890 read_write_operation (distributor, dst_data);
2895 template <
int dim,
int n_components_,
typename Number>
2896 template<
typename VectorType>
2902 Assert (dof_values_initialized==
true,
2903 internal::ExcAccessToUninitializedField());
2907 typename internal::BlockVectorSelector<VectorType,
2909 for (
unsigned int d=0; d<n_components; ++d)
2912 internal::VectorSetter<Number> setter;
2913 read_write_operation (setter, dst_data);
2918 template <
int dim,
int n_components_,
typename Number>
2919 template<
typename VectorType>
2924 const unsigned int first_index)
const 2927 Assert (n_fe_components == 1, ExcNotImplemented());
2928 Assert ((n_fe_components == 1 ?
2929 (first_index+n_components <= dst.size()) :
true),
2930 ExcIndexRange (first_index + n_components_, 0, dst.size()));
2932 Assert (dof_values_initialized==
true,
2933 internal::ExcAccessToUninitializedField());
2935 VectorType *dst_data [n_components];
2936 for (
unsigned int comp=0; comp<n_components; ++comp)
2937 dst_data[comp] = &dst[comp+first_index];
2939 internal::VectorSetter<Number> setter;
2940 read_write_operation (setter, dst_data);
2945 template <
int dim,
int n_components_,
typename Number>
2946 template<
typename VectorType>
2951 const unsigned int first_index)
const 2954 Assert (n_fe_components == 1, ExcNotImplemented());
2955 Assert ((n_fe_components == 1 ?
2956 (first_index+n_components <= dst.size()) :
true),
2957 ExcIndexRange (first_index + n_components_, 0, dst.size()));
2959 Assert (dof_values_initialized==
true,
2960 internal::ExcAccessToUninitializedField());
2962 VectorType *dst_data [n_components];
2963 for (
unsigned int comp=0; comp<n_components; ++comp)
2964 dst_data[comp] = dst[comp+first_index];
2966 internal::VectorSetter<Number> setter;
2967 read_write_operation (setter, dst_data);
2972 template <
int dim,
int n_components_,
typename Number>
2973 template<
typename VectorType>
2980 if (matrix_info == 0)
2982 internal::VectorReader<Number> reader;
2983 read_write_operation (reader, src);
2989 Assert (dof_info != 0, ExcNotInitialized());
2991 ExcNotInitialized());
2993 Assert (dof_info->store_plain_indices ==
true, ExcNotInitialized());
2998 const unsigned int *dof_indices = dof_info->begin_indices_plain(cell);
2999 const unsigned int dofs_per_cell = this->data->
dofs_per_cell;
3001 const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
3002 const bool at_irregular_cell = n_irreg_components_filled > 0;
3006 if (n_fe_components == 1)
3008 const unsigned int n_local_dofs =
3010 for (
unsigned int comp=0; comp<n_components; ++comp)
3011 internal::check_vector_compatibility (*src[comp], *dof_info);
3012 Number *local_src_number [n_components];
3013 for (
unsigned int comp=0; comp<n_components; ++comp)
3014 local_src_number[comp] = &values_dofs[comp][0][0];
3018 if (at_irregular_cell ==
false)
3020 for (
unsigned int j=0; j<n_local_dofs; ++j)
3021 for (
unsigned int comp=0; comp<n_components; ++comp)
3022 local_src_number[comp][j] =
3023 internal::vector_access (*src[comp], dof_indices[j]);
3031 Assert (n_irreg_components_filled > 0, ExcInternalError());
3032 for (
unsigned int ind_local=0; ind_local<n_local_dofs;
3037 for (
unsigned int comp=0; comp<n_components; ++comp)
3038 local_src_number[comp][ind_local] =
3039 internal::vector_access (*src[comp], *dof_indices);
3043 for (
unsigned int comp=0; comp<n_components; ++comp)
3044 local_src_number[comp][ind_local] = 0.;
3056 internal::check_vector_compatibility (*src[0], *dof_info);
3057 Assert (n_fe_components == n_components_, ExcNotImplemented());
3058 const unsigned int n_local_dofs =
3060 Number *local_src_number = &values_dofs[0][0][0];
3061 if (at_irregular_cell ==
false)
3063 for (
unsigned int j=0; j<n_local_dofs; ++j)
3064 local_src_number[j] =
3065 internal::vector_access (*src[0], dof_indices[j]);
3073 Assert (n_irreg_components_filled > 0, ExcInternalError());
3074 for (
unsigned int ind_local=0; ind_local<n_local_dofs; ++dof_indices)
3078 local_src_number[ind_local] =
3079 internal::vector_access (*src[0], *dof_indices);
3083 local_src_number[ind_local] = 0.;
3091 dof_values_initialized =
true;
3100 template <
int dim,
int n_components,
typename Number>
3102 const std::vector<unsigned int> &
3112 template <
int dim,
int n_components,
typename Number>
3118 return &values_dofs[0][0];
3123 template <
int dim,
int n_components,
typename Number>
3130 dof_values_initialized =
true;
3132 return &values_dofs[0][0];
3137 template <
int dim,
int n_components,
typename Number>
3143 Assert (values_quad_initialized || values_quad_submitted,
3144 ExcNotInitialized());
3145 return &values_quad[0][0];
3150 template <
int dim,
int n_components,
typename Number>
3157 values_quad_submitted =
true;
3159 return &values_quad[0][0];
3164 template <
int dim,
int n_components,
typename Number>
3170 Assert (gradients_quad_initialized || gradients_quad_submitted,
3171 ExcNotInitialized());
3172 return &gradients_quad[0][0][0];
3177 template <
int dim,
int n_components,
typename Number>
3184 gradients_quad_submitted =
true;
3186 return &gradients_quad[0][0][0];
3191 template <
int dim,
int n_components,
typename Number>
3197 Assert (hessians_quad_initialized, ExcNotInitialized());
3198 return &hessians_quad[0][0][0];
3203 template <
int dim,
int n_components,
typename Number>
3209 return &hessians_quad[0][0][0];
3214 template <
int dim,
int n_components_,
typename Number>
3222 for (
unsigned int comp=0; comp<n_components; comp++)
3223 return_value[comp] = this->values_dofs[comp][dof];
3224 return return_value;
3229 template <
int dim,
int n_components_,
typename Number>
3235 Assert (this->values_quad_initialized==
true,
3236 internal::ExcAccessToUninitializedField());
3239 for (
unsigned int comp=0; comp<n_components; comp++)
3240 return_value[comp] = this->values_quad[comp][q_point];
3241 return return_value;
3246 template <
int dim,
int n_components_,
typename Number>
3252 Assert (this->gradients_quad_initialized==
true,
3253 internal::ExcAccessToUninitializedField());
3259 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3261 for (
unsigned int comp=0; comp<n_components; comp++)
3262 for (
unsigned int d=0; d<dim; ++d)
3263 grad_out[comp][d] = (this->gradients_quad[comp][d][q_point] *
3264 cartesian_data[0][d]);
3270 this->cell_type == internal::MatrixFreeFunctions::general ?
3271 jacobian[q_point] : jacobian[0];
3272 for (
unsigned int comp=0; comp<n_components; comp++)
3274 for (
unsigned int d=0; d<dim; ++d)
3276 grad_out[comp][d] = (jac[d][0] *
3277 this->gradients_quad[comp][0][q_point]);
3278 for (
unsigned int e=1; e<dim; ++e)
3279 grad_out[comp][d] += (jac[d][e] *
3280 this->gradients_quad[comp][e][q_point]);
3293 template <
typename Number>
3298 const unsigned int q_point,
3301 tmp[0][0] = jac[0][0] * hessians_quad[0][q_point];
3304 template <
typename Number>
3309 const unsigned int q_point,
3312 for (
unsigned int d=0; d<2; ++d)
3314 tmp[0][d] = (jac[d][0] * hessians_quad[0][q_point] +
3315 jac[d][1] * hessians_quad[2][q_point]);
3316 tmp[1][d] = (jac[d][0] * hessians_quad[2][q_point] +
3317 jac[d][1] * hessians_quad[1][q_point]);
3321 template <
typename Number>
3326 const unsigned int q_point,
3329 for (
unsigned int d=0; d<3; ++d)
3331 tmp[0][d] = (jac[d][0] * hessians_quad[0][q_point] +
3332 jac[d][1] * hessians_quad[3][q_point] +
3333 jac[d][2] * hessians_quad[4][q_point]);
3334 tmp[1][d] = (jac[d][0] * hessians_quad[3][q_point] +
3335 jac[d][1] * hessians_quad[1][q_point] +
3336 jac[d][2] * hessians_quad[5][q_point]);
3337 tmp[2][d] = (jac[d][0] * hessians_quad[4][q_point] +
3338 jac[d][1] * hessians_quad[5][q_point] +
3339 jac[d][2] * hessians_quad[2][q_point]);
3346 template <
int dim,
int n_components_,
typename Number>
3352 Assert (this->hessians_quad_initialized==
true,
3353 internal::ExcAccessToUninitializedField());
3359 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3362 for (
unsigned int comp=0; comp<n_components; comp++)
3363 for (
unsigned int d=0; d<dim; ++d)
3365 hessian_out[comp][d][d] = (this->hessians_quad[comp][d][q_point] *
3372 hessian_out[comp][0][1] = (this->hessians_quad[comp][2][q_point] *
3376 hessian_out[comp][0][1] = (this->hessians_quad[comp][3][q_point] *
3378 hessian_out[comp][0][2] = (this->hessians_quad[comp][4][q_point] *
3380 hessian_out[comp][1][2] = (this->hessians_quad[comp][5][q_point] *
3384 Assert (
false, ExcNotImplemented());
3386 for (
unsigned int e=d+1; e<dim; ++e)
3387 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3391 else if (this->cell_type == internal::MatrixFreeFunctions::general)
3394 ExcNotInitialized());
3399 & jac_grad_UT = jacobian_grad_upper[q_point];
3400 for (
unsigned int comp=0; comp<n_components; comp++)
3405 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3409 for (
unsigned int d=0; d<dim; ++d)
3410 for (
unsigned int e=d; e<dim; ++e)
3412 hessian_out[comp][d][e] = jac[d][0] * tmp[0][e];
3413 for (
unsigned int f=1; f<dim; ++f)
3414 hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3418 for (
unsigned int d=0; d<dim; ++d)
3419 for (
unsigned int e=0; e<dim; ++e)
3420 hessian_out[comp][d][d] += (jac_grad[d][e] *
3421 this->gradients_quad[comp][e][q_point]);
3424 for (
unsigned int d=0, count=0; d<dim; ++d)
3425 for (
unsigned int e=d+1; e<dim; ++e, ++count)
3426 for (
unsigned int f=0; f<dim; ++f)
3427 hessian_out[comp][d][e] += (jac_grad_UT[count][f] *
3428 this->gradients_quad[comp][f][q_point]);
3431 for (
unsigned int d=0; d<dim; ++d)
3432 for (
unsigned int e=d+1; e<dim; ++e)
3433 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3440 for (
unsigned int comp=0; comp<n_components; comp++)
3445 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3449 for (
unsigned int d=0; d<dim; ++d)
3450 for (
unsigned int e=d; e<dim; ++e)
3452 hessian_out[comp][d][e] = jac[d][0] * tmp[0][e];
3453 for (
unsigned int f=1; f<dim; ++f)
3454 hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3461 for (
unsigned int d=0; d<dim; ++d)
3462 for (
unsigned int e=d+1; e<dim; ++e)
3463 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3471 template <
int dim,
int n_components_,
typename Number>
3477 Assert (this->hessians_quad_initialized==
true,
3478 internal::ExcAccessToUninitializedField());
3484 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3487 for (
unsigned int comp=0; comp<n_components; comp++)
3488 for (
unsigned int d=0; d<dim; ++d)
3489 hessian_out[comp][d] = (this->hessians_quad[comp][d][q_point] *
3493 else if (this->cell_type == internal::MatrixFreeFunctions::general)
3496 ExcNotInitialized());
3499 for (
unsigned int comp=0; comp<n_components; comp++)
3504 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3509 for (
unsigned int d=0; d<dim; ++d)
3511 hessian_out[comp][d] = jac[d][0] * tmp[0][d];
3512 for (
unsigned int f=1; f<dim; ++f)
3513 hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3516 for (
unsigned int d=0; d<dim; ++d)
3517 for (
unsigned int e=0; e<dim; ++e)
3518 hessian_out[comp][d] += (jac_grad[d][e] *
3519 this->gradients_quad[comp][e][q_point]);
3526 for (
unsigned int comp=0; comp<n_components; comp++)
3531 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3536 for (
unsigned int d=0; d<dim; ++d)
3538 hessian_out[comp][d] = jac[d][0] * tmp[0][d];
3539 for (
unsigned int f=1; f<dim; ++f)
3540 hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3549 template <
int dim,
int n_components_,
typename Number>
3555 Assert (this->hessians_quad_initialized==
true,
3556 internal::ExcAccessToUninitializedField());
3560 = get_hessian_diagonal(q_point);
3561 for (
unsigned int comp=0; comp<n_components; ++comp)
3563 laplacian_out[comp] = hess_diag[comp][0];
3564 for (
unsigned int d=1; d<dim; ++d)
3565 laplacian_out[comp] += hess_diag[comp][d];
3567 return laplacian_out;
3572 template <
int dim,
int n_components_,
typename Number>
3577 const unsigned int dof)
3580 this->dof_values_initialized =
true;
3583 for (
unsigned int comp=0; comp<n_components; comp++)
3584 this->values_dofs[comp][dof] = val_in[comp];
3589 template <
int dim,
int n_components_,
typename Number>
3594 const unsigned int q_point)
3599 this->values_quad_submitted =
true;
3601 if (this->cell_type == internal::MatrixFreeFunctions::general)
3604 for (
unsigned int comp=0; comp<n_components; ++comp)
3605 this->values_quad[comp][q_point] = val_in[comp] * JxW;
3610 for (
unsigned int comp=0; comp<n_components; ++comp)
3611 this->values_quad[comp][q_point] = val_in[comp] * JxW;
3617 template <
int dim,
int n_components_,
typename Number>
3623 const unsigned int q_point)
3628 this->gradients_quad_submitted =
true;
3630 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3633 for (
unsigned int comp=0; comp<n_components; comp++)
3634 for (
unsigned int d=0; d<dim; ++d)
3635 this->gradients_quad[comp][d][q_point] = (grad_in[comp][d] *
3636 cartesian_data[0][d] * JxW);
3641 this->cell_type == internal::MatrixFreeFunctions::general ?
3642 jacobian[q_point] : jacobian[0];
3644 this->cell_type == internal::MatrixFreeFunctions::general ?
3645 J_value[q_point] : J_value[0] * quadrature_weights[q_point];
3646 for (
unsigned int comp=0; comp<n_components; ++comp)
3647 for (
unsigned int d=0; d<dim; ++d)
3650 for (
unsigned int e=1; e<dim; ++e)
3651 new_val += (jac[e][d] * grad_in[comp][e]);
3652 this->gradients_quad[comp][d][q_point] = new_val * JxW;
3659 template <
int dim,
int n_components_,
typename Number>
3667 Assert (this->values_quad_submitted ==
true,
3668 internal::ExcAccessToUninitializedField());
3671 for (
unsigned int comp=0; comp<n_components; ++comp)
3672 return_value[comp] = this->values_quad[comp][0];
3673 const unsigned int n_q_points = this->data->
n_q_points;
3674 for (
unsigned int q=1; q<n_q_points; ++q)
3675 for (
unsigned int comp=0; comp<n_components; ++comp)
3676 return_value[comp] += this->values_quad[comp][q];
3677 return (return_value);
3685 template <
int dim,
int n_components_,
typename Number>
3689 const unsigned int fe_no,
3690 const unsigned int quad_no_in,
3691 const unsigned int fe_degree,
3692 const unsigned int n_q_points)
3695 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
3700 template <
int dim,
int n_components_,
typename Number>
3701 template <
int n_components_other>
3708 const unsigned int first_selected_component,
3712 first_selected_component, other)
3717 template <
int dim,
int n_components_,
typename Number>
3731 template <
int dim,
typename Number>
3735 const unsigned int fe_no,
3736 const unsigned int quad_no_in,
3737 const unsigned int fe_degree,
3738 const unsigned int n_q_points)
3741 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
3746 template <
int dim,
typename Number>
3747 template <
int n_components_other>
3754 const unsigned int first_selected_component,
3758 first_selected_component, other)
3763 template <
int dim,
typename Number>
3773 template <
int dim,
typename Number>
3780 return this->values_dofs[0][dof];
3785 template <
int dim,
typename Number>
3791 Assert (this->values_quad_initialized==
true,
3792 internal::ExcAccessToUninitializedField());
3794 return this->values_quad[0][q_point];
3799 template <
int dim,
typename Number>
3808 Assert (this->gradients_quad_initialized==
true,
3809 internal::ExcAccessToUninitializedField());
3815 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3817 for (
unsigned int d=0; d<dim; ++d)
3818 grad_out[d] = (this->gradients_quad[0][d][q_point] *
3819 this->cartesian_data[0][d]);
3825 this->cell_type == internal::MatrixFreeFunctions::general ?
3826 this->jacobian[q_point] : this->jacobian[0];
3827 for (
unsigned int d=0; d<dim; ++d)
3829 grad_out[d] = (jac[d][0] * this->gradients_quad[0][0][q_point]);
3830 for (
unsigned int e=1; e<dim; ++e)
3831 grad_out[d] += (jac[d][e] * this->gradients_quad[0][e][q_point]);
3839 template <
int dim,
typename Number>
3845 return BaseClass::get_hessian(q_point)[0];
3850 template <
int dim,
typename Number>
3856 return BaseClass::get_hessian_diagonal(q_point)[0];
3861 template <
int dim,
typename Number>
3867 return BaseClass::get_laplacian(q_point)[0];
3872 template <
int dim,
typename Number>
3877 const unsigned int dof)
3880 this->dof_values_initialized =
true;
3883 this->values_dofs[0][dof] = val_in;
3888 template <
int dim,
typename Number>
3893 const unsigned int q_point)
3898 this->values_quad_submitted =
true;
3900 if (this->cell_type == internal::MatrixFreeFunctions::general)
3903 this->values_quad[0][q_point] = val_in * JxW;
3908 this->values_quad[0][q_point] = val_in * JxW;
3914 template <
int dim,
typename Number>
3919 const unsigned int q_point)
3924 this->gradients_quad_submitted =
true;
3926 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3929 for (
unsigned int d=0; d<dim; ++d)
3930 this->gradients_quad[0][d][q_point] = (grad_in[d] *
3931 this->cartesian_data[0][d] *
3938 this->cell_type == internal::MatrixFreeFunctions::general ?
3939 this->jacobian[q_point] : this->jacobian[0];
3941 this->cell_type == internal::MatrixFreeFunctions::general ?
3942 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
3943 for (
unsigned int d=0; d<dim; ++d)
3946 for (
unsigned int e=1; e<dim; ++e)
3947 new_val += jac[e][d] * grad_in[e];
3948 this->gradients_quad[0][d][q_point] = new_val * JxW;
3955 template <
int dim,
typename Number>
3961 return BaseClass::integrate_value()[0];
3970 template <
int dim,
typename Number>
3974 const unsigned int fe_no,
3975 const unsigned int quad_no_in,
3976 const unsigned int fe_degree,
3977 const unsigned int n_q_points)
3980 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
3985 template <
int dim,
typename Number>
3986 template <
int n_components_other>
3993 const unsigned int first_selected_component,
3997 first_selected_component, other)
4002 template <
int dim,
typename Number>
4012 template <
int dim,
typename Number>
4018 return BaseClass::get_gradient (q_point);
4023 template <
int dim,
typename Number>
4029 Assert (this->gradients_quad_initialized==
true,
4030 internal::ExcAccessToUninitializedField());
4036 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4038 divergence = (this->gradients_quad[0][0][q_point] *
4039 this->cartesian_data[0][0]);
4040 for (
unsigned int d=1; d<dim; ++d)
4041 divergence += (this->gradients_quad[d][d][q_point] *
4042 this->cartesian_data[0][d]);
4048 this->cell_type == internal::MatrixFreeFunctions::general ?
4049 this->jacobian[q_point] : this->jacobian[0];
4050 divergence = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4051 for (
unsigned int e=1; e<dim; ++e)
4052 divergence += (jac[0][e] * this->gradients_quad[0][e][q_point]);
4053 for (
unsigned int d=1; d<dim; ++d)
4054 for (
unsigned int e=0; e<dim; ++e)
4055 divergence += (jac[d][e] * this->gradients_quad[d][e][q_point]);
4062 template <
int dim,
typename Number>
4072 for (
unsigned int d=0; d<dim; ++d)
4073 symmetrized[d] = grad[d][d];
4079 symmetrized[2] = grad[0][1] + grad[1][0];
4080 symmetrized[2] *= half;
4083 symmetrized[3] = grad[0][1] + grad[1][0];
4084 symmetrized[3] *= half;
4085 symmetrized[4] = grad[0][2] + grad[2][0];
4086 symmetrized[4] *= half;
4087 symmetrized[5] = grad[1][2] + grad[2][1];
4088 symmetrized[5] *= half;
4091 Assert (
false, ExcNotImplemented());
4098 template <
int dim,
typename Number>
4102 ::get_curl (
const unsigned int q_point)
const 4111 ExcMessage(
"Computing the curl in 1d is not a useful operation"));
4114 curl[0] = grad[1][0] - grad[0][1];
4117 curl[0] = grad[2][1] - grad[1][2];
4118 curl[1] = grad[0][2] - grad[2][0];
4119 curl[2] = grad[1][0] - grad[0][1];
4122 Assert (
false, ExcNotImplemented());
4129 template <
int dim,
typename Number>
4135 Assert (this->hessians_quad_initialized==
true,
4136 internal::ExcAccessToUninitializedField());
4139 return BaseClass::get_hessian_diagonal (q_point);
4144 template <
int dim,
typename Number>
4150 Assert (this->hessians_quad_initialized==
true,
4151 internal::ExcAccessToUninitializedField());
4153 return BaseClass::get_hessian(q_point);
4158 template <
int dim,
typename Number>
4163 const unsigned int q_point)
4165 BaseClass::submit_gradient (grad_in, q_point);
4170 template <
int dim,
typename Number>
4176 const unsigned int q_point)
4178 BaseClass::submit_gradient(grad_in, q_point);
4183 template <
int dim,
typename Number>
4188 const unsigned int q_point)
4193 this->gradients_quad_submitted =
true;
4195 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4198 this->quadrature_weights[q_point] * div_in;
4199 for (
unsigned int d=0; d<dim; ++d)
4201 this->gradients_quad[d][d][q_point] = (fac *
4202 this->cartesian_data[0][d]);
4203 for (
unsigned int e=d+1; e<dim; ++e)
4213 this->cell_type == internal::MatrixFreeFunctions::general ?
4214 this->jacobian[q_point] : this->jacobian[0];
4216 (this->cell_type == internal::MatrixFreeFunctions::general ?
4217 this->J_value[q_point] : this->J_value[0] *
4218 this->quadrature_weights[q_point]) * div_in;
4219 for (
unsigned int d=0; d<dim; ++d)
4221 for (
unsigned int e=0; e<dim; ++e)
4222 this->gradients_quad[d][e][q_point] = jac[d][e] * fac;
4229 template <
int dim,
typename Number>
4235 const unsigned int q_point)
4243 this->gradients_quad_submitted =
true;
4245 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4248 for (
unsigned int d=0; d<dim; ++d)
4249 this->gradients_quad[d][d][q_point] = (sym_grad.access_raw_entry(d) *
4251 this->cartesian_data[0][d]);
4252 for (
unsigned int e=0, counter=dim; e<dim; ++e)
4253 for (
unsigned int d=e+1; d<dim; ++d, ++counter)
4256 this->gradients_quad[e][d][q_point] = (value *
4257 this->cartesian_data[0][d]);
4258 this->gradients_quad[d][e][q_point] = (value *
4259 this->cartesian_data[0][e]);
4266 this->cell_type == internal::MatrixFreeFunctions::general ?
4267 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4269 this->cell_type == internal::MatrixFreeFunctions::general ?
4270 this->jacobian[q_point] : this->jacobian[0];
4272 for (
unsigned int i=0; i<dim; ++i)
4273 weighted[i][i] = sym_grad.access_raw_entry(i) * JxW;
4274 for (
unsigned int i=0, counter=dim; i<dim; ++i)
4275 for (
unsigned int j=i+1; j<dim; ++j, ++counter)
4278 weighted[i][j] = value;
4279 weighted[j][i] = value;
4281 for (
unsigned int comp=0; comp<dim; ++comp)
4282 for (
unsigned int d=0; d<dim; ++d)
4285 for (
unsigned int e=1; e<dim; ++e)
4286 new_val += jac[e][d] * weighted[comp][e];
4287 this->gradients_quad[comp][d][q_point] = new_val;
4294 template <
int dim,
typename Number>
4299 const unsigned int q_point)
4306 ExcMessage(
"Testing by the curl in 1d is not a useful operation"));
4309 grad[1][0] = curl[0];
4310 grad[0][1] = -curl[0];
4313 grad[2][1] = curl[0];
4314 grad[1][2] = -curl[0];
4315 grad[0][2] = curl[1];
4316 grad[2][0] = -curl[1];
4317 grad[1][0] = curl[2];
4318 grad[0][1] = -curl[2];
4321 Assert (
false, ExcNotImplemented());
4323 submit_gradient (grad, q_point);
4330 template <
typename Number>
4334 const unsigned int fe_no,
4335 const unsigned int quad_no_in,
4336 const unsigned int fe_degree,
4337 const unsigned int n_q_points)
4340 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4345 template <
typename Number>
4346 template <
int n_components_other>
4353 const unsigned int first_selected_component,
4357 first_selected_component, other)
4362 template <
typename Number>
4372 template <
typename Number>
4379 return this->values_dofs[0][dof];
4384 template <
typename Number>
4390 Assert (this->values_quad_initialized==
true,
4391 internal::ExcAccessToUninitializedField());
4393 return this->values_quad[0][q_point];
4398 template <
typename Number>
4407 Assert (this->gradients_quad_initialized==
true,
4408 internal::ExcAccessToUninitializedField());
4414 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4416 grad_out[0] = (this->gradients_quad[0][0][q_point] *
4417 this->cartesian_data[0][0]);
4423 this->cell_type == internal::MatrixFreeFunctions::general ?
4424 this->jacobian[q_point] : this->jacobian[0];
4426 grad_out[0] = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4433 template <
typename Number>
4439 return BaseClass::get_hessian(q_point)[0];
4444 template <
typename Number>
4450 return BaseClass::get_hessian_diagonal(q_point)[0];
4455 template <
typename Number>
4461 return BaseClass::get_laplacian(q_point)[0];
4466 template <
typename Number>
4471 const unsigned int dof)
4474 this->dof_values_initialized =
true;
4477 this->values_dofs[0][dof] = val_in;
4482 template <
typename Number>
4487 const unsigned int q_point)
4492 this->values_quad_submitted =
true;
4494 if (this->cell_type == internal::MatrixFreeFunctions::general)
4497 this->values_quad[0][q_point] = val_in * JxW;
4502 this->values_quad[0][q_point] = val_in * JxW;
4508 template <
typename Number>
4513 const unsigned int q_point)
4518 this->gradients_quad_submitted =
true;
4520 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4523 this->gradients_quad[0][0][q_point] = (grad_in[0] *
4524 this->cartesian_data[0][0] *
4531 this->cell_type == internal::MatrixFreeFunctions::general ?
4532 this->jacobian[q_point] : this->jacobian[0];
4534 this->cell_type == internal::MatrixFreeFunctions::general ?
4535 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4537 this->gradients_quad[0][0][q_point] = jac[0][0] * grad_in[0] * JxW;
4543 template <
typename Number>
4549 return BaseClass::integrate_value()[0];
4560 enum EvaluatorVariant
4570 template <EvaluatorVariant variant,
int dim,
int fe_degree,
int n_q_points_1d,
4572 struct EvaluatorTensorProduct
4579 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
4580 struct EvaluatorTensorProduct<evaluate_general,dim,fe_degree,n_q_points_1d,Number>
4589 EvaluatorTensorProduct ()
4592 shape_gradients (0),
4603 shape_values (shape_values.
begin()),
4604 shape_gradients (shape_gradients.
begin()),
4605 shape_hessians (shape_hessians.
begin())
4608 template <
int direction,
bool dof_to_quad,
bool add>
4610 values (
const Number in [],
4613 apply<direction,dof_to_quad,add>(shape_values, in, out);
4616 template <
int direction,
bool dof_to_quad,
bool add>
4618 gradients (
const Number in [],
4621 apply<direction,dof_to_quad,add>(shape_gradients, in, out);
4624 template <
int direction,
bool dof_to_quad,
bool add>
4626 hessians (
const Number in [],
4629 apply<direction,dof_to_quad,add>(shape_hessians, in, out);
4632 template <
int direction,
bool dof_to_quad,
bool add>
4633 static void apply (
const Number *shape_data,
4637 const Number *shape_values;
4638 const Number *shape_gradients;
4639 const Number *shape_hessians;
4646 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
4647 template <
int direction,
bool dof_to_quad,
bool add>
4650 EvaluatorTensorProduct<evaluate_general,dim,fe_degree,n_q_points_1d,Number>
4651 ::apply (
const Number *shape_data,
4656 const int mm = dof_to_quad ? (fe_degree+1) : n_q_points_1d,
4657 nn = dof_to_quad ? n_q_points_1d : (fe_degree+1);
4659 const int n_blocks1 = (dim > 1 ? (direction > 0 ? nn : mm) : 1);
4660 const int n_blocks2 = (dim > 2 ? (direction > 1 ? nn : mm) : 1);
4663 for (
int i2=0; i2<n_blocks2; ++i2)
4665 for (
int i1=0; i1<n_blocks1; ++i1)
4667 for (
int col=0; col<nn; ++col)
4670 if (dof_to_quad ==
true)
4671 val0 = shape_data[col];
4673 val0 = shape_data[col*n_q_points_1d];
4674 Number res0 = val0 * in[0];
4675 for (
int ind=1; ind<mm; ++ind)
4677 if (dof_to_quad ==
true)
4678 val0 = shape_data[ind*n_q_points_1d+col];
4680 val0 = shape_data[col*n_q_points_1d+ind];
4681 res0 += val0 * in[stride*ind];
4684 out[stride*col] = res0;
4686 out[stride*col] += res0;
4704 Assert (
false, ExcNotImplemented());
4722 template <
int dim,
int fe_degree,
typename Number,
int face_direction,
4723 bool dof_to_quad,
bool add>
4726 apply_tensor_product_face (
const Number *shape_data,
4730 const int n_blocks1 = dim > 1 ? (fe_degree+1) : 1;
4731 const int n_blocks2 = dim > 2 ? (fe_degree+1) : 1;
4734 const int mm = dof_to_quad ? (fe_degree+1) : 1,
4735 nn = dof_to_quad ? 1 : (fe_degree+1);
4739 for (
int i2=0; i2<n_blocks2; ++i2)
4741 for (
int i1=0; i1<n_blocks1; ++i1)
4743 if (dof_to_quad ==
true)
4745 Number res0 = shape_data[0] * in[0];
4746 for (
int ind=1; ind<mm; ++ind)
4747 res0 += shape_data[ind] * in[stride*ind];
4755 for (
int col=0; col<nn; ++col)
4757 out[col*stride] = shape_data[col] * in[0];
4759 out[col*stride] += shape_data[col] * in[0];
4765 switch (face_direction)
4790 Assert (
false, ExcNotImplemented());
4793 if (face_direction == 1 && dim == 3)
4799 out -= (fe_degree+1)*(fe_degree+1)-1;
4801 in -= (fe_degree+1)*(fe_degree+1)-1;
4811 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
4812 struct EvaluatorTensorProduct<evaluate_symmetric,dim,fe_degree,n_q_points_1d,Number>
4824 shape_values (shape_values.
begin()),
4825 shape_gradients (shape_gradients.
begin()),
4826 shape_hessians (shape_hessians.
begin())
4829 template <
int direction,
bool dof_to_quad,
bool add>
4831 values (
const Number in [],
4832 Number out[])
const;
4834 template <
int direction,
bool dof_to_quad,
bool add>
4836 gradients (
const Number in [],
4837 Number out[])
const;
4839 template <
int direction,
bool dof_to_quad,
bool add>
4841 hessians (
const Number in [],
4842 Number out[])
const;
4844 const Number *shape_values;
4845 const Number *shape_gradients;
4846 const Number *shape_hessians;
4869 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
4870 template <
int direction,
bool dof_to_quad,
bool add>
4873 EvaluatorTensorProduct<evaluate_symmetric,dim,fe_degree,n_q_points_1d,Number>
4874 ::values (
const Number in [],
4875 Number out [])
const 4878 const int mm = dof_to_quad ? (fe_degree+1) : n_q_points_1d,
4879 nn = dof_to_quad ? n_q_points_1d : (fe_degree+1);
4880 const int n_cols = nn / 2;
4881 const int mid = mm / 2;
4883 const int n_blocks1 = (dim > 1 ? (direction > 0 ? nn : mm) : 1);
4884 const int n_blocks2 = (dim > 2 ? (direction > 1 ? nn : mm) : 1);
4887 for (
int i2=0; i2<n_blocks2; ++i2)
4889 for (
int i1=0; i1<n_blocks1; ++i1)
4891 for (
int col=0; col<n_cols; ++col)
4893 Number val0, val1, in0, in1, res0, res1;
4894 if (dof_to_quad ==
true)
4896 val0 = shape_values[col];
4897 val1 = shape_values[nn-1-col];
4901 val0 = shape_values[col*n_q_points_1d];
4902 val1 = shape_values[(col+1)*n_q_points_1d-1];
4907 in1 = in[stride*(mm-1)];
4912 for (
int ind=1; ind<mid; ++ind)
4914 if (dof_to_quad ==
true)
4916 val0 = shape_values[ind*n_q_points_1d+col];
4917 val1 = shape_values[ind*n_q_points_1d+nn-1-col];
4921 val0 = shape_values[col*n_q_points_1d+ind];
4922 val1 = shape_values[(col+1)*n_q_points_1d-1-ind];
4924 in0 = in[stride*ind];
4925 in1 = in[stride*(mm-1-ind)];
4933 res0 = res1 = Number();
4934 if (dof_to_quad ==
true)
4938 val0 = shape_values[mid*n_q_points_1d+col];
4939 val1 = val0 * in[stride*mid];
4946 if (mm % 2 == 1 && nn % 2 == 0)
4948 val0 = shape_values[col*n_q_points_1d+mid];
4949 val1 = val0 * in[stride*mid];
4956 out[stride*col] = res0;
4957 out[stride*(nn-1-col)] = res1;
4961 out[stride*col] += res0;
4962 out[stride*(nn-1-col)] += res1;
4965 if ( dof_to_quad ==
true && nn%2==1 && mm%2==1 )
4968 out[stride*n_cols] = in[stride*mid];
4970 out[stride*n_cols] += in[stride*mid];
4972 else if (dof_to_quad ==
true && nn%2==1)
4975 Number val0 = shape_values[n_cols];
4978 res0 = in[0] + in[stride*(mm-1)];
4980 for (
int ind=1; ind<mid; ++ind)
4982 val0 = shape_values[ind*n_q_points_1d+n_cols];
4983 Number val1 = in[stride*ind] + in[stride*(mm-1-ind)];
4991 out[stride*n_cols] = res0;
4993 out[stride*n_cols] += res0;
4995 else if (dof_to_quad ==
false && nn%2 == 1)
5000 Number val0 = shape_values[n_cols*n_q_points_1d];
5001 res0 = in[0] + in[stride*(mm-1)];
5003 for (
int ind=1; ind<mid; ++ind)
5005 val0 = shape_values[n_cols*n_q_points_1d+ind];
5006 Number val1 = in[stride*ind] + in[stride*(mm-1-ind)];
5011 res0 += in[stride*mid];
5016 out[stride*n_cols] = res0;
5018 out[stride*n_cols] += res0;
5036 Assert (
false, ExcNotImplemented());
5068 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
5069 template <
int direction,
bool dof_to_quad,
bool add>
5072 EvaluatorTensorProduct<evaluate_symmetric,dim,fe_degree,n_q_points_1d,Number>
5073 ::gradients (
const Number in [],
5074 Number out [])
const 5077 const int mm = dof_to_quad ? (fe_degree+1) : n_q_points_1d,
5078 nn = dof_to_quad ? n_q_points_1d : (fe_degree+1);
5079 const int n_cols = nn / 2;
5080 const int mid = mm / 2;
5082 const int n_blocks1 = (dim > 1 ? (direction > 0 ? nn : mm) : 1);
5083 const int n_blocks2 = (dim > 2 ? (direction > 1 ? nn : mm) : 1);
5086 for (
int i2=0; i2<n_blocks2; ++i2)
5088 for (
int i1=0; i1<n_blocks1; ++i1)
5090 for (
int col=0; col<n_cols; ++col)
5092 Number val0, val1, in0, in1, res0, res1;
5093 if (dof_to_quad ==
true)
5095 val0 = shape_gradients[col];
5096 val1 = shape_gradients[nn-1-col];
5100 val0 = shape_gradients[col*n_q_points_1d];
5101 val1 = shape_gradients[(nn-col-1)*n_q_points_1d];
5106 in1 = in[stride*(mm-1)];
5111 for (
int ind=1; ind<mid; ++ind)
5113 if (dof_to_quad ==
true)
5115 val0 = shape_gradients[ind*n_q_points_1d+col];
5116 val1 = shape_gradients[ind*n_q_points_1d+nn-1-col];
5120 val0 = shape_gradients[col*n_q_points_1d+ind];
5121 val1 = shape_gradients[(nn-col-1)*n_q_points_1d+ind];
5123 in0 = in[stride*ind];
5124 in1 = in[stride*(mm-1-ind)];
5132 res0 = res1 = Number();
5135 if (dof_to_quad ==
true)
5136 val0 = shape_gradients[mid*n_q_points_1d+col];
5138 val0 = shape_gradients[col*n_q_points_1d+mid];
5139 val1 = val0 * in[stride*mid];
5145 out[stride*col] = res0;
5146 out[stride*(nn-1-col)] = res1;
5150 out[stride*col] += res0;
5151 out[stride*(nn-1-col)] += res1;
5157 if (dof_to_quad ==
true)
5158 val0 = shape_gradients[n_cols];
5160 val0 = shape_gradients[n_cols*n_q_points_1d];
5161 res0 = in[0] - in[stride*(mm-1)];
5163 for (
int ind=1; ind<mid; ++ind)
5165 if (dof_to_quad ==
true)
5166 val0 = shape_gradients[ind*n_q_points_1d+n_cols];
5168 val0 = shape_gradients[n_cols*n_q_points_1d+ind];
5169 Number val1 = in[stride*ind] - in[stride*(mm-1-ind)];
5174 out[stride*n_cols] = res0;
5176 out[stride*n_cols] += res0;
5195 Assert (
false, ExcNotImplemented());
5212 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
5213 template <
int direction,
bool dof_to_quad,
bool add>
5216 EvaluatorTensorProduct<evaluate_symmetric,dim,fe_degree,n_q_points_1d,Number>
5217 ::hessians (
const Number in [],
5218 Number out [])
const 5221 const int mm = dof_to_quad ? (fe_degree+1) : n_q_points_1d,
5222 nn = dof_to_quad ? n_q_points_1d : (fe_degree+1);
5223 const int n_cols = nn / 2;
5224 const int mid = mm / 2;
5226 const int n_blocks1 = (dim > 1 ? (direction > 0 ? nn : mm) : 1);
5227 const int n_blocks2 = (dim > 2 ? (direction > 1 ? nn : mm) : 1);
5230 for (
int i2=0; i2<n_blocks2; ++i2)
5232 for (
int i1=0; i1<n_blocks1; ++i1)
5234 for (
int col=0; col<n_cols; ++col)
5236 Number val0, val1, in0, in1, res0, res1;
5237 if (dof_to_quad ==
true)
5239 val0 = shape_hessians[col];
5240 val1 = shape_hessians[nn-1-col];
5244 val0 = shape_hessians[col*n_q_points_1d];
5245 val1 = shape_hessians[(col+1)*n_q_points_1d-1];
5250 in1 = in[stride*(mm-1)];
5255 for (
int ind=1; ind<mid; ++ind)
5257 if (dof_to_quad ==
true)
5259 val0 = shape_hessians[ind*n_q_points_1d+col];
5260 val1 = shape_hessians[ind*n_q_points_1d+nn-1-col];
5264 val0 = shape_hessians[col*n_q_points_1d+ind];
5265 val1 = shape_hessians[(col+1)*n_q_points_1d-1-ind];
5267 in0 = in[stride*ind];
5268 in1 = in[stride*(mm-1-ind)];
5276 res0 = res1 = Number();
5279 if (dof_to_quad ==
true)
5280 val0 = shape_hessians[mid*n_q_points_1d+col];
5282 val0 = shape_hessians[col*n_q_points_1d+mid];
5283 val1 = val0 * in[stride*mid];
5289 out[stride*col] = res0;
5290 out[stride*(nn-1-col)] = res1;
5294 out[stride*col] += res0;
5295 out[stride*(nn-1-col)] += res1;
5301 if (dof_to_quad ==
true)
5302 val0 = shape_hessians[n_cols];
5304 val0 = shape_hessians[n_cols*n_q_points_1d];
5307 res0 = in[0] + in[stride*(mm-1)];
5309 for (
int ind=1; ind<mid; ++ind)
5311 if (dof_to_quad ==
true)
5312 val0 = shape_hessians[ind*n_q_points_1d+n_cols];
5314 val0 = shape_hessians[n_cols*n_q_points_1d+ind];
5315 Number val1 = in[stride*ind] + in[stride*(mm-1-ind)];
5324 if (dof_to_quad ==
true)
5325 val0 = shape_hessians[mid*n_q_points_1d+n_cols];
5327 val0 = shape_hessians[n_cols*n_q_points_1d+mid];
5328 res0 += val0 * in[stride*mid];
5331 out[stride*n_cols] = res0;
5333 out[stride*n_cols] += res0;
5351 Assert (
false, ExcNotImplemented());
5377 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
5378 struct EvaluatorTensorProduct<evaluate_evenodd,dim,fe_degree,n_q_points_1d,Number>
5387 EvaluatorTensorProduct ()
5390 shape_gradients (0),
5402 shape_values (shape_values.
begin()),
5403 shape_gradients (shape_gradients.
begin()),
5404 shape_hessians (shape_hessians.
begin())
5407 template <
int direction,
bool dof_to_quad,
bool add>
5409 values (
const Number in [],
5412 apply<direction,dof_to_quad,add,0>(shape_values, in, out);
5415 template <
int direction,
bool dof_to_quad,
bool add>
5417 gradients (
const Number in [],
5420 apply<direction,dof_to_quad,add,1>(shape_gradients, in, out);
5423 template <
int direction,
bool dof_to_quad,
bool add>
5425 hessians (
const Number in [],
5428 apply<direction,dof_to_quad,add,2>(shape_hessians, in, out);
5431 template <
int direction,
bool dof_to_quad,
bool add,
int type>
5432 static void apply (
const Number *shape_data,
5436 const Number *shape_values;
5437 const Number *shape_gradients;
5438 const Number *shape_hessians;
5443 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
typename Number>
5444 template <
int direction,
bool dof_to_quad,
bool add,
int type>
5447 EvaluatorTensorProduct<evaluate_evenodd,dim,fe_degree,n_q_points_1d,Number>
5448 ::apply (
const Number *shapes,
5454 const int mm = dof_to_quad ? (fe_degree+1) : n_q_points_1d,
5455 nn = dof_to_quad ? n_q_points_1d : (fe_degree+1);
5456 const int n_cols = nn / 2;
5457 const int mid = mm / 2;
5459 const int n_blocks1 = (dim > 1 ? (direction > 0 ? nn : mm) : 1);
5460 const int n_blocks2 = (dim > 2 ? (direction > 1 ? nn : mm) : 1);
5463 const int offset = (n_q_points_1d+1)/2;
5469 for (
int i2=0; i2<n_blocks2; ++i2)
5471 for (
int i1=0; i1<n_blocks1; ++i1)
5473 Number xp[mid>0?mid:1], xm[mid>0?mid:1];
5474 for (
int i=0; i<mid; ++i)
5476 if (dof_to_quad ==
true && type == 1)
5478 xp[i] = in[stride*i] - in[stride*(mm-1-i)];
5479 xm[i] = in[stride*i] + in[stride*(mm-1-i)];
5483 xp[i] = in[stride*i] + in[stride*(mm-1-i)];
5484 xm[i] = in[stride*i] - in[stride*(mm-1-i)];
5487 for (
int col=0; col<n_cols; ++col)
5492 if (dof_to_quad ==
true)
5494 r0 = shapes[col] * xp[0];
5495 r1 = shapes[fe_degree*offset + col] * xm[0];
5499 r0 = shapes[col*offset] * xp[0];
5500 r1 = shapes[(fe_degree-col)*offset] * xm[0];
5502 for (
int ind=1; ind<mid; ++ind)
5504 if (dof_to_quad ==
true)
5506 r0 += shapes[ind*offset+col] * xp[ind];
5507 r1 += shapes[(fe_degree-ind)*offset+col] * xm[ind];
5511 r0 += shapes[col*offset+ind] * xp[ind];
5512 r1 += shapes[(fe_degree-col)*offset+ind] * xm[ind];
5518 if (mm % 2 == 1 && dof_to_quad ==
true)
5521 r1 += shapes[mid*offset+col] * in[stride*mid];
5523 r0 += shapes[mid*offset+col] * in[stride*mid];
5525 else if (mm % 2 == 1 && (nn % 2 == 0 || type > 0))
5526 r0 += shapes[col*offset+mid] * in[stride*mid];
5530 out[stride*col] = r0 + r1;
5531 if (type == 1 && dof_to_quad ==
false)
5532 out[stride*(nn-1-col)] = r1 - r0;
5534 out[stride*(nn-1-col)] = r0 - r1;
5538 out[stride*col] += r0 + r1;
5539 if (type == 1 && dof_to_quad ==
false)
5540 out[stride*(nn-1-col)] += r1 - r0;
5542 out[stride*(nn-1-col)] += r0 - r1;
5545 if ( type == 0 && dof_to_quad ==
true && nn%2==1 && mm%2==1 )
5548 out[stride*n_cols] = in[stride*mid];
5550 out[stride*n_cols] += in[stride*mid];
5552 else if (dof_to_quad ==
true && nn%2==1)
5557 r0 = shapes[n_cols] * xp[0];
5558 for (
int ind=1; ind<mid; ++ind)
5559 r0 += shapes[ind*offset+n_cols] * xp[ind];
5563 if (type != 1 && mm % 2 == 1)
5564 r0 += shapes[mid*offset+n_cols] * in[stride*mid];
5567 out[stride*n_cols] = r0;
5569 out[stride*n_cols] += r0;
5571 else if (dof_to_quad ==
false && nn%2 == 1)
5578 r0 = shapes[n_cols*offset] * xm[0];
5579 for (
int ind=1; ind<mid; ++ind)
5580 r0 += shapes[n_cols*offset+ind] * xm[ind];
5584 r0 = shapes[n_cols*offset] * xp[0];
5585 for (
int ind=1; ind<mid; ++ind)
5586 r0 += shapes[n_cols*offset+ind] * xp[ind];
5592 if (type == 0 && mm % 2 == 1)
5593 r0 += in[stride*mid];
5594 else if (type == 2 && mm % 2 == 1)
5595 r0 += shapes[n_cols*offset+mid] * in[stride*mid];
5598 out[stride*n_cols] = r0;
5600 out[stride*n_cols] += r0;
5618 Assert (
false, ExcNotImplemented());
5632 template <MatrixFreeFunctions::ElementType element,
bool is_
long>
5633 struct EvaluatorSelector {};
5635 template <
bool is_
long>
5636 struct EvaluatorSelector<MatrixFreeFunctions::tensor_general,is_long>
5638 static const EvaluatorVariant variant = evaluate_general;
5642 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric,false>
5644 static const EvaluatorVariant variant = evaluate_symmetric;
5647 template <>
struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric,true>
5649 static const EvaluatorVariant variant = evaluate_evenodd;
5652 template <
bool is_
long>
5653 struct EvaluatorSelector<MatrixFreeFunctions::truncated_tensor,is_long>
5655 static const EvaluatorVariant variant = evaluate_general;
5659 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric_plus_dg0,false>
5661 static const EvaluatorVariant variant = evaluate_symmetric;
5665 struct EvaluatorSelector<MatrixFreeFunctions::tensor_symmetric_plus_dg0,true>
5667 static const EvaluatorVariant variant = evaluate_evenodd;
5670 template <
bool is_
long>
5671 struct EvaluatorSelector<MatrixFreeFunctions::tensor_gausslobatto,is_long>
5673 static const EvaluatorVariant variant = evaluate_evenodd;
5689 template <MatrixFreeFunctions::ElementType type,
int dim,
int fe_degree,
5690 int n_q_points_1d,
int n_components,
typename Number>
5691 struct FEEvaluationImpl
5694 void evaluate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
5699 const bool evaluate_val,
5700 const bool evaluate_grad,
5701 const bool evaluate_lapl);
5704 void integrate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
5708 const bool evaluate_val,
5709 const bool evaluate_grad);
5713 template <MatrixFreeFunctions::ElementType type,
int dim,
int fe_degree,
5714 int n_q_points_1d,
int n_components,
typename Number>
5717 FEEvaluationImpl<type,dim,fe_degree,n_q_points_1d,n_components,Number>
5718 ::evaluate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
5723 const bool evaluate_val,
5724 const bool evaluate_grad,
5725 const bool evaluate_lapl)
5727 if (evaluate_val ==
false && evaluate_grad ==
false && evaluate_lapl ==
false)
5730 const EvaluatorVariant variant =
5731 EvaluatorSelector<type,(fe_degree+n_q_points_1d>4)>::variant;
5732 typedef EvaluatorTensorProduct<variant, dim, fe_degree, n_q_points_1d,
5734 Eval eval (variant == evaluate_evenodd ? shape_info.shape_val_evenodd :
5735 shape_info.shape_values,
5736 variant == evaluate_evenodd ? shape_info.shape_gra_evenodd :
5737 shape_info.shape_gradients,
5738 variant == evaluate_evenodd ? shape_info.shape_hes_evenodd :
5739 shape_info.shape_hessians);
5741 const unsigned int temp_size = Eval::dofs_per_cell > Eval::n_q_points ?
5742 Eval::dofs_per_cell : Eval::n_q_points;
5746 n_components*Eval::dofs_per_cell];
5748 if (type == MatrixFreeFunctions::truncated_tensor)
5750 for (
unsigned int c=0; c<n_components; ++c)
5751 expanded_dof_values[c] = &data_array[c*Eval::dofs_per_cell];
5752 values_dofs = expanded_dof_values;
5754 unsigned int count_p = 0, count_q = 0;
5755 for (
unsigned int i=0; i<(dim>2?fe_degree+1:1); ++i)
5757 for (
unsigned int j=0; j<(dim>1?fe_degree+1-i:1); ++j)
5759 for (
unsigned int k=0; k<fe_degree+1-j-i; ++k, ++count_p, ++count_q)
5760 for (
unsigned int c=0; c<n_components; ++c)
5761 expanded_dof_values[c][count_q] = values_dofs_actual[c][count_p];
5762 for (
unsigned int k=fe_degree+1-j-i; k<fe_degree+1; ++k, ++count_q)
5763 for (
unsigned int c=0; c<n_components; ++c)
5766 for (
unsigned int j=fe_degree+1-i; j<fe_degree+1; ++j)
5767 for (
unsigned int k=0; k<fe_degree+1; ++k, ++count_q)
5768 for (
unsigned int c=0; c<n_components; ++c)
5777 const unsigned int d1 = dim>1?1:0;
5778 const unsigned int d2 = dim>2?2:0;
5779 const unsigned int d3 = dim>2?3:0;
5780 const unsigned int d4 = dim>2?4:0;
5781 const unsigned int d5 = dim>2?5:0;
5786 for (
unsigned int c=0; c<n_components; c++)
5788 if (evaluate_val ==
true)
5789 eval.template values<0,true,false> (values_dofs[c], values_quad[c]);
5790 if (evaluate_grad ==
true)
5791 eval.template gradients<0,true,false>(values_dofs[c], gradients_quad[c][0]);
5792 if (evaluate_lapl ==
true)
5793 eval.template hessians<0,true,false> (values_dofs[c], hessians_quad[c][0]);
5798 for (
unsigned int c=0; c<n_components; c++)
5803 if (evaluate_grad ==
true)
5805 eval.template gradients<0,true,false> (values_dofs[c], temp1);
5806 eval.template values<1,true,false> (temp1, gradients_quad[c][0]);
5808 if (evaluate_lapl ==
true)
5811 if (evaluate_grad ==
false)
5812 eval.template gradients<0,true,false>(values_dofs[c], temp1);
5813 eval.template gradients<1,true,false> (temp1, hessians_quad[c][d1+d1]);
5816 eval.template hessians<0,true,false>(values_dofs[c], temp1);
5817 eval.template values<1,true,false> (temp1, hessians_quad[c][0]);
5821 eval.template values<0,true,false> (values_dofs[c], temp1);
5822 if (evaluate_grad ==
true)
5823 eval.template gradients<1,true,false> (temp1, gradients_quad[c][d1]);
5826 if (evaluate_lapl ==
true)
5827 eval.template hessians<1,true,false> (temp1, hessians_quad[c][d1]);
5830 if (evaluate_val ==
true)
5831 eval.template values<1,true,false> (temp1, values_quad[c]);
5836 for (
unsigned int c=0; c<n_components; c++)
5841 if (evaluate_grad ==
true)
5844 eval.template gradients<0,true,false> (values_dofs[c], temp1);
5845 eval.template values<1,true,false> (temp1, temp2);
5846 eval.template values<2,true,false> (temp2, gradients_quad[c][0]);
5849 if (evaluate_lapl ==
true)
5852 if (evaluate_grad ==
false)
5854 eval.template gradients<0,true,false> (values_dofs[c], temp1);
5855 eval.template values<1,true,false> (temp1, temp2);
5857 eval.template gradients<2,true,false> (temp2, hessians_quad[c][d4]);
5860 eval.template gradients<1,true,false> (temp1, temp2);
5861 eval.template values<2,true,false> (temp2, hessians_quad[c][d3]);
5864 eval.template hessians<0,true,false>(values_dofs[c], temp1);
5865 eval.template values<1,true,false> (temp1, temp2);
5866 eval.template values<2,true,false> (temp2, hessians_quad[c][0]);
5870 eval.template values<0,true,false> (values_dofs[c], temp1);
5871 if (evaluate_grad ==
true)
5873 eval.template gradients<1,true,false>(temp1, temp2);
5874 eval.template values<2,true,false> (temp2, gradients_quad[c][d1]);
5877 if (evaluate_lapl ==
true)
5880 if (evaluate_grad ==
false)
5881 eval.template gradients<1,true,false>(temp1, temp2);
5882 eval.template gradients<2,true,false> (temp2, hessians_quad[c][d5]);
5885 eval.template hessians<1,true,false> (temp1, temp2);
5886 eval.template values<2,true,false> (temp2, hessians_quad[c][d1]);
5890 eval.template values<1,true,false> (temp1, temp2);
5891 if (evaluate_grad ==
true)
5892 eval.template gradients<2,true,false> (temp2, gradients_quad[c][d2]);
5896 if (evaluate_lapl ==
true)
5897 eval.template hessians<2,true,false>(temp2, hessians_quad[c][d2]);
5900 if (evaluate_val ==
true)
5901 eval.template values<2,true,false> (temp2, values_quad[c]);
5911 if (type == MatrixFreeFunctions::tensor_symmetric_plus_dg0 && evaluate_val)
5912 for (
unsigned int c=0; c<n_components; ++c)
5913 for (
unsigned int q=0; q<Eval::n_q_points; ++q)
5914 values_quad[c][q] += values_dofs[c][Eval::dofs_per_cell];
5919 template <MatrixFreeFunctions::ElementType type,
int dim,
int fe_degree,
5920 int n_q_points_1d,
int n_components,
typename Number>
5923 FEEvaluationImpl<type,dim,fe_degree,n_q_points_1d,n_components,Number>
5924 ::integrate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
5928 const bool integrate_val,
5929 const bool integrate_grad)
5931 const EvaluatorVariant variant =
5932 EvaluatorSelector<type,(fe_degree+n_q_points_1d>4)>::variant;
5933 typedef EvaluatorTensorProduct<variant, dim, fe_degree, n_q_points_1d,
5935 Eval eval (variant == evaluate_evenodd ? shape_info.shape_val_evenodd :
5936 shape_info.shape_values,
5937 variant == evaluate_evenodd ? shape_info.shape_gra_evenodd :
5938 shape_info.shape_gradients,
5939 variant == evaluate_evenodd ? shape_info.shape_hes_evenodd :
5940 shape_info.shape_hessians);
5942 const unsigned int temp_size = Eval::dofs_per_cell > Eval::n_q_points ?
5943 Eval::dofs_per_cell : Eval::n_q_points;
5950 n_components*Eval::dofs_per_cell];
5952 if (type == MatrixFreeFunctions::truncated_tensor)
5954 for (
unsigned int c=0; c<n_components; ++c)
5955 expanded_dof_values[c] = &data_array[c*Eval::dofs_per_cell];
5956 values_dofs = expanded_dof_values;
5962 const unsigned int d1 = dim>1?1:0;
5963 const unsigned int d2 = dim>2?2:0;
5968 for (
unsigned int c=0; c<n_components; c++)
5970 if (integrate_val ==
true)
5971 eval.template values<0,false,false> (values_quad[c], values_dofs[c]);
5972 if (integrate_grad ==
true)
5974 if (integrate_val ==
true)
5975 eval.template gradients<0,false,true> (gradients_quad[c][0], values_dofs[c]);
5977 eval.template gradients<0,false,false> (gradients_quad[c][0], values_dofs[c]);
5983 for (
unsigned int c=0; c<n_components; c++)
5985 if (integrate_val ==
true)
5988 eval.template values<0,false,false> (values_quad[c], temp1);
5990 if (integrate_grad ==
true)
5991 eval.template gradients<0,false,true> (gradients_quad[c][0], temp1);
5992 eval.template values<1,false,false>(temp1, values_dofs[c]);
5994 if (integrate_grad ==
true)
5997 eval.template values<0,false,false> (gradients_quad[c][d1], temp1);
5998 if (integrate_val ==
false)
6000 eval.template gradients<1,false,false>(temp1, values_dofs[c]);
6002 eval.template gradients<0,false,false> (gradients_quad[c][0], temp1);
6003 eval.template values<1,false,true> (temp1, values_dofs[c]);
6006 eval.template gradients<1,false,true>(temp1, values_dofs[c]);
6012 for (
unsigned int c=0; c<n_components; c++)
6014 if (integrate_val ==
true)
6017 eval.template values<0,false,false> (values_quad[c], temp1);
6019 if (integrate_grad ==
true)
6020 eval.template gradients<0,false,true> (gradients_quad[c][0], temp1);
6021 eval.template values<1,false,false>(temp1, temp2);
6022 if (integrate_grad ==
true)
6024 eval.template values<0,false,false> (gradients_quad[c][d1], temp1);
6025 eval.template gradients<1,false,true>(temp1, temp2);
6027 eval.template values<2,false,false> (temp2, values_dofs[c]);
6029 else if (integrate_grad ==
true)
6031 eval.template gradients<0,false,false>(gradients_quad[c][0], temp1);
6032 eval.template values<1,false,false> (temp1, temp2);
6033 eval.template values<0,false,false> (gradients_quad[c][d1], temp1);
6034 eval.template gradients<1,false,true>(temp1, temp2);
6035 eval.template values<2,false,false> (temp2, values_dofs[c]);
6037 if (integrate_grad ==
true)
6040 eval.template values<0,false,false> (gradients_quad[c][d2], temp1);
6041 eval.template values<1,false,false> (temp1, temp2);
6042 eval.template gradients<2,false,true> (temp2, values_dofs[c]);
6052 if (type == MatrixFreeFunctions::tensor_symmetric_plus_dg0)
6055 for (
unsigned int c=0; c<n_components; ++c)
6057 values_dofs[c][Eval::dofs_per_cell] = values_quad[c][0];
6058 for (
unsigned int q=1; q<Eval::n_q_points; ++q)
6059 values_dofs[c][Eval::dofs_per_cell] += values_quad[c][q];
6062 for (
unsigned int c=0; c<n_components; ++c)
6066 if (type == MatrixFreeFunctions::truncated_tensor)
6068 unsigned int count_p = 0, count_q = 0;
6069 for (
unsigned int i=0; i<(dim>2?fe_degree+1:1); ++i)
6071 for (
unsigned int j=0; j<(dim>1?fe_degree+1-i:1); ++j)
6073 for (
unsigned int k=0; k<fe_degree+1-j-i; ++k, ++count_p, ++count_q)
6075 for (
unsigned int c=0; c<n_components; ++c)
6076 values_dofs_actual[c][count_p] = expanded_dof_values[c][count_q];
6080 count_q += i*(fe_degree+1);
6088 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename Number>
6089 struct FEEvaluationImpl<MatrixFreeFunctions::tensor_gausslobatto, dim,
6090 fe_degree, n_q_points_1d, n_components, Number>
6093 void evaluate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
6098 const bool evaluate_val,
6099 const bool evaluate_grad,
6100 const bool evaluate_lapl);
6103 void integrate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
6107 const bool integrate_val,
6108 const bool integrate_grad);
6111 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename Number>
6114 FEEvaluationImpl<MatrixFreeFunctions::tensor_gausslobatto, dim,
6115 fe_degree, n_q_points_1d, n_components, Number>
6116 ::evaluate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
6121 const bool evaluate_val,
6122 const bool evaluate_grad,
6123 const bool evaluate_lapl)
6125 typedef EvaluatorTensorProduct<evaluate_evenodd, dim, fe_degree, fe_degree+1,
6127 Eval eval (shape_info.shape_val_evenodd, shape_info.shape_gra_evenodd,
6128 shape_info.shape_hes_evenodd);
6133 const unsigned int d1 = dim>1?1:0;
6134 const unsigned int d2 = dim>2?2:0;
6135 const unsigned int d3 = dim>2?3:0;
6136 const unsigned int d4 = dim>2?4:0;
6137 const unsigned int d5 = dim>2?5:0;
6142 if (evaluate_val ==
true)
6143 std::memcpy (values_quad[0], values_dofs[0],
6144 eval.dofs_per_cell * n_components *
6145 sizeof (values_dofs[0][0]));
6146 for (
unsigned int c=0; c<n_components; c++)
6148 if (evaluate_grad ==
true)
6149 eval.template gradients<0,true,false>(values_dofs[c], gradients_quad[c][0]);
6150 if (evaluate_lapl ==
true)
6151 eval.template hessians<0,true,false> (values_dofs[c], hessians_quad[c][0]);
6156 if (evaluate_val ==
true)
6158 std::memcpy (values_quad[0], values_dofs[0],
6159 Eval::dofs_per_cell * n_components *
6160 sizeof (values_dofs[0][0]));
6162 if (evaluate_grad ==
true)
6163 for (
unsigned int comp=0; comp<n_components; comp++)
6166 eval.template gradients<0,true,false> (values_dofs[comp],
6167 gradients_quad[comp][0]);
6169 eval.template gradients<1,true,false> (values_dofs[comp],
6170 gradients_quad[comp][d1]);
6172 if (evaluate_lapl ==
true)
6173 for (
unsigned int comp=0; comp<n_components; comp++)
6176 eval.template hessians<0,true,false> (values_dofs[comp],
6177 hessians_quad[comp][0]);
6179 eval.template hessians<1,true,false> (values_dofs[comp],
6180 hessians_quad[comp][d1]);
6184 eval.template gradients<0,true,false> (values_dofs[comp], temp1);
6185 eval.template gradients<1,true,false> (temp1, hessians_quad[comp][d1+d1]);
6190 if (evaluate_val ==
true)
6192 std::memcpy (values_quad[0], values_dofs[0],
6193 Eval::dofs_per_cell * n_components *
6194 sizeof (values_dofs[0][0]));
6196 if (evaluate_grad ==
true)
6197 for (
unsigned int comp=0; comp<n_components; comp++)
6200 eval.template gradients<0,true,false> (values_dofs[comp],
6201 gradients_quad[comp][0]);
6203 eval.template gradients<1,true,false> (values_dofs[comp],
6204 gradients_quad[comp][d1]);
6206 eval.template gradients<2,true,false> (values_dofs[comp],
6207 gradients_quad[comp][d2]);
6209 if (evaluate_lapl ==
true)
6210 for (
unsigned int comp=0; comp<n_components; comp++)
6213 eval.template hessians<0,true,false> (values_dofs[comp],
6214 hessians_quad[comp][0]);
6216 eval.template hessians<1,true,false> (values_dofs[comp],
6217 hessians_quad[comp][d1]);
6219 eval.template hessians<2,true,false> (values_dofs[comp],
6220 hessians_quad[comp][d2]);
6224 eval.template gradients<0,true,false> (values_dofs[comp], temp1);
6225 eval.template gradients<1,true,false> (temp1, hessians_quad[comp][d3]);
6227 eval.template gradients<2,true,false> (temp1, hessians_quad[comp][d4]);
6229 eval.template gradients<1,true,false> (values_dofs[comp], temp1);
6230 eval.template gradients<2,true,false> (temp1, hessians_quad[comp][d5]);
6238 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename Number>
6241 FEEvaluationImpl<MatrixFreeFunctions::tensor_gausslobatto, dim,
6242 fe_degree, n_q_points_1d, n_components, Number>
6243 ::integrate (
const MatrixFreeFunctions::ShapeInfo<Number> &shape_info,
6247 const bool integrate_val,
6248 const bool integrate_grad)
6250 typedef EvaluatorTensorProduct<evaluate_evenodd, dim, fe_degree, fe_degree+1,
6252 Eval eval (shape_info.shape_val_evenodd, shape_info.shape_gra_evenodd,
6253 shape_info.shape_hes_evenodd);
6258 const unsigned int d1 = dim>1?1:0;
6259 const unsigned int d2 = dim>2?2:0;
6261 if (integrate_val ==
true)
6262 std::memcpy (values_dofs[0], values_quad[0],
6263 Eval::dofs_per_cell * n_components *
6264 sizeof (values_dofs[0][0]));
6268 for (
unsigned int c=0; c<n_components; c++)
6270 if (integrate_grad ==
true)
6272 if (integrate_val ==
true)
6273 eval.template gradients<0,false,true> (gradients_quad[c][0],
6276 eval.template gradients<0,false,false> (gradients_quad[c][0],
6283 if (integrate_grad ==
true)
6284 for (
unsigned int comp=0; comp<n_components; comp++)
6288 if (integrate_val ==
true)
6289 eval.template gradients<0, false, true> (gradients_quad[comp][0],
6292 eval.template gradients<0, false, false> (gradients_quad[comp][0],
6296 eval.template gradients<1, false, true> (gradients_quad[comp][d1],
6302 if (integrate_grad ==
true)
6303 for (
unsigned int comp=0; comp<n_components; comp++)
6307 if (integrate_val ==
true)
6308 eval.template gradients<0, false, true> (gradients_quad[comp][0],
6311 eval.template gradients<0, false, false> (gradients_quad[comp][0],
6315 eval.template gradients<1, false, true> (gradients_quad[comp][d1],
6319 eval.template gradients<2, false, true> (gradients_quad[comp][d2],
6336 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6341 const unsigned int fe_no,
6342 const unsigned int quad_no)
6344 BaseClass (data_in, fe_no, quad_no, fe_degree, n_q_points),
6347 check_template_arguments(fe_no);
6348 set_data_pointers();
6353 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6361 const unsigned int first_selected_component)
6363 BaseClass (mapping, fe, quadrature, update_flags,
6364 first_selected_component,
6369 set_data_pointers();
6374 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6381 const unsigned int first_selected_component)
6384 first_selected_component,
6389 set_data_pointers();
6394 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6396 template <
int n_components_other>
6401 const unsigned int first_selected_component)
6406 first_selected_component, &other),
6410 set_data_pointers();
6415 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6424 set_data_pointers();
6429 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6437 const unsigned int desired_dofs_per_cell = this->data->
dofs_per_cell;
6440 for (
unsigned int c=0; c<n_components_; ++c)
6442 this->values_dofs[c] = &my_data_array[c*desired_dofs_per_cell];
6443 this->values_quad[c] = &my_data_array[n_components*desired_dofs_per_cell+c*n_q_points];
6444 for (
unsigned int d=0; d<dim; ++d)
6445 this->gradients_quad[c][d] = &my_data_array[n_components*(desired_dofs_per_cell+
6448 (c*dim+d)*n_q_points];
6449 for (
unsigned int d=0; d<(dim*dim+dim)/2; ++d)
6450 this->hessians_quad[c][d] = &my_data_array[n_components*((dim+1)*n_q_points+
6451 desired_dofs_per_cell)
6453 (c*(dim*dim+dim)+d)*n_q_points];
6458 case internal::MatrixFreeFunctions::tensor_symmetric:
6460 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
6461 dim, fe_degree, n_q_points_1d, n_components_,
6464 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
6465 dim, fe_degree, n_q_points_1d, n_components_,
6469 case internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0:
6471 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
6472 dim, fe_degree, n_q_points_1d, n_components_,
6475 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
6476 dim, fe_degree, n_q_points_1d, n_components_,
6480 case internal::MatrixFreeFunctions::tensor_general:
6482 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
6483 dim, fe_degree, n_q_points_1d, n_components_,
6486 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
6487 dim, fe_degree, n_q_points_1d, n_components_,
6491 case internal::MatrixFreeFunctions::tensor_gausslobatto:
6493 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
6494 dim, fe_degree, n_q_points_1d, n_components_,
6497 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
6498 dim, fe_degree, n_q_points_1d, n_components_,
6502 case internal::MatrixFreeFunctions::truncated_tensor:
6504 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
6505 dim, fe_degree, n_q_points_1d, n_components_,
6508 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
6509 dim, fe_degree, n_q_points_1d, n_components_,
6521 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6533 n_q_points != this->data->n_q_points)
6535 std::string message =
6536 "-------------------------------------------------------\n";
6537 message +=
"Illegal arguments in constructor/wrong template arguments!\n";
6538 message +=
" Called --> FEEvaluation<dim,";
6542 message +=
",Number>(data";
6557 proposed_dof_comp = fe_no;
6559 for (
unsigned int no=0; no<this->matrix_info->
n_components(); ++no)
6560 if (this->matrix_info->
get_shape_info(no,0,this->active_fe_index,0).fe_degree
6563 proposed_dof_comp = no;
6567 this->mapping_info->
mapping_data_gen[this->quad_no].n_q_points[this->active_quad_index])
6568 proposed_quad_comp = this->quad_no;
6570 for (
unsigned int no=0; no<this->mapping_info->
mapping_data_gen.size(); ++no)
6571 if (this->mapping_info->
mapping_data_gen[no].n_q_points[this->active_quad_index]
6574 proposed_quad_comp = no;
6581 if (proposed_dof_comp != fe_no)
6582 message +=
"Wrong vector component selection:\n";
6584 message +=
"Wrong quadrature formula selection:\n";
6585 message +=
" Did you mean FEEvaluation<dim,";
6589 message +=
",Number>(data";
6596 std::string correct_pos;
6597 if (proposed_dof_comp != fe_no)
6598 correct_pos =
" ^ ";
6601 if (proposed_quad_comp != this->quad_no)
6602 correct_pos +=
" ^\n";
6604 correct_pos +=
" \n";
6605 message +=
" " + correct_pos;
6609 const unsigned int proposed_n_q_points_1d =
static_cast<unsigned int>(std::pow(1.001*this->data->
n_q_points,1./dim));
6610 message +=
"Wrong template arguments:\n";
6611 message +=
" Did you mean FEEvaluation<dim,";
6615 message +=
",Number>(data";
6622 std::string correct_pos;
6627 if (proposed_n_q_points_1d != n_q_points_1d)
6628 correct_pos +=
" ^\n";
6630 correct_pos +=
" \n";
6631 message +=
" " + correct_pos;
6634 n_q_points == this->data->n_q_points,
6641 n_q_points[this->active_quad_index]);
6643 this->dof_info->dofs_per_cell[this->active_fe_index]);
6650 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6658 ExcNotInitialized());
6664 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
6670 return this->quadrature_points[q];
6672 point[0] = this->quadrature_points[q%n_q_points_1d][0];
6673 point[1] = this->quadrature_points[q/n_q_points_1d][1];
6676 point[0] = this->quadrature_points[q%n_q_points_1d][0];
6677 point[1] = this->quadrature_points[(q/n_q_points_1d)%n_q_points_1d][1];
6678 point[2] = this->quadrature_points[q/(n_q_points_1d*n_q_points_1d)][2];
6681 Assert (
false, ExcNotImplemented());
6687 return this->quadrature_points[q];
6692 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6698 const bool evaluate_grad,
6699 const bool evaluate_lapl)
6701 Assert (this->dof_values_initialized ==
true,
6702 internal::ExcAccessToUninitializedField());
6703 Assert(this->matrix_info != 0 ||
6704 this->mapped_geometry->is_initialized(), ExcNotInitialized());
6708 evaluate_funct (*this->data, &this->values_dofs[0],
6709 this->values_quad, this->gradients_quad, this->hessians_quad,
6710 evaluate_val, evaluate_grad, evaluate_lapl);
6713 if (evaluate_val ==
true)
6714 this->values_quad_initialized =
true;
6715 if (evaluate_grad ==
true)
6716 this->gradients_quad_initialized =
true;
6717 if (evaluate_lapl ==
true)
6718 this->hessians_quad_initialized =
true;
6724 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
6729 ::integrate (
bool integrate_val,
bool integrate_grad)
6731 if (integrate_val ==
true)
6732 Assert (this->values_quad_submitted ==
true,
6733 internal::ExcAccessToUninitializedField());
6734 if (integrate_grad ==
true)
6735 Assert (this->gradients_quad_submitted ==
true,
6736 internal::ExcAccessToUninitializedField());
6737 Assert(this->matrix_info != 0 ||
6738 this->mapped_geometry->is_initialized(), ExcNotInitialized());
6742 integrate_funct (*this->data, this->values_dofs, this->values_quad,
6743 this->gradients_quad, integrate_val, integrate_grad);
6746 this->dof_values_initialized =
true;
6752 #endif // ifndef DOXYGEN 6755 DEAL_II_NAMESPACE_CLOSE
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 distribute_local_to_global(VectorType &dst) const
static const unsigned int invalid_unsigned_int
const Number * constraint_pool_begin(const unsigned int pool_index) const
#define AssertDimension(dim1, dim2)
AlignedVector< std::pair< Tensor< 1, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > cartesian_data
unsigned int get_cell_data_number() const
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info(const unsigned int fe_component=0) const
void integrate(const bool integrate_val, const bool integrate_grad)
Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArray< Number > > > get_hessian(const unsigned int q_point) const
void set_dof_values(VectorType &dst) const
const MatrixFree< dim, Number > * matrix_info
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
unsigned int cell_data_number
void read_write_operation(const VectorOperation &operation, VectorType *vectors[]) const
const unsigned int first_selected_component
value_type get_laplacian(const unsigned int q_point) const
bool quadrature_points_initialized
const unsigned int dofs_per_cell
unsigned int n_components() const
const unsigned int n_fe_components
bool gradients_quad_initialized
::ExceptionBase & ExcMessage(std::string arg1)
#define AssertIndexRange(index, range)
const VectorizedArray< Number > * begin_gradients() const
std::vector< types::global_dof_index > old_style_dof_indices
unsigned int get_cell_data_index(const unsigned int cell_chunk_no) const
const VectorizedArray< Number > * begin_hessians() const
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian
friend class FEEvaluationBase
std::vector< MappingInfoDependent > mapping_data_gen
#define AssertThrow(cond, exc)
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)
const VectorizedArray< Number > * begin_dof_values() const
void evaluate(const bool evaluate_val, const bool evaluate_grad, const bool evaluate_hess=false)
void read_dof_values(const VectorType &src)
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
Number local_element(const size_type local_index) const
void fill_JxW_values(AlignedVector< VectorizedArray< Number > > &JxW_values) const
value_type get_dof_value(const unsigned int dof) const
bool values_quad_submitted
const unsigned int active_fe_index
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
unsigned int global_dof_index
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian_grad
#define Assert(cond, exc)
unsigned int element_multiplicity(const unsigned int index) const
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)
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
bool mapping_initialized() const
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
bool JxW_values_initialized
const VectorizedArray< Number > * quadrature_weights
gradient_type get_hessian_diagonal(const unsigned int q_point) const
bool hessians_quad_initialized
const internal::MatrixFreeFunctions::SizeInfo & get_size_info() const
value_type integrate_value() const
void submit_value(const value_type val_in, const unsigned int q_point)
internal::MatrixFreeFunctions::CellType cell_type
std_cxx11::shared_ptr< internal::MatrixFreeFunctions::ShapeInfo< Number > > stored_shape_info
void read_dof_values_plain(const VectorType &src)
gradient_type get_gradient(const unsigned int q_point) const
const VectorizedArray< Number > * begin_values() const
bool partitioners_are_compatible(const Utilities::MPI::Partitioner &part) const
unsigned int dofs_per_cell
bool second_derivatives_initialized
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info() const
bool indices_initialized() const
void check_template_arguments(const unsigned int fe_no)
std::vector< unsigned int > lexicographic_numbering
bool values_quad_initialized
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
internal::MatrixFreeFunctions::CellType get_cell_type() const
Point< dim, VectorizedArray< Number > > quadrature_point(const unsigned int q_point) const
bool gradients_quad_submitted
std_cxx1x::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
bool dof_values_initialized