17 #ifndef dealii__matrix_free_fe_evaluation_h 18 #define dealii__matrix_free_fe_evaluation_h 21 #include <deal.II/base/array_view.h> 22 #include <deal.II/base/config.h> 23 #include <deal.II/base/exceptions.h> 24 #include <deal.II/base/smartpointer.h> 25 #include <deal.II/base/symmetric_tensor.h> 26 #include <deal.II/base/template_constraints.h> 27 #include <deal.II/base/vectorization.h> 28 #include <deal.II/matrix_free/mapping_data_on_the_fly.h> 29 #include <deal.II/matrix_free/matrix_free.h> 30 #include <deal.II/matrix_free/shape_info.h> 31 #include <deal.II/matrix_free/evaluation_kernels.h> 32 #include <deal.II/matrix_free/tensor_product_kernels.h> 35 DEAL_II_NAMESPACE_OPEN
44 template <
typename>
class Vector;
52 template <
int dim,
int fe_degree,
int n_q_points_1d = fe_degree+1,
53 int n_components_ = 1,
typename Number =
double >
class FEEvaluation;
79 template <
int dim,
int n_components_,
typename Number>
83 typedef Number number_type;
86 static const unsigned int dimension = dim;
87 static const unsigned int n_components = n_components_;
120 template <
typename DoFHandlerType,
bool level_dof_access>
151 internal::MatrixFreeFunctions::CellType
get_cell_type()
const;
203 template <
typename VectorType>
205 const unsigned int first_index = 0);
235 template <
typename VectorType>
237 const unsigned int first_index = 0);
262 template<
typename VectorType>
264 const unsigned int first_index = 0)
const;
292 template<
typename VectorType>
294 const unsigned int first_index = 0)
const;
328 const unsigned int dof);
355 const unsigned int q_point);
381 const unsigned int q_point);
580 const std::vector<unsigned int> &
604 const unsigned int fe_no,
606 const unsigned int fe_degree,
607 const unsigned int n_q_points);
643 template <
int n_components_other>
673 template<
typename VectorType,
typename VectorOperation>
675 VectorType *vectors[])
const;
684 template<
typename VectorType>
787 const internal::MatrixFreeFunctions::DoFInfo *
dof_info;
914 std_cxx1x::shared_ptr<internal::MatrixFreeFunctions::MappingDataOnTheFly<dim,Number> >
mapped_geometry;
945 template <
int,
int,
int,
int,
typename>
friend class FEEvaluation;
957 template <
int dim,
int n_components_,
typename Number>
961 typedef Number number_type;
964 static const unsigned int dimension = dim;
965 static const unsigned int n_components = n_components_;
977 const unsigned int fe_no,
979 const unsigned int fe_degree,
980 const unsigned int n_q_points);
986 template <
int n_components_other>
1016 template <
int dim,
typename Number>
1020 typedef Number number_type;
1023 static const unsigned int dimension = dim;
1042 const unsigned int dof);
1061 const unsigned int q_point);
1079 const unsigned int q_point);
1121 const unsigned int fe_no,
1123 const unsigned int fe_degree,
1124 const unsigned int n_q_points);
1130 template <
int n_components_other>
1160 template <
int dim,
typename Number>
1164 typedef Number number_type;
1167 static const unsigned int dimension = dim;
1168 static const unsigned int n_components = dim;
1197 get_curl (
const unsigned int q_point)
const;
1223 const unsigned int q_point);
1234 const unsigned int q_point);
1245 const unsigned int q_point);
1256 const unsigned int q_point);
1263 const unsigned int q_point);
1274 const unsigned int fe_no,
1276 const unsigned int dofs_per_cell,
1277 const unsigned int n_q_points);
1283 template <
int n_components_other>
1312 template <
typename Number>
1316 typedef Number number_type;
1319 static const unsigned int dimension = 1;
1338 const unsigned int dof);
1357 const unsigned int q_point);
1375 const unsigned int q_point);
1417 const unsigned int fe_no,
1419 const unsigned int fe_degree,
1420 const unsigned int n_q_points);
1426 template <
int n_components_other>
1944 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
1950 typedef Number number_type;
1951 typedef typename BaseClass::value_type value_type;
1952 typedef typename BaseClass::gradient_type gradient_type;
1953 static const unsigned int dimension = dim;
1954 static const unsigned int n_components = n_components_;
1965 const unsigned int fe_no = 0,
1966 const unsigned int quad_no = 0);
2020 template <
int n_components_other>
2049 void evaluate (
const bool evaluate_values,
2050 const bool evaluate_gradients,
2051 const bool evaluate_hessians =
false);
2060 void integrate (
const bool integrate_values,
2061 const bool integrate_gradients);
2100 const bool evaluate_val,
2101 const bool evaluate_grad,
2102 const bool evaluate_lapl);
2112 const bool evaluate_val,
2113 const bool evaluate_grad);
2120 namespace MatrixFreeFunctions
2124 template <
int dim,
int degree>
2125 struct DGP_dofs_per_cell
2128 static const unsigned int value =
2129 (DGP_dofs_per_cell<dim-1,degree>::value * (degree+dim)) / dim;
2133 template <
int degree>
2134 struct DGP_dofs_per_cell<1,degree>
2136 static const unsigned int value = degree+1;
2150 template <
int dim,
int n_components_,
typename Number>
2154 const unsigned int fe_no_in,
2155 const unsigned int quad_no_in,
2156 const unsigned int fe_degree,
2157 const unsigned int n_q_points)
2159 scratch_data_array (data_in.acquire_scratch_data()),
2160 quad_no (quad_no_in),
2161 n_fe_components (data_in.get_dof_info(fe_no_in).
n_components),
2162 active_fe_index (fe_degree !=
numbers::invalid_unsigned_int ?
2163 data_in.get_dof_info(fe_no_in).fe_index_from_degree(fe_degree)
2166 active_quad_index (fe_degree !=
numbers::invalid_unsigned_int ?
2167 data_in.get_mapping_info().
2168 mapping_data_gen[quad_no_in].
2169 quad_index_from_n_q_points(n_q_points)
2172 matrix_info (&data_in),
2173 dof_info (&data_in.get_dof_info(fe_no_in)),
2174 mapping_info (&data_in.get_mapping_info()),
2175 data (&data_in.get_shape_info
2176 (fe_no_in, quad_no_in, active_fe_index,
2177 active_quad_index)),
2181 quadrature_weights (mapping_info->mapping_data_gen[quad_no].
2182 quadrature_weights[active_quad_index].begin()),
2183 quadrature_points (0),
2185 jacobian_grad_upper(0),
2186 cell (
numbers::invalid_unsigned_int),
2187 cell_type (
internal::MatrixFreeFunctions::undefined),
2188 cell_data_number (
numbers::invalid_unsigned_int),
2189 dof_values_initialized (false),
2190 values_quad_initialized (false),
2191 gradients_quad_initialized(false),
2192 hessians_quad_initialized (false),
2193 values_quad_submitted (false),
2194 gradients_quad_submitted (false),
2195 first_selected_component (0)
2197 set_data_pointers();
2203 dof_info->dofs_per_cell[active_fe_index]/n_fe_components);
2206 Assert (n_fe_components == 1 ||
2209 ExcMessage (
"The underlying FE is vector-valued. In this case, the " 2210 "template argument n_components must be a the same " 2211 "as the number of underlying vector components."));
2220 template <
int dim,
int n_components_,
typename Number>
2221 template <
int n_components_other>
2228 const unsigned int first_selected_component,
2232 quad_no (
numbers::invalid_unsigned_int),
2233 n_fe_components (n_components_),
2234 active_fe_index (
numbers::invalid_unsigned_int),
2235 active_quad_index (
numbers::invalid_unsigned_int),
2240 data (new
internal::MatrixFreeFunctions::ShapeInfo<Number>(quadrature, fe, fe.component_to_base_index(first_selected_component).first)),
2244 quadrature_weights (0),
2245 quadrature_points (0),
2247 jacobian_grad_upper(0),
2249 cell_type (
internal::MatrixFreeFunctions::general),
2250 cell_data_number (
numbers::invalid_unsigned_int),
2251 dof_values_initialized (false),
2252 values_quad_initialized (false),
2253 gradients_quad_initialized(false),
2254 hessians_quad_initialized (false),
2255 values_quad_submitted (false),
2256 gradients_quad_submitted (false),
2259 first_selected_component (fe.component_to_base_index(first_selected_component).second)
2261 const unsigned int base_element_number =
2263 set_data_pointers();
2271 MappingDataOnTheFly<dim,Number>(mapping, quadrature,
2273 jacobian = mapped_geometry->get_inverse_jacobians().begin();
2274 J_value = mapped_geometry->get_JxW_values().begin();
2275 quadrature_points = mapped_geometry->get_quadrature_points().begin();
2279 ExcMessage(
"The underlying element must at least contain as many " 2280 "components as requested by this class"));
2281 (void) base_element_number;
2286 template <
int dim,
int n_components_,
typename Number>
2291 scratch_data_array (other.matrix_info == 0 ?
2293 other.matrix_info->acquire_scratch_data()),
2294 quad_no (other.quad_no),
2295 n_fe_components (other.n_fe_components),
2296 active_fe_index (other.active_fe_index),
2297 active_quad_index (other.active_quad_index),
2298 matrix_info (other.matrix_info),
2299 dof_info (other.dof_info),
2300 mapping_info (other.mapping_info),
2301 data (other.matrix_info == 0 ?
2302 new
internal::MatrixFreeFunctions::ShapeInfo<Number>(*other.data) :
2307 quadrature_weights (mapping_info != 0 ?
2308 mapping_info->mapping_data_gen[quad_no].
2309 quadrature_weights[active_quad_index].begin()
2312 quadrature_points (0),
2314 jacobian_grad_upper(0),
2315 cell (
numbers::invalid_unsigned_int),
2316 cell_type (
internal::MatrixFreeFunctions::general),
2317 cell_data_number (
numbers::invalid_unsigned_int),
2318 dof_values_initialized (false),
2319 values_quad_initialized (false),
2320 gradients_quad_initialized(false),
2321 hessians_quad_initialized (false),
2322 values_quad_submitted (false),
2323 gradients_quad_submitted (false),
2324 first_selected_component (other.first_selected_component)
2326 set_data_pointers();
2331 mapped_geometry.reset
2333 MappingDataOnTheFly<dim,Number>(other.
mapped_geometry->get_fe_values().get_mapping(),
2336 jacobian = mapped_geometry->get_inverse_jacobians().begin();
2337 J_value = mapped_geometry->get_JxW_values().begin();
2338 quadrature_points = mapped_geometry->get_quadrature_points().begin();
2345 template <
int dim,
int n_components_,
typename Number>
2358 if (matrix_info == 0)
2361 delete scratch_data_array;
2381 set_data_pointers();
2386 quadrature_weights = mapping_info != 0 ?
2388 quadrature_weights[active_quad_index].begin()
2391 quadrature_points = 0;
2393 jacobian_grad_upper = 0;
2395 cell_type = internal::MatrixFreeFunctions::general;
2401 mapped_geometry.reset
2403 MappingDataOnTheFly<dim,Number>(other.
mapped_geometry->get_fe_values().get_mapping(),
2406 jacobian = mapped_geometry->get_inverse_jacobians().begin();
2407 J_value = mapped_geometry->get_JxW_values().begin();
2408 quadrature_points = mapped_geometry->get_quadrature_points().begin();
2417 template <
int dim,
int n_components_,
typename Number>
2421 if (matrix_info != 0)
2432 delete scratch_data_array;
2439 template <
int dim,
int n_components_,
typename Number>
2447 const unsigned int tensor_dofs_per_cell =
2448 Utilities::fixed_power<dim>(this->data->
fe_degree+1);
2449 const unsigned int dofs_per_cell = this->data->
dofs_per_cell;
2450 const unsigned int n_quadrature_points = this->data->
n_q_points;
2452 const unsigned int shift = std::max(tensor_dofs_per_cell+1, dofs_per_cell)*
2453 (n_components_+2) + 2 * n_quadrature_points;
2454 const unsigned int allocated_size = shift + n_components_ * dofs_per_cell
2455 + (n_components_*(dim*dim+2*dim+1)*n_quadrature_points);
2456 scratch_data_array->resize_fast(allocated_size);
2459 for (
unsigned int c=0; c<n_components_; ++c)
2461 this->values_dofs[c] = scratch_data_array->begin() + c*dofs_per_cell;
2462 this->values_quad[c] = scratch_data_array->begin() +
2464 for (
unsigned int d=0;
d<dim; ++
d)
2465 this->gradients_quad[c][d] = scratch_data_array->begin() +
2467 (c*dim+d)*n_quadrature_points;
2468 for (
unsigned int d=0;
d<(dim*dim+dim)/2; ++
d)
2469 this->hessians_quad[c][d] = scratch_data_array->begin() +
2470 n_components*((dim+1)*n_quadrature_points + dofs_per_cell) +
2471 (c*(dim*dim+dim)+d)*n_quadrature_points;
2473 scratch_data = scratch_data_array->begin() + n_components_ * dofs_per_cell
2474 + (n_components_*(dim*dim+2*dim+1)*n_quadrature_points);
2479 template <
int dim,
int n_components_,
typename Number>
2484 Assert (mapped_geometry == 0,
2485 ExcMessage(
"FEEvaluation was initialized without a matrix-free object." 2486 " Integer indexing is not possible"));
2487 if (mapped_geometry != 0)
2493 dof_info->cell_active_fe_index[cell_in] : 0),
2502 mapping_data_gen[quad_no].rowstart_q_points.size());
2504 rowstart_q_points[cell];
2506 quadrature_points.size());
2511 if (cell_type == internal::MatrixFreeFunctions::cartesian)
2513 cartesian_data = &mapping_info->
cartesian_data[cell_data_number].first;
2514 J_value = &mapping_info->
cartesian_data[cell_data_number].second;
2516 else if (cell_type == internal::MatrixFreeFunctions::affine)
2518 jacobian = &mapping_info->
affine_data[cell_data_number].first;
2519 J_value = &mapping_info->
affine_data[cell_data_number].second;
2523 const unsigned int rowstart = mapping_info->
2524 mapping_data_gen[quad_no].rowstart_jacobians[cell_data_number];
2526 mapping_data_gen[quad_no].jacobians.size());
2532 mapping_data_gen[quad_no].JxW_values.size());
2534 JxW_values[rowstart]);
2539 mapping_data_gen[quad_no].jacobians_grad_diag.size());
2541 jacobians_grad_diag[rowstart];
2543 mapping_data_gen[quad_no].jacobians_grad_upper.size());
2545 jacobians_grad_upper[rowstart];
2550 dof_values_initialized =
false;
2551 values_quad_initialized =
false;
2552 gradients_quad_initialized =
false;
2553 hessians_quad_initialized =
false;
2559 template <
int dim,
int n_components_,
typename Number>
2560 template <
typename DoFHandlerType,
bool level_dof_access>
2567 ExcMessage(
"Cannot use initialization from cell iterator if " 2568 "initialized from MatrixFree object. Use variant for " 2569 "on the fly computation with arguments as for FEValues " 2573 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2574 if (level_dof_access)
2575 cell->get_mg_dof_indices(local_dof_indices);
2577 cell->get_dof_indices(local_dof_indices);
2582 template <
int dim,
int n_components_,
typename Number>
2589 ExcMessage(
"Cannot use initialization from cell iterator if " 2590 "initialized from MatrixFree object. Use variant for " 2591 "on the fly computation with arguments as for FEValues " 2594 mapped_geometry->reinit(cell);
2599 template <
int dim,
int n_components_,
typename Number>
2606 return cell_data_number;
2611 template <
int dim,
int n_components_,
typename Number>
2613 internal::MatrixFreeFunctions::CellType
2622 template <
int dim,
int n_components_,
typename Number>
2633 template <
int dim,
int n_components_,
typename Number>
2641 if (this->cell_type == internal::MatrixFreeFunctions::cartesian ||
2642 this->cell_type == internal::MatrixFreeFunctions::affine)
2646 for (
unsigned int q=0; q<this->data->
n_q_points; ++q)
2647 JxW_values[q] = J * this->quadrature_weights[q];
2650 for (
unsigned int q=0; q<data->
n_q_points; ++q)
2651 JxW_values[q] = this->J_value[q];
2656 template <
int dim,
int n_components_,
typename Number>
2662 if (this->cell_type == internal::MatrixFreeFunctions::cartesian ||
2663 this->cell_type == internal::MatrixFreeFunctions::affine)
2666 return this->J_value[0] * this->quadrature_weights[q_point];
2669 return this->J_value[q_point];
2677 template <
typename VectorType>
2679 typename VectorType::value_type &
2680 vector_access (VectorType &vec,
2681 const unsigned int entry)
2689 template <
typename VectorType>
2691 typename VectorType::value_type
2692 vector_access (
const VectorType &vec,
2693 const unsigned int entry)
2703 template <
typename Number>
2707 const unsigned int entry)
2717 template <
typename Number>
2721 const unsigned int entry)
2731 template <
typename VectorType>
2733 void check_vector_compatibility (
const VectorType &vec,
2734 const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2740 dof_info.vector_partitioner->size());
2743 template <
typename Number>
2746 const internal::MatrixFreeFunctions::DoFInfo &dof_info)
2751 ExcMessage(
"The parallel layout of the given vector is not " 2752 "compatible with the parallel partitioning in MatrixFree. " 2753 "Use MatrixFree::initialize_dof_vector to get a " 2754 "compatible vector."));
2758 template <
typename Number>
2761 template <
typename VectorType>
2762 void process_dof (
const unsigned int index,
2766 res = vector_access (const_cast<const VectorType &>(vec), index);
2769 template <
typename VectorType>
2774 res =
const_cast<const VectorType &
>(vec)(index);
2777 void pre_constraints (
const Number &,
2783 template <
typename VectorType>
2784 void process_constraint (
const unsigned int index,
2785 const Number weight,
2789 res += weight * vector_access (const_cast<const VectorType &>(vec), index);
2792 void post_constraints (
const Number &sum,
2793 Number &write_pos)
const 2798 void process_empty (Number &res)
const 2805 template <
typename Number>
2806 struct VectorDistributorLocalToGlobal
2808 template <
typename VectorType>
2809 void process_dof (
const unsigned int index,
2813 vector_access (vec, index) += res;
2816 template <
typename VectorType>
2824 void pre_constraints (
const Number &input,
2830 template <
typename VectorType>
2831 void process_constraint (
const unsigned int index,
2832 const Number weight,
2836 vector_access (vec, index) += weight * res;
2839 void post_constraints (
const Number &,
2844 void process_empty (Number &)
const 2851 template <
typename Number>
2854 template <
typename VectorType>
2855 void process_dof (
const unsigned int index,
2859 vector_access (vec, index) = res;
2862 template <
typename VectorType>
2870 void pre_constraints (
const Number &,
2875 template <
typename VectorType>
2876 void process_constraint (
const unsigned int,
2883 void post_constraints (
const Number &,
2888 void process_empty (Number &)
const 2896 template <
typename VectorType,
bool>
2897 struct BlockVectorSelector {};
2899 template <
typename VectorType>
2900 struct BlockVectorSelector<VectorType,true>
2902 typedef typename VectorType::BlockType BaseVectorType;
2904 static BaseVectorType *get_vector_component (VectorType &vec,
2905 const unsigned int component)
2908 return &vec.block(component);
2912 template <
typename VectorType>
2913 struct BlockVectorSelector<VectorType,false>
2915 typedef VectorType BaseVectorType;
2917 static BaseVectorType *get_vector_component (VectorType &vec,
2924 template <
typename VectorType>
2925 struct BlockVectorSelector<
std::vector<VectorType>,false>
2927 typedef VectorType BaseVectorType;
2929 static BaseVectorType *get_vector_component (std::vector<VectorType> &vec,
2930 const unsigned int component)
2933 return &vec[component];
2937 template <
typename VectorType>
2938 struct BlockVectorSelector<
std::vector<VectorType *>,false>
2940 typedef VectorType BaseVectorType;
2942 static BaseVectorType *get_vector_component (std::vector<VectorType *> &vec,
2943 const unsigned int component)
2946 return vec[component];
2953 template <
int dim,
int n_components_,
typename Number>
2954 template<
typename VectorType,
typename VectorOperation>
2959 VectorType *src[])
const 2970 if (matrix_info == 0)
2974 unsigned int index = first_selected_component * this->data->
dofs_per_cell;
2977 for (
unsigned int i=0; i<this->data->
dofs_per_cell; ++i, ++index)
2980 *src[0], values_dofs[comp][i][0]);
2981 for (
unsigned int v=1; v<VectorizedArray<Number>::n_array_elements; ++v)
2982 operation.process_empty(values_dofs[comp][i][v]);
2996 const unsigned int *dof_indices = dof_info->begin_indices(cell);
2997 const std::pair<unsigned short,unsigned short> *indicators =
2998 dof_info->begin_indicators(cell);
2999 const std::pair<unsigned short,unsigned short> *indicators_end =
3000 dof_info->end_indicators(cell);
3001 unsigned int ind_local = 0;
3002 const unsigned int dofs_per_cell = this->data->
dofs_per_cell;
3004 const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
3005 const bool at_irregular_cell = n_irreg_components_filled > 0;
3009 if (n_fe_components == 1)
3011 const unsigned int n_local_dofs =
3014 internal::check_vector_compatibility (*src[comp], *dof_info);
3018 const_cast<Number *>(&values_dofs[comp][0][0]);
3022 if (at_irregular_cell ==
false)
3025 if (indicators != indicators_end)
3027 for ( ; indicators != indicators_end; ++indicators)
3030 for (
unsigned int j=0; j<indicators->first; ++j)
3032 operation.process_dof (dof_indices[j], *src[comp],
3033 local_data[comp][ind_local+j]);
3035 ind_local += indicators->first;
3036 dof_indices += indicators->first;
3042 operation.pre_constraints (local_data[comp][ind_local],
3045 const Number *data_val =
3047 const Number *end_pool =
3049 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3051 operation.process_constraint (*dof_indices, *data_val,
3052 *src[comp], value[comp]);
3055 operation.post_constraints (value[comp],
3056 local_data[comp][ind_local]);
3062 for (; ind_local < n_local_dofs; ++dof_indices, ++ind_local)
3065 operation.process_dof (*dof_indices, *src[comp],
3066 local_data[comp][ind_local]);
3074 static_cast<int>(n_local_dofs));
3075 for (
unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
3076 for (
unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
3078 operation.process_dof (dof_indices[j+v], *src[comp],
3079 local_data[comp][j+v]);
3089 for ( ; indicators != indicators_end; ++indicators)
3091 for (
unsigned int j=0; j<indicators->first; ++j)
3096 operation.process_dof (dof_indices[j], *src[comp],
3097 local_data[comp][ind_local]);
3102 >= n_irreg_components_filled)
3105 operation.process_empty (local_data[comp][ind_local]);
3109 dof_indices += indicators->first;
3115 operation.pre_constraints (local_data[comp][ind_local],
3118 const Number *data_val =
3120 const Number *end_pool =
3123 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3125 operation.process_constraint (*dof_indices, *data_val,
3126 *src[comp], value[comp]);
3129 operation.post_constraints (value[comp],
3130 local_data[comp][ind_local]);
3133 >= n_irreg_components_filled)
3136 operation.process_empty (local_data[comp][ind_local]);
3140 for (; ind_local<n_local_dofs; ++dof_indices)
3142 Assert (dof_indices != dof_info->end_indices(cell),
3148 operation.process_dof (*dof_indices, *src[comp],
3149 local_data[comp][ind_local]);
3152 >= n_irreg_components_filled)
3155 operation.process_empty(local_data[comp][ind_local]);
3167 internal::check_vector_compatibility (*src[0], *dof_info);
3169 const unsigned int n_local_dofs =
3171 Number *local_data =
3172 const_cast<Number *
>(&values_dofs[0][0][0]);
3173 if (at_irregular_cell ==
false)
3176 if (indicators != indicators_end)
3178 for ( ; indicators != indicators_end; ++indicators)
3181 for (
unsigned int j=0; j<indicators->first; ++j)
3182 operation.process_dof (dof_indices[j], *src[0],
3183 local_data[ind_local+j]);
3184 ind_local += indicators->first;
3185 dof_indices += indicators->first;
3190 operation.pre_constraints (local_data[ind_local], value);
3192 const Number *data_val =
3194 const Number *end_pool =
3197 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3198 operation.process_constraint (*dof_indices, *data_val,
3201 operation.post_constraints (value, local_data[ind_local]);
3206 for (; ind_local<n_local_dofs; ++dof_indices, ++ind_local)
3207 operation.process_dof (*dof_indices, *src[0],
3208 local_data[ind_local]);
3209 Assert (dof_indices == dof_info->end_indices(cell),
3217 static_cast<int>(n_local_dofs));
3218 for (
unsigned int j=0; j<n_local_dofs; j+=VectorizedArray<Number>::n_array_elements)
3219 for (
unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
3220 operation.process_dof (dof_indices[j+v], *src[0],
3231 for ( ; indicators != indicators_end; ++indicators)
3233 for (
unsigned int j=0; j<indicators->first; ++j)
3237 operation.process_dof (dof_indices[j], *src[0],
3238 local_data[ind_local]);
3243 >= n_irreg_components_filled)
3245 operation.process_empty (local_data[ind_local]);
3249 dof_indices += indicators->first;
3254 operation.pre_constraints (local_data[ind_local], value);
3256 const Number *data_val =
3258 const Number *end_pool =
3261 for ( ; data_val != end_pool; ++data_val, ++dof_indices)
3262 operation.process_constraint (*dof_indices, *data_val,
3265 operation.post_constraints (value, local_data[ind_local]);
3268 >= n_irreg_components_filled)
3270 operation.process_empty (local_data[ind_local]);
3274 for (; ind_local<n_local_dofs; ++dof_indices)
3276 Assert (dof_indices != dof_info->end_indices(cell),
3281 operation.process_dof (*dof_indices, *src[0],
3282 local_data[ind_local]);
3285 >= n_irreg_components_filled)
3287 operation.process_empty (local_data[ind_local]);
3297 template <
int dim,
int n_components_,
typename Number>
3298 template<
typename VectorType>
3303 const unsigned int first_index)
3307 typename internal::BlockVectorSelector<VectorType,
3310 src_data[d] = internal::BlockVectorSelector<VectorType,
IsBlockVector<VectorType>::value>::get_vector_component(const_cast<VectorType &>(src), d+first_index);
3312 internal::VectorReader<Number> reader;
3313 read_write_operation (reader, src_data);
3316 dof_values_initialized =
true;
3322 template <
int dim,
int n_components_,
typename Number>
3323 template<
typename VectorType>
3328 const unsigned int first_index)
3332 const typename internal::BlockVectorSelector<VectorType,
3335 src_data[d] = internal::BlockVectorSelector<VectorType,
IsBlockVector<VectorType>::value>::get_vector_component(const_cast<VectorType &>(src), d+first_index);
3337 read_dof_values_plain (src_data);
3342 template <
int dim,
int n_components_,
typename Number>
3343 template<
typename VectorType>
3348 const unsigned int first_index)
const 3350 Assert (dof_values_initialized==
true,
3355 typename internal::BlockVectorSelector<VectorType,
3360 internal::VectorDistributorLocalToGlobal<Number> distributor;
3361 read_write_operation (distributor, dst_data);
3366 template <
int dim,
int n_components_,
typename Number>
3367 template<
typename VectorType>
3372 const unsigned int first_index)
const 3374 Assert (dof_values_initialized==
true,
3379 typename internal::BlockVectorSelector<VectorType,
3384 internal::VectorSetter<Number> setter;
3385 read_write_operation (setter, dst_data);
3390 template <
int dim,
int n_components_,
typename Number>
3391 template<
typename VectorType>
3398 if (matrix_info == 0)
3400 internal::VectorReader<Number> reader;
3401 read_write_operation (reader, src);
3416 const unsigned int *dof_indices = dof_info->begin_indices_plain(cell);
3417 const unsigned int dofs_per_cell = this->data->
dofs_per_cell;
3419 const unsigned int n_irreg_components_filled = dof_info->row_starts[cell][2];
3420 const bool at_irregular_cell = n_irreg_components_filled > 0;
3424 if (n_fe_components == 1)
3426 const unsigned int n_local_dofs =
3429 internal::check_vector_compatibility (*src[comp], *dof_info);
3432 local_src_number[comp] = &values_dofs[comp][0][0];
3436 if (at_irregular_cell ==
false)
3438 for (
unsigned int j=0; j<n_local_dofs; ++j)
3440 local_src_number[comp][j] =
3441 internal::vector_access (*src[comp], dof_indices[j]);
3450 for (
unsigned int ind_local=0; ind_local<n_local_dofs;
3456 local_src_number[comp][ind_local] =
3457 internal::vector_access (*src[comp], *dof_indices);
3462 local_src_number[comp][ind_local] = 0.;
3474 internal::check_vector_compatibility (*src[0], *dof_info);
3476 const unsigned int n_local_dofs =
3478 Number *local_src_number = &values_dofs[0][0][0];
3479 if (at_irregular_cell ==
false)
3481 for (
unsigned int j=0; j<n_local_dofs; ++j)
3482 local_src_number[j] =
3483 internal::vector_access (*src[0], dof_indices[j]);
3492 for (
unsigned int ind_local=0; ind_local<n_local_dofs; ++dof_indices)
3496 local_src_number[ind_local] =
3497 internal::vector_access (*src[0], *dof_indices);
3501 local_src_number[ind_local] = 0.;
3509 dof_values_initialized =
true;
3518 template <
int dim,
int n_components,
typename Number>
3520 const std::vector<unsigned int> &
3524 return data->lexicographic_numbering;
3529 template <
int dim,
int n_components,
typename Number>
3542 template <
int dim,
int n_components,
typename Number>
3553 template <
int dim,
int n_components,
typename Number>
3567 template <
int dim,
int n_components,
typename Number>
3580 template <
int dim,
int n_components,
typename Number>
3594 template <
int dim,
int n_components,
typename Number>
3607 template <
int dim,
int n_components,
typename Number>
3621 template <
int dim,
int n_components,
typename Number>
3633 template <
int dim,
int n_components,
typename Number>
3644 template <
int dim,
int n_components_,
typename Number>
3653 return_value[comp] = this->values_dofs[comp][dof];
3654 return return_value;
3659 template <
int dim,
int n_components_,
typename Number>
3665 Assert (this->values_quad_initialized==
true,
3670 return_value[comp] = this->values_quad[comp][q_point];
3671 return return_value;
3676 template <
int dim,
int n_components_,
typename Number>
3682 Assert (this->gradients_quad_initialized==
true,
3689 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3692 for (
unsigned int d=0;
d<dim; ++
d)
3693 grad_out[comp][d] = (this->gradients_quad[comp][d][q_point] *
3694 cartesian_data[0][d]);
3700 this->cell_type == internal::MatrixFreeFunctions::general ?
3701 jacobian[q_point] : jacobian[0];
3704 for (
unsigned int d=0;
d<dim; ++
d)
3706 grad_out[comp][
d] = (jac[
d][0] *
3707 this->gradients_quad[comp][0][q_point]);
3708 for (
unsigned int e=1;
e<dim; ++
e)
3709 grad_out[comp][d] += (jac[d][e] *
3710 this->gradients_quad[comp][e][q_point]);
3723 template <
typename Number>
3728 const unsigned int q_point,
3731 tmp[0][0] = jac[0][0] * hessians_quad[0][q_point];
3734 template <
typename Number>
3739 const unsigned int q_point,
3742 for (
unsigned int d=0;
d<2; ++
d)
3744 tmp[0][
d] = (jac[
d][0] * hessians_quad[0][q_point] +
3745 jac[
d][1] * hessians_quad[2][q_point]);
3746 tmp[1][
d] = (jac[
d][0] * hessians_quad[2][q_point] +
3747 jac[
d][1] * hessians_quad[1][q_point]);
3751 template <
typename Number>
3756 const unsigned int q_point,
3759 for (
unsigned int d=0;
d<3; ++
d)
3761 tmp[0][
d] = (jac[
d][0] * hessians_quad[0][q_point] +
3762 jac[
d][1] * hessians_quad[3][q_point] +
3763 jac[
d][2] * hessians_quad[4][q_point]);
3764 tmp[1][
d] = (jac[
d][0] * hessians_quad[3][q_point] +
3765 jac[
d][1] * hessians_quad[1][q_point] +
3766 jac[
d][2] * hessians_quad[5][q_point]);
3767 tmp[2][
d] = (jac[
d][0] * hessians_quad[4][q_point] +
3768 jac[
d][1] * hessians_quad[5][q_point] +
3769 jac[
d][2] * hessians_quad[2][q_point]);
3776 template <
int dim,
int n_components_,
typename Number>
3782 Assert (this->hessians_quad_initialized==
true,
3789 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3793 for (
unsigned int d=0;
d<dim; ++
d)
3795 hessian_out[comp][
d][
d] = (this->hessians_quad[comp][
d][q_point] *
3802 hessian_out[comp][0][1] = (this->hessians_quad[comp][2][q_point] *
3806 hessian_out[comp][0][1] = (this->hessians_quad[comp][3][q_point] *
3808 hessian_out[comp][0][2] = (this->hessians_quad[comp][4][q_point] *
3810 hessian_out[comp][1][2] = (this->hessians_quad[comp][5][q_point] *
3816 for (
unsigned int e=d+1;
e<dim; ++
e)
3817 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3821 else if (this->cell_type == internal::MatrixFreeFunctions::general)
3829 & jac_grad_UT = jacobian_grad_upper[q_point];
3835 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3839 for (
unsigned int d=0;
d<dim; ++
d)
3840 for (
unsigned int e=d;
e<dim; ++
e)
3842 hessian_out[comp][
d][
e] = jac[
d][0] * tmp[0][
e];
3843 for (
unsigned int f=1; f<dim; ++f)
3844 hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3848 for (
unsigned int d=0;
d<dim; ++
d)
3849 for (
unsigned int e=0;
e<dim; ++
e)
3850 hessian_out[comp][d][d] += (jac_grad[d][e] *
3851 this->gradients_quad[comp][e][q_point]);
3854 for (
unsigned int d=0, count=0;
d<dim; ++
d)
3855 for (
unsigned int e=d+1;
e<dim; ++
e, ++count)
3856 for (
unsigned int f=0; f<dim; ++f)
3857 hessian_out[comp][d][e] += (jac_grad_UT[count][f] *
3858 this->gradients_quad[comp][f][q_point]);
3861 for (
unsigned int d=0;
d<dim; ++
d)
3862 for (
unsigned int e=d+1;
e<dim; ++
e)
3863 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3875 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3879 for (
unsigned int d=0;
d<dim; ++
d)
3880 for (
unsigned int e=d;
e<dim; ++
e)
3882 hessian_out[comp][
d][
e] = jac[
d][0] * tmp[0][
e];
3883 for (
unsigned int f=1; f<dim; ++f)
3884 hessian_out[comp][d][e] += jac[d][f] * tmp[f][e];
3891 for (
unsigned int d=0;
d<dim; ++
d)
3892 for (
unsigned int e=d+1;
e<dim; ++
e)
3893 hessian_out[comp][e][d] = hessian_out[comp][d][e];
3901 template <
int dim,
int n_components_,
typename Number>
3907 Assert (this->hessians_quad_initialized==
true,
3914 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
3918 for (
unsigned int d=0;
d<dim; ++
d)
3919 hessian_out[comp][d] = (this->hessians_quad[comp][d][q_point] *
3923 else if (this->cell_type == internal::MatrixFreeFunctions::general)
3934 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3939 for (
unsigned int d=0;
d<dim; ++
d)
3941 hessian_out[comp][
d] = jac[
d][0] * tmp[0][
d];
3942 for (
unsigned int f=1; f<dim; ++f)
3943 hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3946 for (
unsigned int d=0;
d<dim; ++
d)
3947 for (
unsigned int e=0;
e<dim; ++
e)
3948 hessian_out[comp][d] += (jac_grad[d][e] *
3949 this->gradients_quad[comp][e][q_point]);
3961 internal::hessian_unit_times_jac (jac, this->hessians_quad[comp],
3966 for (
unsigned int d=0;
d<dim; ++
d)
3968 hessian_out[comp][
d] = jac[
d][0] * tmp[0][
d];
3969 for (
unsigned int f=1; f<dim; ++f)
3970 hessian_out[comp][d] += jac[d][f] * tmp[f][d];
3979 template <
int dim,
int n_components_,
typename Number>
3985 Assert (this->hessians_quad_initialized==
true,
3990 = get_hessian_diagonal(q_point);
3993 laplacian_out[comp] = hess_diag[comp][0];
3994 for (
unsigned int d=1;
d<dim; ++
d)
3995 laplacian_out[comp] += hess_diag[comp][d];
3997 return laplacian_out;
4002 template <
int dim,
int n_components_,
typename Number>
4007 const unsigned int dof)
4010 this->dof_values_initialized =
true;
4014 this->values_dofs[comp][dof] = val_in[comp];
4019 template <
int dim,
int n_components_,
typename Number>
4024 const unsigned int q_point)
4029 this->values_quad_submitted =
true;
4031 if (this->cell_type == internal::MatrixFreeFunctions::general)
4035 this->values_quad[comp][q_point] = val_in[comp] * JxW;
4041 this->values_quad[comp][q_point] = val_in[comp] * JxW;
4047 template <
int dim,
int n_components_,
typename Number>
4053 const unsigned int q_point)
4058 this->gradients_quad_submitted =
true;
4060 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4064 for (
unsigned int d=0;
d<dim; ++
d)
4065 this->gradients_quad[comp][d][q_point] = (grad_in[comp][d] *
4066 cartesian_data[0][d] * JxW);
4071 this->cell_type == internal::MatrixFreeFunctions::general ?
4072 jacobian[q_point] : jacobian[0];
4074 this->cell_type == internal::MatrixFreeFunctions::general ?
4075 J_value[q_point] : J_value[0] * quadrature_weights[q_point];
4077 for (
unsigned int d=0;
d<dim; ++
d)
4080 for (
unsigned int e=1;
e<dim; ++
e)
4081 new_val += (jac[e][d] * grad_in[comp][e]);
4082 this->gradients_quad[comp][
d][q_point] = new_val * JxW;
4089 template <
int dim,
int n_components_,
typename Number>
4097 Assert (this->values_quad_submitted ==
true,
4102 return_value[comp] = this->values_quad[comp][0];
4103 const unsigned int n_q_points = this->data->
n_q_points;
4104 for (
unsigned int q=1; q<n_q_points; ++q)
4106 return_value[comp] += this->values_quad[comp][q];
4107 return (return_value);
4115 template <
int dim,
int n_components_,
typename Number>
4119 const unsigned int fe_no,
4120 const unsigned int quad_no_in,
4121 const unsigned int fe_degree,
4122 const unsigned int n_q_points)
4125 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4130 template <
int dim,
int n_components_,
typename Number>
4131 template <
int n_components_other>
4138 const unsigned int first_selected_component,
4141 FEEvaluationBase <dim,n_components_,Number>(mapping, fe, quadrature, update_flags,
4142 first_selected_component, other)
4147 template <
int dim,
int n_components_,
typename Number>
4157 template <
int dim,
int n_components_,
typename Number>
4172 template <
int dim,
typename Number>
4176 const unsigned int fe_no,
4177 const unsigned int quad_no_in,
4178 const unsigned int fe_degree,
4179 const unsigned int n_q_points)
4182 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4187 template <
int dim,
typename Number>
4188 template <
int n_components_other>
4195 const unsigned int first_selected_component,
4199 first_selected_component, other)
4204 template <
int dim,
typename Number>
4214 template <
int dim,
typename Number>
4226 template <
int dim,
typename Number>
4233 return this->values_dofs[0][dof];
4238 template <
int dim,
typename Number>
4244 Assert (this->values_quad_initialized==
true,
4247 return this->values_quad[0][q_point];
4252 template <
int dim,
typename Number>
4261 Assert (this->gradients_quad_initialized==
true,
4268 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4270 for (
unsigned int d=0;
d<dim; ++
d)
4271 grad_out[d] = (this->gradients_quad[0][d][q_point] *
4272 this->cartesian_data[0][d]);
4278 this->cell_type == internal::MatrixFreeFunctions::general ?
4279 this->jacobian[q_point] : this->jacobian[0];
4280 for (
unsigned int d=0;
d<dim; ++
d)
4282 grad_out[
d] = (jac[
d][0] * this->gradients_quad[0][0][q_point]);
4283 for (
unsigned int e=1;
e<dim; ++
e)
4284 grad_out[d] += (jac[d][e] * this->gradients_quad[0][e][q_point]);
4292 template <
int dim,
typename Number>
4298 return BaseClass::get_hessian(q_point)[0];
4303 template <
int dim,
typename Number>
4309 return BaseClass::get_hessian_diagonal(q_point)[0];
4314 template <
int dim,
typename Number>
4320 return BaseClass::get_laplacian(q_point)[0];
4325 template <
int dim,
typename Number>
4330 const unsigned int dof)
4333 this->dof_values_initialized =
true;
4336 this->values_dofs[0][dof] = val_in;
4341 template <
int dim,
typename Number>
4346 const unsigned int q_point)
4351 this->values_quad_submitted =
true;
4353 if (this->cell_type == internal::MatrixFreeFunctions::general)
4356 this->values_quad[0][q_point] = val_in * JxW;
4361 this->values_quad[0][q_point] = val_in * JxW;
4367 template <
int dim,
typename Number>
4372 const unsigned int q_point)
4377 this->gradients_quad_submitted =
true;
4379 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4382 for (
unsigned int d=0;
d<dim; ++
d)
4383 this->gradients_quad[0][d][q_point] = (grad_in[d] *
4384 this->cartesian_data[0][d] *
4391 this->cell_type == internal::MatrixFreeFunctions::general ?
4392 this->jacobian[q_point] : this->jacobian[0];
4394 this->cell_type == internal::MatrixFreeFunctions::general ?
4395 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4396 for (
unsigned int d=0;
d<dim; ++
d)
4399 for (
unsigned int e=1;
e<dim; ++
e)
4400 new_val += jac[e][d] * grad_in[e];
4401 this->gradients_quad[0][
d][q_point] = new_val * JxW;
4408 template <
int dim,
typename Number>
4414 return BaseClass::integrate_value()[0];
4423 template <
int dim,
typename Number>
4427 const unsigned int fe_no,
4428 const unsigned int quad_no_in,
4429 const unsigned int fe_degree,
4430 const unsigned int n_q_points)
4433 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4438 template <
int dim,
typename Number>
4439 template <
int n_components_other>
4446 const unsigned int first_selected_component,
4450 first_selected_component, other)
4455 template <
int dim,
typename Number>
4465 template <
int dim,
typename Number>
4477 template <
int dim,
typename Number>
4483 return BaseClass::get_gradient (q_point);
4488 template <
int dim,
typename Number>
4494 Assert (this->gradients_quad_initialized==
true,
4501 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4503 divergence = (this->gradients_quad[0][0][q_point] *
4504 this->cartesian_data[0][0]);
4505 for (
unsigned int d=1;
d<dim; ++
d)
4506 divergence += (this->gradients_quad[d][d][q_point] *
4507 this->cartesian_data[0][d]);
4513 this->cell_type == internal::MatrixFreeFunctions::general ?
4514 this->jacobian[q_point] : this->jacobian[0];
4515 divergence = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4516 for (
unsigned int e=1;
e<dim; ++
e)
4517 divergence += (jac[0][e] * this->gradients_quad[0][e][q_point]);
4518 for (
unsigned int d=1;
d<dim; ++
d)
4519 for (
unsigned int e=0;
e<dim; ++
e)
4520 divergence += (jac[d][e] * this->gradients_quad[d][e][q_point]);
4527 template <
int dim,
typename Number>
4537 for (
unsigned int d=0;
d<dim; ++
d)
4538 symmetrized[d] = grad[d][d];
4544 symmetrized[2] = grad[0][1] + grad[1][0];
4545 symmetrized[2] *= half;
4548 symmetrized[3] = grad[0][1] + grad[1][0];
4549 symmetrized[3] *= half;
4550 symmetrized[4] = grad[0][2] + grad[2][0];
4551 symmetrized[4] *= half;
4552 symmetrized[5] = grad[1][2] + grad[2][1];
4553 symmetrized[5] *= half;
4563 template <
int dim,
typename Number>
4567 ::get_curl (
const unsigned int q_point)
const 4576 ExcMessage(
"Computing the curl in 1d is not a useful operation"));
4579 curl[0] = grad[1][0] - grad[0][1];
4582 curl[0] = grad[2][1] - grad[1][2];
4583 curl[1] = grad[0][2] - grad[2][0];
4584 curl[2] = grad[1][0] - grad[0][1];
4594 template <
int dim,
typename Number>
4600 Assert (this->hessians_quad_initialized==
true,
4604 return BaseClass::get_hessian_diagonal (q_point);
4609 template <
int dim,
typename Number>
4615 Assert (this->hessians_quad_initialized==
true,
4618 return BaseClass::get_hessian(q_point);
4623 template <
int dim,
typename Number>
4628 const unsigned int q_point)
4630 BaseClass::submit_gradient (grad_in, q_point);
4635 template <
int dim,
typename Number>
4641 const unsigned int q_point)
4643 BaseClass::submit_gradient(grad_in, q_point);
4648 template <
int dim,
typename Number>
4653 const unsigned int q_point)
4658 this->gradients_quad_submitted =
true;
4660 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4663 this->quadrature_weights[q_point] * div_in;
4664 for (
unsigned int d=0;
d<dim; ++
d)
4666 this->gradients_quad[
d][
d][q_point] = (fac *
4667 this->cartesian_data[0][
d]);
4668 for (
unsigned int e=d+1;
e<dim; ++
e)
4678 this->cell_type == internal::MatrixFreeFunctions::general ?
4679 this->jacobian[q_point] : this->jacobian[0];
4681 (this->cell_type == internal::MatrixFreeFunctions::general ?
4682 this->J_value[q_point] : this->J_value[0] *
4683 this->quadrature_weights[q_point]) * div_in;
4684 for (
unsigned int d=0;
d<dim; ++
d)
4686 for (
unsigned int e=0;
e<dim; ++
e)
4687 this->gradients_quad[d][e][q_point] = jac[d][e] * fac;
4694 template <
int dim,
typename Number>
4700 const unsigned int q_point)
4708 this->gradients_quad_submitted =
true;
4710 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4713 for (
unsigned int d=0;
d<dim; ++
d)
4714 this->gradients_quad[d][d][q_point] = (sym_grad.access_raw_entry(d) *
4716 this->cartesian_data[0][
d]);
4717 for (
unsigned int e=0, counter=dim;
e<dim; ++
e)
4718 for (
unsigned int d=e+1;
d<dim; ++
d, ++counter)
4721 this->gradients_quad[
e][
d][q_point] = (value *
4722 this->cartesian_data[0][
d]);
4723 this->gradients_quad[
d][
e][q_point] = (value *
4724 this->cartesian_data[0][
e]);
4731 this->cell_type == internal::MatrixFreeFunctions::general ?
4732 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
4734 this->cell_type == internal::MatrixFreeFunctions::general ?
4735 this->jacobian[q_point] : this->jacobian[0];
4737 for (
unsigned int i=0; i<dim; ++i)
4738 weighted[i][i] = sym_grad.access_raw_entry(i) * JxW;
4739 for (
unsigned int i=0, counter=dim; i<dim; ++i)
4740 for (
unsigned int j=i+1; j<dim; ++j, ++counter)
4743 weighted[i][j] = value;
4744 weighted[j][i] = value;
4746 for (
unsigned int comp=0; comp<dim; ++comp)
4747 for (
unsigned int d=0;
d<dim; ++
d)
4750 for (
unsigned int e=1;
e<dim; ++
e)
4751 new_val += jac[e][d] * weighted[comp][e];
4752 this->gradients_quad[comp][
d][q_point] = new_val;
4759 template <
int dim,
typename Number>
4764 const unsigned int q_point)
4771 ExcMessage(
"Testing by the curl in 1d is not a useful operation"));
4774 grad[1][0] = curl[0];
4775 grad[0][1] = -curl[0];
4778 grad[2][1] = curl[0];
4779 grad[1][2] = -curl[0];
4780 grad[0][2] = curl[1];
4781 grad[2][0] = -curl[1];
4782 grad[1][0] = curl[2];
4783 grad[0][1] = -curl[2];
4788 submit_gradient (grad, q_point);
4795 template <
typename Number>
4799 const unsigned int fe_no,
4800 const unsigned int quad_no_in,
4801 const unsigned int fe_degree,
4802 const unsigned int n_q_points)
4805 (data_in, fe_no, quad_no_in, fe_degree, n_q_points)
4810 template <
typename Number>
4811 template <
int n_components_other>
4818 const unsigned int first_selected_component,
4822 first_selected_component, other)
4827 template <
typename Number>
4837 template <
typename Number>
4849 template <
typename Number>
4856 return this->values_dofs[0][dof];
4861 template <
typename Number>
4867 Assert (this->values_quad_initialized==
true,
4870 return this->values_quad[0][q_point];
4875 template <
typename Number>
4884 Assert (this->gradients_quad_initialized==
true,
4891 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
4893 grad_out[0] = (this->gradients_quad[0][0][q_point] *
4894 this->cartesian_data[0][0]);
4900 this->cell_type == internal::MatrixFreeFunctions::general ?
4901 this->jacobian[q_point] : this->jacobian[0];
4903 grad_out[0] = (jac[0][0] * this->gradients_quad[0][0][q_point]);
4910 template <
typename Number>
4916 return BaseClass::get_hessian(q_point)[0];
4921 template <
typename Number>
4927 return BaseClass::get_hessian_diagonal(q_point)[0];
4932 template <
typename Number>
4938 return BaseClass::get_laplacian(q_point)[0];
4943 template <
typename Number>
4948 const unsigned int dof)
4951 this->dof_values_initialized =
true;
4954 this->values_dofs[0][dof] = val_in;
4959 template <
typename Number>
4964 const unsigned int q_point)
4969 this->values_quad_submitted =
true;
4971 if (this->cell_type == internal::MatrixFreeFunctions::general)
4974 this->values_quad[0][q_point] = val_in * JxW;
4979 this->values_quad[0][q_point] = val_in * JxW;
4985 template <
typename Number>
4990 const unsigned int q_point)
4995 this->gradients_quad_submitted =
true;
4997 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
5000 this->gradients_quad[0][0][q_point] = (grad_in[0] *
5001 this->cartesian_data[0][0] *
5008 this->cell_type == internal::MatrixFreeFunctions::general ?
5009 this->jacobian[q_point] : this->jacobian[0];
5011 this->cell_type == internal::MatrixFreeFunctions::general ?
5012 this->J_value[q_point] : this->J_value[0] * this->quadrature_weights[q_point];
5014 this->gradients_quad[0][0][q_point] = jac[0][0] * grad_in[0] * JxW;
5020 template <
typename Number>
5026 return BaseClass::integrate_value()[0];
5035 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5040 const unsigned int fe_no,
5041 const unsigned int quad_no)
5043 BaseClass (data_in, fe_no, quad_no, fe_degree, static_n_q_points),
5044 dofs_per_cell (this->data->dofs_per_cell),
5045 n_q_points (this->data->n_q_points)
5047 check_template_arguments(fe_no, 0);
5052 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5060 const unsigned int first_selected_component)
5062 BaseClass (mapping, fe, quadrature, update_flags,
5063 first_selected_component,
5065 dofs_per_cell (this->data->dofs_per_cell),
5066 n_q_points (this->data->n_q_points)
5073 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5080 const unsigned int first_selected_component)
5082 BaseClass (
StaticMappingQ1<dim>::mapping, fe, quadrature, update_flags,
5083 first_selected_component,
5085 dofs_per_cell (this->data->dofs_per_cell),
5086 n_q_points (this->data->n_q_points)
5093 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5095 template <
int n_components_other>
5100 const unsigned int first_selected_component)
5102 BaseClass (other.mapped_geometry->get_fe_values().get_mapping(),
5103 fe, other.mapped_geometry->get_quadrature(),
5104 other.mapped_geometry->get_fe_values().get_update_flags(),
5105 first_selected_component, &other),
5106 dofs_per_cell (this->data->dofs_per_cell),
5107 n_q_points (this->data->n_q_points)
5114 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5121 dofs_per_cell (this->data->dofs_per_cell),
5122 n_q_points (this->data->n_q_points)
5129 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5136 BaseClass::operator=(other);
5143 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5149 const unsigned int first_selected_component)
5151 const unsigned int fe_degree_templ = fe_degree != -1 ? fe_degree : 0;
5152 const unsigned int n_q_points_1d_templ = n_q_points_1d > 0 ? n_q_points_1d : 1;
5153 if (fe_degree == -1)
5155 if (this->data->
element_type == internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0)
5157 evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5158 dim, -1, 0, n_components_, Number>::evaluate;
5159 integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5160 dim, -1, 0, n_components_, Number>::integrate;
5162 else if (this->data->
element_type == internal::MatrixFreeFunctions::truncated_tensor)
5164 evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5165 dim, -1, 0, n_components_, Number>::evaluate;
5166 integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5167 dim, -1, 0, n_components_, Number>::integrate;
5171 evaluate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5172 dim, -1, 0, n_components_, Number>::evaluate;
5173 integrate_funct = internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5174 dim, -1, 0, n_components_, Number>::integrate;
5180 case internal::MatrixFreeFunctions::tensor_symmetric:
5182 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
5183 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5186 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric,
5187 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5191 case internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0:
5193 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5194 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5197 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_symmetric_plus_dg0,
5198 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5202 case internal::MatrixFreeFunctions::tensor_general:
5204 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5205 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5208 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_general,
5209 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5213 case internal::MatrixFreeFunctions::tensor_gausslobatto:
5215 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
5216 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5219 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::tensor_gausslobatto,
5220 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5224 case internal::MatrixFreeFunctions::truncated_tensor:
5226 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5227 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5230 internal::FEEvaluationImpl<internal::MatrixFreeFunctions::truncated_tensor,
5231 dim, fe_degree_templ, n_q_points_1d_templ, n_components_,
5240 (void)first_selected_component;
5245 if ((fe_degree != -1 && static_cast<unsigned int>(fe_degree) != this->data->
fe_degree)
5247 (fe_degree != -1 && static_n_q_points != this->data->
n_q_points))
5249 std::string message =
5250 "-------------------------------------------------------\n";
5251 message +=
"Illegal arguments in constructor/wrong template arguments!\n";
5252 message +=
" Called --> FEEvaluation<dim,";
5256 message +=
",Number>(data";
5270 if (static_cast<unsigned int>(fe_degree) == this->data->
fe_degree)
5271 proposed_dof_comp = fe_no;
5273 for (
unsigned int no=0; no<this->matrix_info->
n_components(); ++no)
5274 if (this->matrix_info->
get_shape_info(no,0,this->active_fe_index,0).fe_degree
5275 ==
static_cast<unsigned int>(fe_degree))
5277 proposed_dof_comp = no;
5280 if (static_n_q_points ==
5281 this->mapping_info->
mapping_data_gen[this->quad_no].n_q_points[this->active_quad_index])
5282 proposed_quad_comp = this->quad_no;
5284 for (
unsigned int no=0; no<this->mapping_info->
mapping_data_gen.size(); ++no)
5285 if (this->mapping_info->
mapping_data_gen[no].n_q_points[this->active_quad_index]
5286 == static_n_q_points)
5288 proposed_quad_comp = no;
5295 if (proposed_dof_comp != fe_no)
5296 message +=
"Wrong vector component selection:\n";
5298 message +=
"Wrong quadrature formula selection:\n";
5299 message +=
" Did you mean FEEvaluation<dim,";
5303 message +=
",Number>(data";
5310 std::string correct_pos;
5311 if (proposed_dof_comp != fe_no)
5312 correct_pos =
" ^ ";
5315 if (proposed_quad_comp != this->quad_no)
5316 correct_pos +=
" ^\n";
5318 correct_pos +=
" \n";
5319 message +=
" " + correct_pos;
5323 const unsigned int proposed_n_q_points_1d =
static_cast<unsigned int>(std::pow(1.001*this->data->
n_q_points,1./dim));
5324 message +=
"Wrong template arguments:\n";
5325 message +=
" Did you mean FEEvaluation<dim,";
5329 message +=
",Number>(data";
5336 std::string correct_pos;
5337 if (this->data->
fe_degree != static_cast<unsigned int>(fe_degree))
5341 if (proposed_n_q_points_1d != n_q_points_1d)
5342 correct_pos +=
" ^\n";
5344 correct_pos +=
" \n";
5345 message +=
" " + correct_pos;
5347 Assert (static_cast<unsigned int>(fe_degree) == this->data->
fe_degree &&
5348 static_n_q_points == this->data->n_q_points,
5355 n_q_points[this->active_quad_index]);
5357 this->dof_info->dofs_per_cell[this->active_fe_index]);
5364 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5379 if (this->cell_type == internal::MatrixFreeFunctions::cartesian)
5385 return this->quadrature_points[q];
5387 point[0] = this->quadrature_points[q%n_q_points_1d][0];
5388 point[1] = this->quadrature_points[q/n_q_points_1d][1];
5391 point[0] = this->quadrature_points[q%n_q_points_1d][0];
5392 point[1] = this->quadrature_points[(q/n_q_points_1d)%n_q_points_1d][1];
5393 point[2] = this->quadrature_points[q/(n_q_points_1d*n_q_points_1d)][2];
5402 return this->quadrature_points[q];
5407 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5413 const bool evaluate_grad,
5414 const bool evaluate_lapl)
5416 Assert (this->dof_values_initialized ==
true,
5418 Assert(this->matrix_info != 0 ||
5423 evaluate_funct (*this->data, &this->values_dofs[0], this->values_quad,
5424 this->gradients_quad, this->hessians_quad, this->scratch_data,
5425 evaluate_val, evaluate_grad, evaluate_lapl);
5428 if (evaluate_val ==
true)
5429 this->values_quad_initialized =
true;
5430 if (evaluate_grad ==
true)
5431 this->gradients_quad_initialized =
true;
5432 if (evaluate_lapl ==
true)
5433 this->hessians_quad_initialized =
true;
5439 template <
int dim,
int fe_degree,
int n_q_points_1d,
int n_components_,
5444 ::integrate (
bool integrate_val,
bool integrate_grad)
5446 if (integrate_val ==
true)
5447 Assert (this->values_quad_submitted ==
true,
5449 if (integrate_grad ==
true)
5450 Assert (this->gradients_quad_submitted ==
true,
5452 Assert(this->matrix_info != 0 ||
5457 integrate_funct (*this->data, this->values_dofs, this->values_quad,
5458 this->gradients_quad, this->scratch_data,
5459 integrate_val, integrate_grad);
5462 this->dof_values_initialized =
true;
5468 #endif // ifndef DOXYGEN 5471 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
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
Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArray< Number > > > get_hessian(const unsigned int q_point) const
void check_template_arguments(const unsigned int fe_no, const unsigned int first_selected_component)
const MatrixFree< dim, Number > * matrix_info
FEEvaluation & operator=(const FEEvaluation &other)
void reinit(const unsigned int cell)
const Tensor< 1,(dim >1?dim *(dim-1)/2:1), Tensor< 1, dim, VectorizedArray< Number > > > * jacobian_grad_upper
AlignedVector< std::pair< Tensor< 2, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > affine_data
unsigned int cell_data_number
void read_write_operation(const VectorOperation &operation, VectorType *vectors[]) const
static ::ExceptionBase & ExcAccessToUninitializedField()
const unsigned int first_selected_component
value_type get_laplacian(const unsigned int q_point) const
bool quadrature_points_initialized
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
const unsigned int dofs_per_cell
unsigned int n_components() const
VectorizedArray< Number > JxW(const unsigned int q_point) const
void evaluate(const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians=false)
const unsigned int n_fe_components
void read_dof_values(const VectorType &src, const unsigned int first_index=0)
void integrate(const bool integrate_values, const bool integrate_gradients)
void read_dof_values_plain(const VectorType &src, const unsigned int first_index=0)
bool gradients_quad_initialized
#define AssertIndexRange(index, range)
const VectorizedArray< Number > * begin_gradients() const
AlignedVector< VectorizedArray< Number > > * scratch_data_array
std::vector< types::global_dof_index > old_style_dof_indices
unsigned int get_cell_data_index(const unsigned int cell_chunk_no) const
const VectorizedArray< Number > * begin_hessians() const
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian
friend class FEEvaluationBase
static ::ExceptionBase & ExcNotInitialized()
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)
SymmetricTensor< 2, dim, VectorizedArray< Number > > get_symmetric_gradient(const unsigned int q_point) const
void submit_symmetric_gradient(const SymmetricTensor< 2, dim, VectorizedArray< Number > > grad_in, const unsigned int q_point)
const VectorizedArray< Number > * begin_dof_values() const
VectorizedArray< Number > * scratch_data
const Tensor< 1, dim, VectorizedArray< Number > > * cartesian_data
void submit_dof_value(const value_type val_in, const unsigned int dof)
const internal::MatrixFreeFunctions::ShapeInfo< Number > * data
VectorizedArray< Number > * values_quad[n_components]
AlignedVector< VectorizedArray< Number > > * acquire_scratch_data() const
FEEvaluationBase & operator=(const FEEvaluationBase &other)
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int n_q_points
void fill_JxW_values(AlignedVector< VectorizedArray< Number > > &JxW_values) const
value_type get_dof_value(const unsigned int dof) const
void submit_divergence(const VectorizedArray< Number > div_in, const unsigned int q_point)
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
T sum(const T &t, const MPI_Comm &mpi_communicator)
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
void submit_curl(const Tensor< 1, dim==2?1:dim, VectorizedArray< Number > > curl_in, const unsigned int q_point)
bool mapping_initialized() const
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
bool JxW_values_initialized
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0) const
const VectorizedArray< Number > * quadrature_weights
gradient_type get_hessian_diagonal(const unsigned int q_point) const
bool hessians_quad_initialized
Tensor< 1,(dim==2?1:dim), VectorizedArray< Number > > get_curl(const unsigned int q_point) const
const internal::MatrixFreeFunctions::SizeInfo & get_size_info() const
FEEvaluationAccess & operator=(const FEEvaluationAccess &other)
void set_dof_values(VectorType &dst, const unsigned int first_index=0) const
void release_scratch_data(const AlignedVector< VectorizedArray< Number > > *memory) const
value_type integrate_value() const
void submit_value(const value_type val_in, const unsigned int q_point)
internal::MatrixFreeFunctions::CellType cell_type
void(* evaluate_funct)(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], VectorizedArray< Number > *hessians_quad[][(dim *(dim+1))/2], VectorizedArray< Number > *scratch_data, const bool evaluate_val, const bool evaluate_grad, const bool evaluate_lapl)
VectorizedArray< Number > get_divergence(const unsigned int q_point) const
gradient_type get_gradient(const unsigned int q_point) const
VectorizedArray< Number > * values_dofs[n_components]
ArrayView< VectorizedArray< Number > > get_scratch_data() const
const VectorizedArray< Number > * begin_values() const
VectorizedArray< Number > * hessians_quad[n_components][(dim *(dim+1))/2]
bool partitioners_are_compatible(const Utilities::MPI::Partitioner &part) const
unsigned int dofs_per_cell
bool second_derivatives_initialized
Number local_element(const size_type local_index) const
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info() const
static ::ExceptionBase & ExcNotImplemented()
bool indices_initialized() const
std::vector< unsigned int > lexicographic_numbering
bool values_quad_initialized
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
const std::vector< unsigned int > & get_internal_dof_numbering() const
value_type get_value(const unsigned int q_point) const
const VectorizedArray< Number > * J_value
CellType get_cell_type(const unsigned int cell_chunk_no) const
Point< 3 > point(const gp_Pnt &p)
VectorizedArray< Number > * gradients_quad[n_components][dim]
internal::MatrixFreeFunctions::CellType get_cell_type() const
Point< dim, VectorizedArray< Number > > quadrature_point(const unsigned int q_point) const
void(* integrate_funct)(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], VectorizedArray< Number > *scratch_data, const bool evaluate_val, const bool evaluate_grad)
bool gradients_quad_submitted
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > make_vectorized_array(const Number &u)
std_cxx1x::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
static ::ExceptionBase & ExcInternalError()
bool dof_values_initialized