From: Daniel Arndt Date: Thu, 29 Jun 2023 16:21:27 +0000 (-0400) Subject: Remove deprecated functions in fe_interface_values.h X-Git-Tag: relicensing~839^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=b0977a3349181bc2e22bd780bc5526237191ff2a;p=dealii.git Remove deprecated functions in fe_interface_values.h --- diff --git a/doc/news/changes/incompatibilities/20230629DanielArndt-2 b/doc/news/changes/incompatibilities/20230629DanielArndt-2 new file mode 100644 index 0000000000..a5f3eda053 --- /dev/null +++ b/doc/news/changes/incompatibilities/20230629DanielArndt-2 @@ -0,0 +1,3 @@ +Removed: Deprecated functions in the FEInterfaceValues nampespace have been removed. +
+(Daniel Arndt, 2023/06/29) diff --git a/include/deal.II/fe/fe_interface_values.h b/include/deal.II/fe/fe_interface_values.h index e032f72d35..dc6684e229 100644 --- a/include/deal.II/fe/fe_interface_values.h +++ b/include/deal.II/fe/fe_interface_values.h @@ -198,16 +198,6 @@ namespace FEInterfaceViews jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_values() function instead. - */ - DEAL_II_DEPRECATED - value_type - jump(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the jump of the gradient $\jump{nabla u}$ on the interface for * the shape @@ -222,16 +212,6 @@ namespace FEInterfaceViews jump_in_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_gradients() function instead. - */ - DEAL_II_DEPRECATED - gradient_type - jump_gradient(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the jump in the gradient $\jump{\nabla u}=\nabla u_{\text{cell0}} * - \nabla u_{\text{cell1}}$ on the interface for the shape function @p @@ -246,16 +226,6 @@ namespace FEInterfaceViews jump_in_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_hessians() function instead. - */ - DEAL_II_DEPRECATED - hessian_type - jump_hessian(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the jump in the third derivative $\jump{\nabla^3 u} = \nabla^3 * u_{\text{cell0}} - \nabla^3 u_{\text{cell1}}$ on the interface for the @@ -270,16 +240,6 @@ namespace FEInterfaceViews jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_third_derivatives() function instead. - */ - DEAL_II_DEPRECATED - third_derivative_type - jump_3rd_derivative(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** @} */ /** @@ -301,26 +261,6 @@ namespace FEInterfaceViews average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the average_of_values() function instead. - */ - DEAL_II_DEPRECATED - value_type - average_value(const unsigned int interface_dof_index, - const unsigned int q_point) const; - - /** - * The same as above. - * - * @deprecated Use the average_of_values() function instead. - */ - DEAL_II_DEPRECATED - value_type - average(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the average of the gradient $\average{\nabla u}$ on the interface * for the shape @@ -335,16 +275,6 @@ namespace FEInterfaceViews average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the average_of_gradients() function instead. - */ - DEAL_II_DEPRECATED - gradient_type - average_gradient(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the average of the Hessian $\average{\nabla^2 u} = * \frac{1}{2}\nabla^2 u_{\text{cell0}} + \frac{1}{2} \nabla^2 @@ -360,16 +290,6 @@ namespace FEInterfaceViews average_of_hessians(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the average_of_hessians() function instead. - */ - DEAL_II_DEPRECATED - hessian_type - average_hessian(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** @} */ /** @@ -830,16 +750,6 @@ namespace FEInterfaceViews jump_in_values(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_values() function instead. - */ - DEAL_II_DEPRECATED - value_type - jump(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the jump of the gradient (a tensor of rank 2) $\jump{\nabla * \mathbf{u}}$ on the interface for the shape @@ -899,16 +809,6 @@ namespace FEInterfaceViews jump_in_third_derivatives(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_third_derivatives() function instead. - */ - DEAL_II_DEPRECATED - third_derivative_type - jump_3rd_derivative(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** @} */ /** @@ -929,16 +829,6 @@ namespace FEInterfaceViews average_of_values(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the average_of_values() function instead. - */ - DEAL_II_DEPRECATED - value_type - average(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the average of the gradient (a tensor of rank 2) $\average{\nabla * \mathbf{u}}$ on the interface for the shape @@ -952,16 +842,6 @@ namespace FEInterfaceViews average_of_gradients(const unsigned int interface_dof_index, const unsigned int q_point) const; - /** - * The same as above. - * - * @deprecated Use the average_of_gradients() function instead. - */ - DEAL_II_DEPRECATED - gradient_type - average_gradient(const unsigned int interface_dof_index, - const unsigned int q_point) const; - /** * Return the average of the Hessian $\average{\nabla^2 u} = * \frac{1}{2}\nabla^2 u_{\text{cell0}} + \frac{1}{2} \nabla^2 @@ -1903,17 +1783,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_shape_values() function instead. - */ - DEAL_II_DEPRECATED - double - jump(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * Return the jump in the gradient $\jump{\nabla u}=\nabla u_{\text{cell0}} - * \nabla u_{\text{cell1}}$ on the interface for the shape function @p @@ -1932,17 +1801,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_shape_gradients() function instead. - */ - DEAL_II_DEPRECATED - Tensor<1, spacedim> - jump_gradient(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * Return the jump in the Hessian $\jump{\nabla^2 u} = \nabla^2 * u_{\text{cell0}} - \nabla^2 u_{\text{cell1}}$ on the interface for the @@ -1962,17 +1820,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_shape_hessians() function instead. - */ - DEAL_II_DEPRECATED - Tensor<2, spacedim> - jump_hessian(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * Return the jump in the third derivative $\jump{\nabla^3 u} = \nabla^3 * u_{\text{cell0}} - \nabla^3 u_{\text{cell1}}$ on the interface for the @@ -1991,17 +1838,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the jump_in_shape_3rd_derivatives() function instead. - */ - DEAL_II_DEPRECATED - Tensor<3, spacedim> - jump_3rd_derivative(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * @} */ @@ -2029,17 +1865,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the average_of_shape_values() function instead. - */ - DEAL_II_DEPRECATED - double - average(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * Return the average of the gradient $\average{\nabla u} = \frac{1}{2}\nabla * u_{\text{cell0}} + \frac{1}{2} \nabla u_{\text{cell1}}$ on the interface @@ -2058,17 +1883,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the average_of_shape_gradients() function instead. - */ - DEAL_II_DEPRECATED - Tensor<1, spacedim> - average_gradient(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * Return the average of the Hessian $\average{\nabla^2 u} = * \frac{1}{2}\nabla^2 u_{\text{cell0}} + \frac{1}{2} \nabla^2 @@ -2088,17 +1902,6 @@ public: const unsigned int q_point, const unsigned int component = 0) const; - /** - * The same as above. - * - * @deprecated Use the average_of_shape_hessians() function instead. - */ - DEAL_II_DEPRECATED - Tensor<2, spacedim> - average_hessian(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component = 0) const; - /** * @} */ @@ -3043,17 +2846,6 @@ FEInterfaceValues::jump_in_shape_values( -template -double -FEInterfaceValues::jump(const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return jump_in_shape_values(interface_dof_index, q_point, component); -} - - - template double FEInterfaceValues::average_of_shape_values( @@ -3084,18 +2876,6 @@ FEInterfaceValues::average_of_shape_values( -template -double -FEInterfaceValues::average( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return average_of_shape_values(interface_dof_index, q_point, component); -} - - - template Tensor<1, spacedim> FEInterfaceValues::average_of_shape_gradients( @@ -3126,18 +2906,6 @@ FEInterfaceValues::average_of_shape_gradients( -template -Tensor<1, spacedim> -FEInterfaceValues::average_gradient( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return average_of_shape_gradients(interface_dof_index, q_point, component); -} - - - template Tensor<2, spacedim> FEInterfaceValues::average_of_shape_hessians( @@ -3168,18 +2936,6 @@ FEInterfaceValues::average_of_shape_hessians( -template -Tensor<2, spacedim> -FEInterfaceValues::average_hessian( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return average_of_shape_hessians(interface_dof_index, q_point, component); -} - - - template Tensor<1, spacedim> FEInterfaceValues::jump_in_shape_gradients( @@ -3210,18 +2966,6 @@ FEInterfaceValues::jump_in_shape_gradients( -template -Tensor<1, spacedim> -FEInterfaceValues::jump_gradient( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return jump_in_shape_gradients(interface_dof_index, q_point, component); -} - - - template Tensor<2, spacedim> FEInterfaceValues::jump_in_shape_hessians( @@ -3252,18 +2996,6 @@ FEInterfaceValues::jump_in_shape_hessians( -template -Tensor<2, spacedim> -FEInterfaceValues::jump_hessian( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return jump_in_shape_hessians(interface_dof_index, q_point, component); -} - - - template Tensor<3, spacedim> FEInterfaceValues::jump_in_shape_3rd_derivatives( @@ -3294,18 +3026,6 @@ FEInterfaceValues::jump_in_shape_3rd_derivatives( -template -Tensor<3, spacedim> -FEInterfaceValues::jump_3rd_derivative( - const unsigned int interface_dof_index, - const unsigned int q_point, - const unsigned int component) const -{ - return jump_in_shape_3rd_derivatives(interface_dof_index, q_point, component); -} - - - template template void @@ -3540,16 +3260,6 @@ namespace FEInterfaceViews - template - typename Scalar::value_type - Scalar::jump(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_values(interface_dof_index, q_point); - } - - - template typename Scalar::value_type Scalar::average_of_values( @@ -3580,16 +3290,6 @@ namespace FEInterfaceViews - template - typename Scalar::value_type - Scalar::average(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return average_of_values(interface_dof_index, q_point); - } - - - template typename Scalar::gradient_type Scalar::average_of_gradients( @@ -3618,16 +3318,6 @@ namespace FEInterfaceViews } - template - typename Scalar::gradient_type - Scalar::average_gradient( - const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return average_of_gradients(interface_dof_index, q_point); - } - - template typename Scalar::gradient_type @@ -3658,16 +3348,6 @@ namespace FEInterfaceViews - template - typename Scalar::gradient_type - Scalar::jump_gradient(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_gradients(interface_dof_index, q_point); - } - - - template typename Scalar::hessian_type Scalar::average_of_hessians( @@ -3697,16 +3377,6 @@ namespace FEInterfaceViews - template - typename Scalar::hessian_type - Scalar::average_hessian(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return average_of_hessians(interface_dof_index, q_point); - } - - - template typename Scalar::third_derivative_type Scalar::jump_in_third_derivatives( @@ -3735,17 +3405,6 @@ namespace FEInterfaceViews - template - typename Scalar::third_derivative_type - Scalar::jump_3rd_derivative( - const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_third_derivatives(interface_dof_index, q_point); - } - - - template typename Scalar::hessian_type Scalar::jump_in_hessians( @@ -3775,16 +3434,6 @@ namespace FEInterfaceViews - template - typename Scalar::hessian_type - Scalar::jump_hessian(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_hessians(interface_dof_index, q_point); - } - - - template template void @@ -4176,16 +3825,6 @@ namespace FEInterfaceViews - template - typename Vector::value_type - Vector::jump(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_values(interface_dof_index, q_point); - } - - - template typename Vector::value_type Vector::average_of_values( @@ -4216,16 +3855,6 @@ namespace FEInterfaceViews - template - typename Vector::value_type - Vector::average(const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return average_of_values(interface_dof_index, q_point); - } - - - template typename Vector::gradient_type Vector::average_of_gradients( @@ -4255,17 +3884,6 @@ namespace FEInterfaceViews - template - typename Vector::gradient_type - Vector::average_gradient( - const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return average_of_gradients(interface_dof_index, q_point); - } - - - template typename Vector::gradient_type Vector::jump_in_gradients( @@ -4411,17 +4029,6 @@ namespace FEInterfaceViews - template - typename Vector::third_derivative_type - Vector::jump_3rd_derivative( - const unsigned int interface_dof_index, - const unsigned int q_point) const - { - return jump_in_third_derivatives(interface_dof_index, q_point); - } - - - template template void diff --git a/tests/feinterface/stokes.cc b/tests/feinterface/stokes.cc index 883c28e279..8ba022a673 100644 --- a/tests/feinterface/stokes.cc +++ b/tests/feinterface/stokes.cc @@ -665,28 +665,28 @@ namespace StokesTests ( // - nu {\nabla u}n . [v] (consistency) -nu * - (fe_fv[velocities].average_gradient(j, point) * + (fe_fv[velocities].average_of_gradients(j, point) * normals[point]) * - fe_fv[velocities].jump(i, point) + fe_fv[velocities].jump_in_values(i, point) // - nu [u] . {\nabla v}n (symmetry) // NIPG: use + - - nu * fe_fv[velocities].jump(j, point) * - (fe_fv[velocities].average_gradient(i, point) * + - nu * fe_fv[velocities].jump_in_values(j, point) * + (fe_fv[velocities].average_of_gradients(i, point) * normals[point]) // nu sigma [u].[v] (penalty) + nu * penalty * - scalar_product(fe_fv[velocities].jump(j, point), - fe_fv[velocities].jump(i, point)) + scalar_product(fe_fv[velocities].jump_in_values(j, point), + fe_fv[velocities].jump_in_values(i, point)) // {p} ([v].n) - + fe_fv[pressure].average(j, point) * - scalar_product(fe_fv[velocities].jump(i, point), + + fe_fv[pressure].average_of_values(j, point) * + scalar_product(fe_fv[velocities].jump_in_values(i, point), normals[point]) // {q} ([u].n) - + fe_fv[pressure].average(i, point) * - scalar_product(fe_fv[velocities].jump(j, point), + + fe_fv[pressure].average_of_values(i, point) * + scalar_product(fe_fv[velocities].jump_in_values(j, point), normals[point])) * JxW[point]; } diff --git a/tests/meshworker/scratch_data_08.cc b/tests/meshworker/scratch_data_08.cc index 6ee4288b10..09b17d7577 100644 --- a/tests/meshworker/scratch_data_08.cc +++ b/tests/meshworker/scratch_data_08.cc @@ -209,12 +209,13 @@ test() for (unsigned int i = 0; i < n_dofs; ++i) for (unsigned int j = 0; j < n_dofs; ++j) { - face_matrix(i, j) += - (-fev.jump_in_shape_gradients(i, q) * n[q] * fev.average(j, q) - - fev.average(i, q) * fev.jump_in_shape_gradients(j, q) * n[q] + - gh * fev.jump_in_shape_values(i, q) * - fev.jump_in_shape_values(j, q)) * - JxW[q]; + face_matrix(i, j) += (-fev.jump_in_shape_gradients(i, q) * n[q] * + fev.average_of_shape_values(j, q) - + fev.average_of_shape_values(i, q) * + fev.jump_in_shape_gradients(j, q) * n[q] + + gh * fev.jump_in_shape_values(i, q) * + fev.jump_in_shape_values(j, q)) * + JxW[q]; } };