/
face_values.jl
290 lines (239 loc) · 12 KB
/
face_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
# Defines FaceScalarValues and FaceVectorValues and common methods
"""
FaceScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
FaceVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
A `FaceValues` object facilitates the process of evaluating values of shape functions, gradients of shape functions,
values of nodal functions, gradients and divergences of nodal functions etc. on the faces of finite elements. There are
two different types of `FaceValues`: `FaceScalarValues` and `FaceVectorValues`. As the names suggest,
`FaceScalarValues` utilizes scalar shape functions and `FaceVectorValues` utilizes vectorial shape functions.
For a scalar field, the `FaceScalarValues` type should be used. For vector field, both subtypes can be used.
!!! note
The quadrature rule for the face should be given with one dimension lower.
I.e. for a 3D case, the quadrature rule should be in 2D.
**Arguments:**
* `T`: an optional argument to determine the type the internal data is stored as.
* `quad_rule`: an instance of a [`QuadratureRule`](@ref)
* `func_interpol`: an instance of an [`Interpolation`](@ref) used to interpolate the approximated function
* `geom_interpol`: an optional instance of an [`Interpolation`](@ref) which is used to interpolate the geometry
**Common methods:**
* [`reinit!`](@ref)
* [`getnquadpoints`](@ref)
* [`getdetJdV`](@ref)
* [`shape_value`](@ref)
* [`shape_gradient`](@ref)
* [`shape_symmetric_gradient`](@ref)
* [`shape_divergence`](@ref)
* [`function_value`](@ref)
* [`function_gradient`](@ref)
* [`function_symmetric_gradient`](@ref)
* [`function_divergence`](@ref)
* [`spatial_coordinate`](@ref)
"""
FaceValues, FaceScalarValues, FaceVectorValues
# FaceScalarValues
struct FaceScalarValues{dim,T<:Real,refshape<:AbstractRefShape} <: FaceValues{dim,T,refshape}
N::Array{T,3}
dNdx::Array{Vec{dim,T},3}
dNdξ::Array{Vec{dim,T},3}
detJdV::Matrix{T}
normals::Vector{Vec{dim,T}}
M::Array{T,3}
dMdξ::Array{Vec{dim,T},3}
# 'Any' is 'dim-1' here -- this is deliberately abstractly typed. Only qr.weights is
# accessed in performance critical code so this doesn't seem to be a problem in
# practice since qr.weights is correctly inferred as Vector{T}, and T is a parameter
# of the struct.
qr::QuadratureRule{<:Any,refshape,T}
current_face::ScalarWrapper{Int}
# The following fields are deliberately abstract -- they are never used in
# performance critical code, just stored here for convenience.
func_interp::Interpolation{dim,refshape}
geo_interp::Interpolation{dim,refshape}
end
function FaceScalarValues(quad_rule::QuadratureRule, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol)
FaceScalarValues(Float64, quad_rule, func_interpol, geom_interpol)
end
function FaceScalarValues(::Type{T}, quad_rule::QuadratureRule{dim_qr,shape}, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol) where {dim_qr,T,shape<:AbstractRefShape}
@assert getdim(func_interpol) == getdim(geom_interpol)
@assert getrefshape(func_interpol) == getrefshape(geom_interpol) == shape
n_qpoints = length(getweights(quad_rule))
dim = dim_qr + 1
face_quad_rule = create_face_quad_rule(quad_rule, func_interpol)
n_faces = length(face_quad_rule)
# Normals
normals = zeros(Vec{dim,T}, n_qpoints)
# Function interpolation
n_func_basefuncs = getnbasefunctions(func_interpol)
N = fill(zero(T) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
dNdx = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
dNdξ = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
# Geometry interpolation
n_geom_basefuncs = getnbasefunctions(geom_interpol)
M = fill(zero(T) * T(NaN), n_geom_basefuncs, n_qpoints, n_faces)
dMdξ = fill(zero(Vec{dim,T}) * T(NaN), n_geom_basefuncs, n_qpoints, n_faces)
for face in 1:n_faces, (qp, ξ) in enumerate(face_quad_rule[face].points)
for i in 1:n_func_basefuncs
dNdξ[i, qp, face], N[i, qp, face] = gradient(ξ -> value(func_interpol, i, ξ), ξ, :all)
end
for i in 1:n_geom_basefuncs
dMdξ[i, qp, face], M[i, qp, face] = gradient(ξ -> value(geom_interpol, i, ξ), ξ, :all)
end
end
detJdV = fill(T(NaN), n_qpoints, n_faces)
FaceScalarValues{dim,T,shape}(N, dNdx, dNdξ, detJdV, normals, M, dMdξ, quad_rule, ScalarWrapper(0), func_interpol, geom_interpol)
end
# FaceVectorValues
struct FaceVectorValues{dim,T<:Real,refshape<:AbstractRefShape,M} <: FaceValues{dim,T,refshape}
N::Array{Vec{dim,T},3}
dNdx::Array{Tensor{2,dim,T,M},3}
dNdξ::Array{Tensor{2,dim,T,M},3}
detJdV::Matrix{T}
normals::Vector{Vec{dim,T}}
M::Array{T,3}
dMdξ::Array{Vec{dim,T},3}
# 'Any' is 'dim-1' here -- this is deliberately abstractly typed. Only qr.weights is
# accessed in performance critical code so this doesn't seem to be a problem in
# practice since qr.weights is correctly inferred as Vector{T}, and T is a parameter
# of the struct.
qr::QuadratureRule{<:Any,refshape,T}
current_face::ScalarWrapper{Int}
# The following fields are deliberately abstract -- they are never used in
# performance critical code, just stored here for convenience.
func_interp::Interpolation{dim,refshape}
geo_interp::Interpolation{dim,refshape}
end
function FaceVectorValues(quad_rule::QuadratureRule, func_interpol::Interpolation, geom_interpol::Interpolation=func_interpol)
FaceVectorValues(Float64, quad_rule, func_interpol, geom_interpol)
end
function FaceVectorValues(::Type{T}, quad_rule::QuadratureRule{dim_qr,shape}, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol) where {dim_qr,T,shape<:AbstractRefShape}
@assert getdim(func_interpol) == getdim(geom_interpol)
@assert getrefshape(func_interpol) == getrefshape(geom_interpol) == shape
n_qpoints = length(getweights(quad_rule))
dim = dim_qr + 1
face_quad_rule = create_face_quad_rule(quad_rule, func_interpol)
n_faces = length(face_quad_rule)
# Normals
normals = zeros(Vec{dim,T}, n_qpoints)
# Function interpolation
n_func_basefuncs = getnbasefunctions(func_interpol) * dim
N = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
dNdx = fill(zero(Tensor{2,dim,T}) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
dNdξ = fill(zero(Tensor{2,dim,T}) * T(NaN), n_func_basefuncs, n_qpoints, n_faces)
# Geometry interpolation
n_geom_basefuncs = getnbasefunctions(geom_interpol)
M = fill(zero(T) * T(NaN), n_geom_basefuncs, n_qpoints, n_faces)
dMdξ = fill(zero(Vec{dim,T}) * T(NaN), n_geom_basefuncs, n_qpoints, n_faces)
for face in 1:n_faces, (qp, ξ) in enumerate(face_quad_rule[face].points)
basefunc_count = 1
for basefunc in 1:getnbasefunctions(func_interpol)
dNdξ_temp, N_temp = gradient(ξ -> value(func_interpol, basefunc, ξ), ξ, :all)
for comp in 1:dim
N_comp = zeros(T, dim)
N_comp[comp] = N_temp
N[basefunc_count, qp, face] = Vec{dim,T}((N_comp...,))
dN_comp = zeros(T, dim, dim)
dN_comp[comp, :] = dNdξ_temp
dNdξ[basefunc_count, qp, face] = Tensor{2,dim,T}((dN_comp...,))
basefunc_count += 1
end
end
for basefunc in 1:n_geom_basefuncs
dMdξ[basefunc, qp, face], M[basefunc, qp, face] = gradient(ξ -> value(geom_interpol, basefunc, ξ), ξ, :all)
end
end
detJdV = fill(T(NaN), n_qpoints, n_faces)
MM = Tensors.n_components(Tensors.get_base(eltype(dNdx)))
FaceVectorValues{dim,T,shape,MM}(N, dNdx, dNdξ, detJdV, normals, M, dMdξ, quad_rule, ScalarWrapper(0), func_interpol, geom_interpol)
end
function reinit!(fv::FaceValues{dim}, x::AbstractVector{Vec{dim,T}}, face::Int) where {dim,T}
n_geom_basefuncs = getngeobasefunctions(fv)
n_func_basefuncs = getnbasefunctions(fv)
length(x) == n_geom_basefuncs || throw_incompatible_coord_length(length(x), n_geom_basefuncs)
@boundscheck checkface(fv, face)
fv.current_face[] = face
cb = getcurrentface(fv)
@inbounds for i in 1:length(fv.qr.weights)
w = fv.qr.weights[i]
fefv_J = zero(Tensor{2,dim})
for j in 1:n_geom_basefuncs
fefv_J += x[j] ⊗ fv.dMdξ[j, i, cb]
end
weight_norm = weighted_normal(fefv_J, fv, cb)
fv.normals[i] = weight_norm / norm(weight_norm)
detJ = norm(weight_norm)
detJ > 0.0 || throw_detJ_not_pos(detJ)
fv.detJdV[i, cb] = detJ * w
Jinv = inv(fefv_J)
for j in 1:n_func_basefuncs
fv.dNdx[j, i, cb] = fv.dNdξ[j, i, cb] ⋅ Jinv
end
end
end
"""
getcurrentface(fv::FaceValues)
Return the current active face of the `FaceValues` object (from last `reinit!`).
"""
getcurrentface(fv::FaceValues) = fv.current_face[]
"""
getnormal(fv::FaceValues, qp::Int)
Return the normal at the quadrature point `qp` for the active face of the
`FaceValues` object(from last `reinit!`).
"""
getnormal(fv::FaceValues, qp::Int) = fv.normals[qp]
"""
BCValues(func_interpol::Interpolation, geom_interpol::Interpolation, boundary_type::Union{Type{<:BoundaryIndex}})
`BCValues` stores the shape values at all faces/edges/vertices (depending on `boundary_type`) for the geomatric interpolation (`geom_interpol`),
for each dof-position determined by the `func_interpol`. Used mainly by the `ConstrainHandler`.
"""
struct BCValues{T}
M::Array{T,3}
nqp::Array{Int}
current_entity::ScalarWrapper{Int}
end
BCValues(func_interpol::Interpolation, geom_interpol::Interpolation, boundary_type::Type{<:BoundaryIndex} = Ferrite.FaceIndex) =
BCValues(Float64, func_interpol, geom_interpol, boundary_type)
function BCValues(::Type{T}, func_interpol::Interpolation{dim,refshape}, geom_interpol::Interpolation{dim,refshape}, boundary_type::Type{<:BoundaryIndex} = Ferrite.FaceIndex) where {T,dim,refshape}
# set up quadrature rules for each boundary entity with dof-positions
# (determined by func_interpol) as the quadrature points
interpolation_coords = reference_coordinates(func_interpol)
qrs = QuadratureRule{dim,refshape,T}[]
for boundarydofs in boundarydof_indices(boundary_type)(func_interpol)
dofcoords = Vec{dim,T}[]
for boundarydof in boundarydofs
push!(dofcoords, interpolation_coords[boundarydof])
end
qrf = QuadratureRule{dim,refshape,T}(fill(T(NaN), length(dofcoords)), dofcoords) # weights will not be used
push!(qrs, qrf)
end
n_boundary_entities = length(qrs)
n_qpoints = n_boundary_entities == 0 ? 0 : maximum(qr->length(getweights(qr)), qrs) # Bound number of qps correctly.
n_geom_basefuncs = getnbasefunctions(geom_interpol)
M = fill(zero(T) * T(NaN), n_geom_basefuncs, n_qpoints, n_boundary_entities)
nqp = zeros(Int,n_boundary_entities)
for n_boundary_entity in 1:n_boundary_entities
for (qp, ξ) in enumerate(qrs[n_boundary_entity].points), i in 1:n_geom_basefuncs
M[i, qp, n_boundary_entity] = value(geom_interpol, i, ξ)
end
nqp[n_boundary_entity] = length(qrs[n_boundary_entity].points)
end
BCValues{T}(M, nqp, ScalarWrapper(0))
end
getnquadpoints(bcv::BCValues) = bcv.nqp[bcv.current_entity.x]
function spatial_coordinate(bcv::BCValues, q_point::Int, xh::AbstractVector{Vec{dim,T}}) where {dim,T}
n_base_funcs = size(bcv.M, 1)
length(xh) == n_base_funcs || throw_incompatible_coord_length(length(xh), n_base_funcs)
x = zero(Vec{dim,T})
face = bcv.current_entity[]
@inbounds for i in 1:n_base_funcs
x += bcv.M[i,q_point,face] * xh[i] # geometric_value(fe_v, q_point, i) * xh[i]
end
return x
end
nfaces(fv) = size(fv.N, 3)
function checkface(fv::FaceValues, face::Int)
0 < face <= nfaces(fv) || error("Face index out of range.")
return nothing
end