-
Notifications
You must be signed in to change notification settings - Fork 85
/
cell_values.jl
170 lines (144 loc) · 7.12 KB
/
cell_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
# Defines CellScalarValues and CellVectorValues and common methods
"""
CellScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
CellVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
A `CellValues` 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. in the finite element cell. There are
two different types of `CellValues`: `CellScalarValues` and `CellVectorValues`. As the names suggest, `CellScalarValues`
utilizes scalar shape functions and `CellVectorValues` utilizes vectorial shape functions. For a scalar field, the
`CellScalarValues` type should be used. For vector field, both subtypes can be used.
**Arguments:**
* `T`: an optional argument (default to `Float64`) 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 a [`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)
"""
CellValues, CellScalarValues, CellVectorValues
# CellScalarValues
struct CellScalarValues{dim,T<:Real,refshape<:AbstractRefShape} <: CellValues{dim,T,refshape}
N::Matrix{T}
dNdx::Matrix{Vec{dim,T}}
dNdξ::Matrix{Vec{dim,T}}
detJdV::Vector{T}
M::Matrix{T}
dMdξ::Matrix{Vec{dim,T}}
qr::QuadratureRule{dim,refshape,T}
# 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 CellScalarValues(quad_rule::QuadratureRule, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol)
CellScalarValues(Float64, quad_rule, func_interpol, geom_interpol)
end
function CellScalarValues(::Type{T}, quad_rule::QuadratureRule{dim,shape}, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol) where {dim,T,shape<:AbstractRefShape}
@assert getdim(func_interpol) == getdim(geom_interpol)
@assert getrefshape(func_interpol) == getrefshape(geom_interpol) == shape
n_qpoints = length(getweights(quad_rule))
# Function interpolation
n_func_basefuncs = getnbasefunctions(func_interpol)
N = fill(zero(T) * T(NaN), n_func_basefuncs, n_qpoints)
dNdx = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints)
dNdξ = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints)
# Geometry interpolation
n_geom_basefuncs = getnbasefunctions(geom_interpol)
M = fill(zero(T) * T(NaN), n_geom_basefuncs, n_qpoints)
dMdξ = fill(zero(Vec{dim,T}) * T(NaN), n_geom_basefuncs, n_qpoints)
for (qp, ξ) in enumerate(quad_rule.points)
for i in 1:n_func_basefuncs
dNdξ[i, qp], N[i, qp] = gradient(ξ -> value(func_interpol, i, ξ), ξ, :all)
end
for i in 1:n_geom_basefuncs
dMdξ[i, qp], M[i, qp] = gradient(ξ -> value(geom_interpol, i, ξ), ξ, :all)
end
end
detJdV = fill(T(NaN), n_qpoints)
CellScalarValues{dim,T,shape}(N, dNdx, dNdξ, detJdV, M, dMdξ, quad_rule, func_interpol, geom_interpol)
end
# CellVectorValues
struct CellVectorValues{dim,T<:Real,refshape<:AbstractRefShape,M} <: CellValues{dim,T,refshape}
N::Matrix{Vec{dim,T}}
dNdx::Matrix{Tensor{2,dim,T,M}}
dNdξ::Matrix{Tensor{2,dim,T,M}}
detJdV::Vector{T}
M::Matrix{T}
dMdξ::Matrix{Vec{dim,T}}
qr::QuadratureRule{dim,refshape,T}
# 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 CellVectorValues(quad_rule::QuadratureRule, func_interpol::Interpolation, geom_interpol::Interpolation=func_interpol)
CellVectorValues(Float64, quad_rule, func_interpol, geom_interpol)
end
function CellVectorValues(::Type{T}, quad_rule::QuadratureRule{dim,shape}, func_interpol::Interpolation,
geom_interpol::Interpolation=func_interpol) where {dim,T,shape<:AbstractRefShape}
@assert getdim(func_interpol) == getdim(geom_interpol)
@assert getrefshape(func_interpol) == getrefshape(geom_interpol) == shape
n_qpoints = length(getweights(quad_rule))
# Function interpolation
n_func_basefuncs = getnbasefunctions(func_interpol) * dim
N = fill(zero(Vec{dim,T}) * T(NaN), n_func_basefuncs, n_qpoints)
dNdx = fill(zero(Tensor{2,dim,T}) * T(NaN), n_func_basefuncs, n_qpoints)
dNdξ = fill(zero(Tensor{2,dim,T}) * T(NaN), n_func_basefuncs, n_qpoints)
# Geometry interpolation
n_geom_basefuncs = getnbasefunctions(geom_interpol)
M = fill(zero(T) * T(NaN), n_geom_basefuncs, n_qpoints)
dMdξ = fill(zero(Vec{dim,T}) * T(NaN), n_geom_basefuncs, n_qpoints)
for (qp, ξ) in enumerate(quad_rule.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] = Vec{dim,T}((N_comp...,))
dN_comp = zeros(T, dim, dim)
dN_comp[comp, :] = dNdξ_temp
dNdξ[basefunc_count, qp] = Tensor{2,dim,T}((dN_comp...,))
basefunc_count += 1
end
end
for basefunc in 1:n_geom_basefuncs
dMdξ[basefunc, qp], M[basefunc, qp] = gradient(ξ -> value(geom_interpol, basefunc, ξ), ξ, :all)
end
end
detJdV = fill(T(NaN), n_qpoints)
MM = Tensors.n_components(Tensors.get_base(eltype(dNdx)))
CellVectorValues{dim,T,shape,MM}(N, dNdx, dNdξ, detJdV, M, dMdξ, quad_rule, func_interpol, geom_interpol)
end
function reinit!(cv::CellValues{dim}, x::AbstractVector{Vec{dim,T}}) where {dim,T}
n_geom_basefuncs = getngeobasefunctions(cv)
n_func_basefuncs = getnbasefunctions(cv)
length(x) == n_geom_basefuncs || throw_incompatible_coord_length(length(x), n_geom_basefuncs)
@inbounds for i in 1:length(cv.qr.weights)
w = cv.qr.weights[i]
fecv_J = zero(Tensor{2,dim})
for j in 1:n_geom_basefuncs
fecv_J += x[j] ⊗ cv.dMdξ[j, i]
end
detJ = det(fecv_J)
detJ > 0.0 || throw_detJ_not_pos(detJ)
cv.detJdV[i] = detJ * w
Jinv = inv(fecv_J)
for j in 1:n_func_basefuncs
cv.dNdx[j, i] = cv.dNdξ[j, i] ⋅ Jinv
end
end
end