Table<2, Tensor<3, dim>> shape_3rd_derivatives;
};
+ /**
+ * Correct the shape Hessians by subtracting the terms corresponding to the
+ * Jacobian pushed forward gradient.
+ *
+ * Before the correction, the Hessians would be given by
+ * @f[
+ * D_{ijk} = \frac{d^2\phi_i}{d \hat x_J d \hat x_K} (J_{jJ})^{-1}
+ * (J_{kK})^{-1},
+ * @f]
+ * where $J_{iI}=\frac{d x_i}{d \hat x_I}$. After the correction, the
+ * correct Hessians would be given by
+ * @f[
+ * \frac{d^2 \phi_i}{d x_j d x_k} = D_{ijk} - H_{mjk} \frac{d \phi_i}{d x_m},
+ * @f]
+ * where $H_{ijk}$ is the Jacobian pushed-forward derivative.
+ */
+ void
+ correct_hessians(
+ internal::FEValuesImplementation::FiniteElementRelatedData<dim, spacedim>
+ &output_data,
+ const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
+ & mapping_data,
+ const unsigned int n_q_points) const;
+
/**
* Correct the shape third derivatives by subtracting the terms
* corresponding to the Jacobian pushed forward gradient and second
const unsigned int n_q_points,
const unsigned int dof) const;
+
/**
* The polynomial space. Its type is given by the template parameter
* PolynomialType.
mapping_internal,
make_array_view(output_data.shape_hessians, k));
- for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
- for (unsigned int i = 0; i < quadrature.size(); ++i)
- output_data.shape_hessians[k][i] -=
- output_data.shape_gradients[k][i] *
- mapping_data.jacobian_pushed_forward_grads[i];
+ correct_hessians(output_data, mapping_data, quadrature.size());
}
if (flags & update_3rd_derivatives &&
mapping_internal,
make_array_view(output_data.shape_hessians, k));
- for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
- for (unsigned int i = 0; i < quadrature.size(); ++i)
- for (unsigned int j = 0; j < spacedim; ++j)
- output_data.shape_hessians[k][i] -=
- mapping_data.jacobian_pushed_forward_grads[i][j] *
- output_data.shape_gradients[k][i][j];
+ correct_hessians(output_data, mapping_data, quadrature.size());
}
if (flags & update_3rd_derivatives)
mapping_internal,
make_array_view(output_data.shape_hessians, k));
- for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
- for (unsigned int i = 0; i < quadrature.size(); ++i)
- for (unsigned int j = 0; j < spacedim; ++j)
- output_data.shape_hessians[k][i] -=
- mapping_data.jacobian_pushed_forward_grads[i][j] *
- output_data.shape_gradients[k][i][j];
+ correct_hessians(output_data, mapping_data, quadrature.size());
}
if (flags & update_3rd_derivatives)
+template <class PolynomialType, int dim, int spacedim>
+inline void
+FE_Poly<PolynomialType, dim, spacedim>::correct_hessians(
+ internal::FEValuesImplementation::FiniteElementRelatedData<dim, spacedim>
+ &output_data,
+ const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
+ & mapping_data,
+ const unsigned int n_q_points) const
+{
+ for (unsigned int dof = 0; dof < this->dofs_per_cell; ++dof)
+ for (unsigned int i = 0; i < n_q_points; ++i)
+ for (unsigned int j = 0; j < spacedim; ++j)
+ output_data.shape_hessians[dof][i] -=
+ mapping_data.jacobian_pushed_forward_grads[i][j] *
+ output_data.shape_gradients[dof][i][j];
+}
+
+
+
template <class PolynomialType, int dim, int spacedim>
inline void
FE_Poly<PolynomialType, dim, spacedim>::correct_third_derivatives(
mapping_internal,
make_array_view(output_data.shape_hessians, k));
- for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
- for (unsigned int i = 0; i < quadrature.size(); ++i)
- for (unsigned int j = 0; j < 2; ++j)
- output_data.shape_hessians[k][i] -=
- mapping_data.jacobian_pushed_forward_grads[i][j] *
- output_data.shape_gradients[k][i][j];
+ correct_hessians(output_data, mapping_data, quadrature.size());
}
if (fe_data.update_each & update_3rd_derivatives &&
mapping_internal,
make_array_view(output_data.shape_hessians, k));
- for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
- for (unsigned int i = 0; i < quadrature.size(); ++i)
- for (unsigned int j = 0; j < 3; ++j)
- output_data.shape_hessians[k][i] -=
- mapping_data.jacobian_pushed_forward_grads[i][j] *
- output_data.shape_gradients[k][i][j];
+ correct_hessians(output_data, mapping_data, quadrature.size());
}
if (fe_data.update_each & update_3rd_derivatives &&