Skip to content

Commit a4fdff4

Browse files
authored
Add shape_symmetric_gradient for PointValues (#1325)
Defines generic `shape_symmetric_gradient` for `AbstractValues`
1 parent efa551d commit a4fdff4

6 files changed

Lines changed: 13 additions & 4 deletions

File tree

CHANGELOG.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,7 @@ The format is based on [Keep a Changelog](https://keepachangelog.com/en/1.0.0/),
66
and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0.html).
77

88
## [Next] - xxxx-xx-xx
9+
- Make `shape_symmetric_gradient` work for `PointValues` ([#1325])
910

1011
## [v1.4.0] - 2026-04-20
1112

src/FEValues/CellValues.jl

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -100,7 +100,6 @@ shape_hessian_type(cv::CellValues) = shape_hessian_type(cv.fun_values)
100100
@propagate_inbounds shape_value(cv::CellValues, q_point::Int, i::Int) = shape_value(cv.fun_values, q_point, i)
101101
@propagate_inbounds shape_gradient(cv::CellValues, q_point::Int, i::Int) = shape_gradient(cv.fun_values, q_point, i)
102102
@propagate_inbounds shape_hessian(cv::CellValues, q_point::Int, i::Int) = shape_hessian(cv.fun_values, q_point, i)
103-
@propagate_inbounds shape_symmetric_gradient(cv::CellValues, q_point::Int, i::Int) = shape_symmetric_gradient(cv.fun_values, q_point, i)
104103

105104
# Access quadrature rule values
106105
get_quadrature_rule(cv::CellValues) = cv.qr

src/FEValues/FacetValues.jl

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -93,7 +93,6 @@ get_fun_values(fv::FacetValues) = @inbounds fv.fun_values[getcurrentfacet(fv)]
9393
@propagate_inbounds shape_value(fv::FacetValues, q_point::Int, i::Int) = shape_value(get_fun_values(fv), q_point, i)
9494
@propagate_inbounds shape_gradient(fv::FacetValues, q_point::Int, i::Int) = shape_gradient(get_fun_values(fv), q_point, i)
9595
@propagate_inbounds shape_hessian(fv::FacetValues, q_point::Int, i::Int) = shape_hessian(get_fun_values(fv), q_point, i)
96-
@propagate_inbounds shape_symmetric_gradient(fv::FacetValues, q_point::Int, i::Int) = shape_symmetric_gradient(get_fun_values(fv), q_point, i)
9796

9897
"""
9998
getcurrentfacet(fv::FacetValues)

src/FEValues/FunctionValues.jl

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -109,7 +109,6 @@ getnquadpoints(funvals::FunctionValues) = size(funvals.Nx, 2)
109109
@propagate_inbounds shape_value(funvals::FunctionValues, q_point::Int, base_func::Int) = funvals.Nx[base_func, q_point]
110110
@propagate_inbounds shape_gradient(funvals::FunctionValues, q_point::Int, base_func::Int) = funvals.dNdx[base_func, q_point]
111111
@propagate_inbounds shape_hessian(funvals::FunctionValues{2}, q_point::Int, base_func::Int) = funvals.d2Ndx2[base_func, q_point]
112-
@propagate_inbounds shape_symmetric_gradient(funvals::FunctionValues, q_point::Int, base_func::Int) = symmetric(shape_gradient(funvals, q_point, base_func))
113112

114113
function_interpolation(funvals::FunctionValues) = funvals.ip
115114
function_difforder(::FunctionValues{DiffOrder}) where {DiffOrder} = DiffOrder

src/FEValues/common_values.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -135,7 +135,7 @@ shape_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)
135135
Return the symmetric gradient of shape function `base_function` evaluated in
136136
quadrature point `q_point`.
137137
"""
138-
function shape_symmetric_gradient end
138+
@propagate_inbounds shape_symmetric_gradient(cv::AbstractValues, q_point::Int, base_func::Int) = symmetric(shape_gradient(cv, q_point, base_func))
139139

140140
"""
141141
shape_divergence(fe_v::AbstractValues, q_point::Int, base_function::Int)

test/test_pointevaluation.jl

Lines changed: 11 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -513,8 +513,19 @@ end
513513
@test function_value(pvv, uv) function_value(cvv, 1, uv)
514514
@test function_gradient(pvv, uv) function_gradient(cvv, 1, uv)
515515
@test function_symmetric_gradient(pvv, uv) function_symmetric_gradient(cvv, 1, uv)
516+
517+
@test shape_gradient(pvv, 1, 1) shape_gradient(cvv, 1, 1)
518+
@test shape_divergence(pvv, 1, 1) shape_divergence(cvv, 1, 1)
519+
@test shape_curl(pvv, 1, 1) shape_curl(cvv, 1, 1)
520+
@test shape_symmetric_gradient(pvv, 1, 1) shape_symmetric_gradient(cvv, 1, 1)
521+
516522
reinit!(pvv, x, ξ₂)
517523
@test function_value(pvv, uv) function_value(cvv, 2, uv)
518524
@test function_gradient(pvv, uv) function_gradient(cvv, 2, uv)
519525
@test function_symmetric_gradient(pvv, uv) function_symmetric_gradient(cvv, 2, uv)
526+
527+
@test shape_gradient(pvv, 1, 1) shape_gradient(cvv, 2, 1)
528+
@test shape_divergence(pvv, 1, 1) shape_divergence(cvv, 2, 1)
529+
@test shape_curl(pvv, 1, 1) shape_curl(cvv, 2, 1)
530+
@test shape_symmetric_gradient(pvv, 1, 1) shape_symmetric_gradient(cvv, 2, 1)
520531
end

0 commit comments

Comments
 (0)