-
Notifications
You must be signed in to change notification settings - Fork 93
/
Copy pathcommon_values.jl
312 lines (254 loc) · 12.9 KB
/
common_values.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
# Common methods for all `AbstractValues` objects
using Base: @propagate_inbounds
@noinline throw_detJ_not_pos(detJ) = throw(ArgumentError("det(J) is not positive: det(J) = $(detJ)"))
function checkquadpoint(fe_v::AbstractValues, qp::Int)
0 < qp <= getnquadpoints(fe_v) || error("quadrature point out of range")
return nothing
end
@noinline function throw_incompatible_dof_length(length_ue, n_base_funcs)
throw(ArgumentError(
"the number of base functions ($(n_base_funcs)) does not match the length " *
"of the vector ($(length_ue)). Perhaps you passed the global vector, " *
"or forgot to pass a dof_range?"
))
end
@noinline function throw_incompatible_coord_length(length_x, n_base_funcs)
throw(ArgumentError(
"the number of (geometric) base functions ($(n_base_funcs)) does not match " *
"the number of coordinates in the vector ($(length_x)). Perhaps you forgot to " *
"use an appropriate geometric interpolation when creating FE values? See " *
"https://github.com/Ferrite-FEM/Ferrite.jl/issues/265 for more details."
))
end
"""
reinit!(cv::CellValues, cell::AbstractCell, x::Vector)
reinit!(cv::CellValues, x::Vector)
reinit!(fv::FaceValues, cell::AbstractCell, x::Vector, face::Int)
reinit!(fv::FaceValues, x::Vector, face::Int)
Update the `CellValues`/`FaceValues` object for a cell or face with coordinates `x`.
The derivatives of the shape functions, and the new integration weights are computed.
For interpolations with non-identity mappings, the current `cell` is also required.
"""
reinit!
"""
getnquadpoints(fe_v::AbstractValues)
Return the number of quadrature points. For `FaceValues`,
this is the number for the current face.
"""
function getnquadpoints end
"""
getdetJdV(fe_v::AbstractValues, q_point::Int)
Return the product between the determinant of the Jacobian and the quadrature
point weight for the given quadrature point: ``\\det(J(\\mathbf{x})) w_q``.
This value is typically used when one integrates a function on a
finite element cell or face as
``\\int\\limits_\\Omega f(\\mathbf{x}) d \\Omega \\approx \\sum\\limits_{q = 1}^{n_q} f(\\mathbf{x}_q) \\det(J(\\mathbf{x})) w_q``
``\\int\\limits_\\Gamma f(\\mathbf{x}) d \\Gamma \\approx \\sum\\limits_{q = 1}^{n_q} f(\\mathbf{x}_q) \\det(J(\\mathbf{x})) w_q``
"""
function getdetJdV end
"""
shape_value(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the value of shape function `base_function` evaluated in
quadrature point `q_point`.
"""
shape_value(fe_v::AbstractValues, q_point::Int, base_function::Int)
"""
geometric_value(fe_v::AbstractValues, q_point, base_function::Int)
Return the value of the geometric shape function `base_function` evaluated in
quadrature point `q_point`.
"""
geometric_value(fe_v::AbstractValues, q_point::Int, base_function::Int)
"""
shape_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the gradient of shape function `base_function` evaluated in
quadrature point `q_point`.
"""
shape_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)
"""
shape_symmetric_gradient(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the symmetric gradient of shape function `base_function` evaluated in
quadrature point `q_point`.
"""
function shape_symmetric_gradient end
"""
shape_divergence(fe_v::AbstractValues, q_point::Int, base_function::Int)
Return the divergence of shape function `base_function` evaluated in
quadrature point `q_point`.
"""
@propagate_inbounds function shape_divergence(cv::AbstractValues, q_point::Int, base_func::Int)
return divergence_from_gradient(shape_gradient(cv, q_point, base_func))
end
divergence_from_gradient(grad::Vec) = sum(grad)
divergence_from_gradient(grad::Tensor{2}) = tr(grad)
function shape_curl(cv::AbstractValues, q_point::Int, base_func::Int)
return curl_from_gradient(shape_gradient(cv, q_point, base_func))
end
curl_from_gradient(∇v::SecondOrderTensor{3}) = Vec{3}((∇v[3,2] - ∇v[2,3], ∇v[1,3] - ∇v[3,1], ∇v[2,1] - ∇v[1,2]))
curl_from_gradient(∇v::SecondOrderTensor{2}) = Vec{1}((∇v[2,1] - ∇v[1,2],)) # Alternatively define as Vec{3}((0,0,v))
"""
function_value(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])
Compute the value of the function in a quadrature point. `u` is a vector with values
for the degrees of freedom. For a scalar valued function, `u` contains scalars.
For a vector valued function, `u` can be a vector of scalars (for use of `VectorValues`)
or `u` can be a vector of `Vec`s (for use with ScalarValues).
The value of a scalar valued function is computed as ``u(\\mathbf{x}) = \\sum\\limits_{i = 1}^n N_i (\\mathbf{x}) u_i``
where ``u_i`` are the value of ``u`` in the nodes. For a vector valued function the value is calculated as
``\\mathbf{u}(\\mathbf{x}) = \\sum\\limits_{i = 1}^n N_i (\\mathbf{x}) \\mathbf{u}_i`` where ``\\mathbf{u}_i`` are the
nodal values of ``\\mathbf{u}``.
"""
function function_value(fe_v::AbstractValues, q_point::Int, u::AbstractVector, dof_range = eachindex(u))
n_base_funcs = getnbasefunctions(fe_v)
length(dof_range) == n_base_funcs || throw_incompatible_dof_length(length(dof_range), n_base_funcs)
@boundscheck checkbounds(u, dof_range)
@boundscheck checkquadpoint(fe_v, q_point)
val = function_value_init(fe_v, u)
@inbounds for (i, j) in pairs(dof_range)
val += shape_value(fe_v, q_point, i) * u[j]
end
return val
end
"""
shape_value_type(fe_v::AbstractValues)
Return the type of `shape_value(fe_v, q_point, base_function)`
"""
function shape_value_type(fe_v::AbstractValues)
# Default fallback
return typeof(shape_value(fe_v, 1, 1))
end
function_value_init(cv::AbstractValues, ::AbstractVector{T}) where {T} = zero(shape_value_type(cv)) * zero(T)
"""
function_gradient(fe_v::AbstractValues{dim}, q_point::Int, u::AbstractVector, [dof_range])
Compute the gradient of the function in a quadrature point. `u` is a vector with values
for the degrees of freedom. For a scalar valued function, `u` contains scalars.
For a vector valued function, `u` can be a vector of scalars (for use of `VectorValues`)
or `u` can be a vector of `Vec`s (for use with ScalarValues).
The gradient of a scalar function or a vector valued function with use of `VectorValues` is computed as
``\\mathbf{\\nabla} u(\\mathbf{x}) = \\sum\\limits_{i = 1}^n \\mathbf{\\nabla} N_i (\\mathbf{x}) u_i`` or
``\\mathbf{\\nabla} u(\\mathbf{x}) = \\sum\\limits_{i = 1}^n \\mathbf{\\nabla} \\mathbf{N}_i (\\mathbf{x}) u_i`` respectively,
where ``u_i`` are the nodal values of the function.
For a vector valued function with use of `ScalarValues` the gradient is computed as
``\\mathbf{\\nabla} \\mathbf{u}(\\mathbf{x}) = \\sum\\limits_{i = 1}^n \\mathbf{u}_i \\otimes \\mathbf{\\nabla} N_i (\\mathbf{x})``
where ``\\mathbf{u}_i`` are the nodal values of ``\\mathbf{u}``.
"""
function function_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector, dof_range = eachindex(u))
n_base_funcs = getnbasefunctions(fe_v)
length(dof_range) == n_base_funcs || throw_incompatible_dof_length(length(dof_range), n_base_funcs)
@boundscheck checkbounds(u, dof_range)
@boundscheck checkquadpoint(fe_v, q_point)
grad = function_gradient_init(fe_v, u)
@inbounds for (i, j) in pairs(dof_range)
grad += shape_gradient(fe_v, q_point, i) * u[j]
end
return grad
end
# TODO: Deprecate this, nobody is using this in practice...
function function_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector{<:Vec})
n_base_funcs = getnbasefunctions(fe_v)
length(u) == n_base_funcs || throw_incompatible_dof_length(length(u), n_base_funcs)
@boundscheck checkquadpoint(fe_v, q_point)
grad = function_gradient_init(fe_v, u)
@inbounds for i in 1:n_base_funcs
grad += u[i] ⊗ shape_gradient(fe_v, q_point, i)
end
return grad
end
"""
shape_gradient_type(fe_v::AbstractValues)
Return the type of `shape_gradient(fe_v, q_point, base_function)`
"""
function shape_gradient_type(fe_v::AbstractValues)
# Default fallback
return typeof(shape_gradient(fe_v, 1, 1))
end
function function_gradient_init(cv::AbstractValues, ::AbstractVector{T}) where {T}
return zero(shape_gradient_type(cv)) * zero(T)
end
function function_gradient_init(cv::AbstractValues, ::AbstractVector{T}) where {T <: AbstractVector}
return zero(T) ⊗ zero(shape_gradient_type(cv))
end
"""
function_symmetric_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])
Compute the symmetric gradient of the function, see [`function_gradient`](@ref).
Return a `SymmetricTensor`.
The symmetric gradient of a scalar function is computed as
``\\left[ \\mathbf{\\nabla} \\mathbf{u}(\\mathbf{x_q}) \\right]^\\text{sym} = \\sum\\limits_{i = 1}^n \\frac{1}{2} \\left[ \\mathbf{\\nabla} N_i (\\mathbf{x}_q) \\otimes \\mathbf{u}_i + \\mathbf{u}_i \\otimes \\mathbf{\\nabla} N_i (\\mathbf{x}_q) \\right]``
where ``\\mathbf{u}_i`` are the nodal values of the function.
"""
function function_symmetric_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector, dof_range)
grad = function_gradient(fe_v, q_point, u, dof_range)
return symmetric(grad)
end
function function_symmetric_gradient(fe_v::AbstractValues, q_point::Int, u::AbstractVector)
grad = function_gradient(fe_v, q_point, u)
return symmetric(grad)
end
"""
function_divergence(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])
Compute the divergence of the vector valued function in a quadrature point.
The divergence of a vector valued functions in the quadrature point ``\\mathbf{x}_q)`` is computed as
``\\mathbf{\\nabla} \\cdot \\mathbf{u}(\\mathbf{x_q}) = \\sum\\limits_{i = 1}^n \\mathbf{\\nabla} N_i (\\mathbf{x_q}) \\cdot \\mathbf{u}_i``
where ``\\mathbf{u}_i`` are the nodal values of the function.
"""
function function_divergence(fe_v::AbstractValues, q_point::Int, u::AbstractVector, dof_range = eachindex(u))
return divergence_from_gradient(function_gradient(fe_v, q_point, u, dof_range))
end
# TODO: Deprecate this, nobody is using this in practice...
function function_divergence(fe_v::AbstractValues, q_point::Int, u::AbstractVector{<:Vec})
n_base_funcs = getnbasefunctions(fe_v)
length(u) == n_base_funcs || throw_incompatible_dof_length(length(u), n_base_funcs)
@boundscheck checkquadpoint(fe_v, q_point)
diverg = zero(eltype(eltype(u)))
@inbounds for i in 1:n_base_funcs
diverg += shape_gradient(fe_v, q_point, i) ⋅ u[i]
end
return diverg
end
"""
function_curl(fe_v::AbstractValues, q_point::Int, u::AbstractVector, [dof_range])
Compute the curl of the vector valued function in a quadrature point.
The curl of a vector valued functions in the quadrature point ``\\mathbf{x}_q)`` is computed as
``\\mathbf{\\nabla} \\times \\mathbf{u}(\\mathbf{x_q}) = \\sum\\limits_{i = 1}^n \\mathbf{\\nabla} N_i \\times (\\mathbf{x_q}) \\cdot \\mathbf{u}_i``
where ``\\mathbf{u}_i`` are the nodal values of the function.
"""
function_curl(fe_v::AbstractValues, q_point::Int, u::AbstractVector, dof_range = eachindex(u)) =
curl_from_gradient(function_gradient(fe_v, q_point, u, dof_range))
# TODO: Deprecate this, nobody is using this in practice...
function_curl(fe_v::AbstractValues, q_point::Int, u::AbstractVector{<:Vec}) =
curl_from_gradient(function_gradient(fe_v, q_point, u))
"""
spatial_coordinate(fe_v::AbstractValues, q_point::Int, x::AbstractVector)
Compute the spatial coordinate in a quadrature point. `x` contains the nodal
coordinates of the cell.
The coordinate is computed, using the geometric interpolation, as
``\\mathbf{x} = \\sum\\limits_{i = 1}^n M_i (\\mathbf{x}) \\mathbf{\\hat{x}}_i``
"""
function spatial_coordinate(fe_v::AbstractValues, q_point::Int, x::AbstractVector{<:Vec})
n_base_funcs = getngeobasefunctions(fe_v)
length(x) == n_base_funcs || throw_incompatible_coord_length(length(x), n_base_funcs)
@boundscheck checkquadpoint(fe_v, q_point)
vec = zero(eltype(x))
@inbounds for i in 1:n_base_funcs
vec += geometric_value(fe_v, q_point, i) * x[i]
end
return vec
end
# Utility functions used by GeometryMapping, FunctionValues
_copy_or_nothing(x) = copy(x)
_copy_or_nothing(::Nothing) = nothing
function shape_values!(values::AbstractMatrix, ip, qr_points::Vector{<:Vec})
for (qp, ξ) in pairs(qr_points)
shape_values!(@view(values[:, qp]), ip, ξ)
end
end
function shape_gradients_and_values!(gradients::AbstractMatrix, values::AbstractMatrix, ip, qr_points::Vector{<:Vec})
for (qp, ξ) in pairs(qr_points)
shape_gradients_and_values!(@view(gradients[:, qp]), @view(values[:, qp]), ip, ξ)
end
end
#= PR798
function shape_hessians_gradients_and_values!(hessians::AbstractMatrix, gradients::AbstractMatrix, values::AbstractMatrix, ip, qr_points::Vector{<:Vec})
for (qp, ξ) in pairs(qr_points)
shape_hessians_gradients_and_values!(@view(hessians[:, qp]), @view(gradients[:, qp]), @view(values[:, qp]), ip, ξ)
end
end
=#