-
Notifications
You must be signed in to change notification settings - Fork 100
/
Copy pathMonomialBases.jl
420 lines (349 loc) · 9.07 KB
/
MonomialBases.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
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
"""
struct MonomialBasis{D,T} <: Field
Type representing a basis of multivariate scalar-valued, vector-valued, or
tensor-valued, iso- or aniso-tropic monomials. The fields
of this `struct` are not public
This type fully implements the [`Field`](@ref) interface, with up to second order
derivatives.
"""
struct MonomialBasis{D,T} <: Field
orders::NTuple{D,Int}
terms::Vector{CartesianIndex{D}}
function MonomialBasis{D}(
::Type{T}, orders::NTuple{D,Int}, terms::Vector{CartesianIndex{D}}) where {D,T}
new{D,T}(orders,terms)
end
end
"""
MonomialBasis{D}(::Type{T}, orders::Tuple [, filter::Function]) where {D,T}
This version of the constructor allows to pass a tuple `orders` containing the
polynomial order to be used in each of the `D` dimensions in order to construct
and anisotropic tensor-product space.
"""
function MonomialBasis{D}(
::Type{T}, orders::NTuple{D,Int}, filter::Function=_q_filter) where {D,T}
terms = _define_terms(filter, orders)
MonomialBasis{D}(T,orders,terms)
end
"""
MonomialBasis{D}(::Type{T}, order::Int [, filter::Function]) where {D,T}
Returns an instance of `MonomialBasis` representing a multivariate polynomial basis
in `D` dimensions, of polynomial degree `order`, whose value is represented by the type `T`.
The type `T` is typically `<:Number`, e.g., `Float64` for scalar-valued functions and `VectorValue{D,Float64}`
for vector-valued ones.
# Filter function
The `filter` function is used to select which terms of the tensor product space
of order `order` in `D` dimensions are to be used. If the filter is not provided, the full tensor-product
space is used by default leading to a multivariate polynomial space of type Q.
The signature of the filter function is
(e,order) -> Bool
where `e` is a tuple of `D` integers containing the exponents of a multivariate monomial. The following filters
are used to select well known polynomial spaces
- Q space: `(e,order) -> true`
- P space: `(e,order) -> sum(e) <= order`
- "Serendipity" space: `(e,order) -> sum( [ i for i in e if i>1 ] ) <= order`
"""
function MonomialBasis{D}(
::Type{T}, order::Int, filter::Function=_q_filter) where {D,T}
orders = tfill(order,Val{D}())
MonomialBasis{D}(T,orders,filter)
end
# API
"""
get_exponents(b::MonomialBasis)
Get a vector of tuples with the exponents of all the terms in the
monomial basis.
# Examples
```jldoctest
using Gridap.Polynomials
b = MonomialBasis{2}(Float64,2)
exponents = get_exponents(b)
println(exponents)
# output
Tuple{Int,Int}[(0, 0), (1, 0), (2, 0), (0, 1), (1, 1), (2, 1), (0, 2), (1, 2), (2, 2)]
```
"""
function get_exponents(b::MonomialBasis)
indexbase = 1
[Tuple(t) .- indexbase for t in b.terms]
end
"""
get_order(b::MonomialBasis)
"""
function get_order(b::MonomialBasis)
maximum(b.orders)
end
"""
get_orders(b::MonomialBasis)
"""
function get_orders(b::MonomialBasis)
b.orders
end
"""
"""
get_value_type(::MonomialBasis{D,T}) where {D,T} = T
get_value_type(::Type{MonomialBasis{D,T}}) where {D,T} = T
# Field implementation
function field_cache(f::MonomialBasis{D,T},x) where {D,T}
@assert D == length(eltype(x)) "Incorrect number of point components"
np = length(x)
ndof = length(f.terms)*num_components(T)
n = 1 + _maximum(f.orders)
r = CachedArray(zeros(T,(np,ndof)))
v = CachedArray(zeros(T,(ndof,)))
c = CachedArray(zeros(eltype(T),(D,n)))
(r, v, c)
end
function evaluate_field!(cache,f::MonomialBasis{D,T},x) where {D,T}
r, v, c = cache
np = length(x)
ndof = length(f.terms)*num_components(T)
n = 1 + _maximum(f.orders)
setsize!(r,(np,ndof))
setsize!(v,(ndof,))
setsize!(c,(D,n))
for i in 1:np
@inbounds xi = x[i]
_evaluate_nd!(v,xi,f.orders,f.terms,c)
for j in 1:ndof
@inbounds r[i,j] = v[j]
end
end
r.array
end
function gradient_cache(f::MonomialBasis{D,V},x) where {D,V}
@assert D == length(eltype(x)) "Incorrect number of point components"
np = length(x)
ndof = length(f.terms)*num_components(V)
xi = testitem(x)
T = gradient_type(V,xi)
n = 1 + _maximum(f.orders)
r = CachedArray(zeros(T,(np,ndof)))
v = CachedArray(zeros(T,(ndof,)))
c = CachedArray(zeros(eltype(T),(D,n)))
g = CachedArray(zeros(eltype(T),(D,n)))
(r, v, c, g)
end
function evaluate_gradient!(cache,f::MonomialBasis{D,T},x) where {D,T}
r, v, c, g = cache
np = length(x)
ndof = length(f.terms) * num_components(T)
n = 1 + _maximum(f.orders)
setsize!(r,(np,ndof))
setsize!(v,(ndof,))
setsize!(c,(D,n))
setsize!(g,(D,n))
for i in 1:np
@inbounds xi = x[i]
_gradient_nd!(v,xi,f.orders,f.terms,c,g,T)
for j in 1:ndof
@inbounds r[i,j] = v[j]
end
end
r.array
end
function hessian_cache(f::MonomialBasis{D,V},x) where {D,V}
@assert D == length(eltype(x)) "Incorrect number of point components"
np = length(x)
ndof = length(f.terms)*num_components(V)
xi = testitem(x)
T = gradient_type(gradient_type(V,xi),xi)
n = 1 + _maximum(f.orders)
r = CachedArray(zeros(T,(np,ndof)))
v = CachedArray(zeros(T,(ndof,)))
c = CachedArray(zeros(eltype(T),(D,n)))
g = CachedArray(zeros(eltype(T),(D,n)))
h = CachedArray(zeros(eltype(T),(D,n)))
(r, v, c, g, h)
end
function evaluate_hessian!(cache,f::MonomialBasis{D,T},x) where {D,T}
r, v, c, g, h = cache
np = length(x)
ndof = length(f.terms) * num_components(T)
n = 1 + _maximum(f.orders)
setsize!(r,(np,ndof))
setsize!(v,(ndof,))
setsize!(c,(D,n))
setsize!(g,(D,n))
setsize!(h,(D,n))
for i in 1:np
@inbounds xi = x[i]
_hessian_nd!(v,xi,f.orders,f.terms,c,g,h,T)
for j in 1:ndof
@inbounds r[i,j] = v[j]
end
end
r.array
end
# Helpers
_q_filter(e,o) = true
function _define_terms(filter,orders)
t = orders .+ 1
g = (0 .* orders) .+ 1
cis = CartesianIndices(t)
co = CartesianIndex(g)
maxorder = _maximum(orders)
[ ci for ci in cis if filter(Int[Tuple(ci-co)...],maxorder) ]
end
function _evaluate_1d!(v::AbstractMatrix{T},x,order,d) where T
n = order + 1
z = one(T)
@inbounds v[d,1] = z
for i in 2:n
@inbounds v[d,i] = x[d]^(i-1)
end
end
function _gradient_1d!(v::AbstractMatrix{T},x,order,d) where T
n = order + 1
z = zero(T)
@inbounds v[d,1] = z
for i in 2:n
@inbounds v[d,i] = (i-1)*x[d]^(i-2)
end
end
function _hessian_1d!(v::AbstractMatrix{T},x,order,d) where T
n = order + 1
z = zero(T)
@inbounds v[d,1] = z
if n>1
@inbounds v[d,2] = z
end
for i in 3:n
@inbounds v[d,i] = (i-1)*(i-2)*x[d]^(i-3)
end
end
function _evaluate_nd!(
v::AbstractVector{V},
x,
orders,
terms::AbstractVector{CartesianIndex{D}},
c::AbstractMatrix{T}) where {V,T,D}
dim = D
for d in 1:dim
_evaluate_1d!(c,x,orders[d],d)
end
o = one(T)
k = 1
for ci in terms
s = o
for d in 1:dim
@inbounds s *= c[d,ci[d]]
end
k = _set_value!(v,s,k)
end
end
@inline function _set_value!(v::AbstractVector{V},s::T,k) where {V,T}
m = zero(mutable(V))
z = zero(T)
js = eachindex(m)
for j in js
for i in js
@inbounds m[i] = z
end
m[j] = s
v[k] = m
k += 1
end
k
end
@inline function _set_value!(v::AbstractVector{<:Real},s,k)
@inbounds v[k] = s
k+1
end
function _gradient_nd!(
v::AbstractVector{G},
x,
orders,
terms::AbstractVector{CartesianIndex{D}},
c::AbstractMatrix{T},
g::AbstractMatrix{T},
::Type{V}) where {G,T,D,V}
dim = D
for d in 1:dim
_evaluate_1d!(c,x,orders[d],d)
_gradient_1d!(g,x,orders[d],d)
end
z = zero(mutable(VectorValue{D,T}))
o = one(T)
k = 1
for ci in terms
s = z
for i in eachindex(s)
@inbounds s[i] = o
end
for q in 1:dim
for d in 1:dim
if d != q
@inbounds s[q] *= c[d,ci[d]]
else
@inbounds s[q] *= g[d,ci[d]]
end
end
end
k = _set_gradient!(v,s,k,V)
end
end
@inline function _set_gradient!(
v::AbstractVector{G},s,k,::Type{<:Real}) where G
@inbounds v[k] = s
k+1
end
@inline function _set_gradient!(
v::AbstractVector{G},s,k,::Type{V}) where {V,G}
T = eltype(s)
m = zero(mutable(G))
w = zero(V)
z = zero(T)
for j in CartesianIndices(w)
for i in CartesianIndices(m)
@inbounds m[i] = z
end
for i in CartesianIndices(s)
@inbounds m[i,j] = s[i]
end
@inbounds v[k] = m
k += 1
end
k
end
function _hessian_nd!(
v::AbstractVector{G},
x,
orders,
terms::AbstractVector{CartesianIndex{D}},
c::AbstractMatrix{T},
g::AbstractMatrix{T},
h::AbstractMatrix{T},
::Type{V}) where {G,T,D,V}
dim = D
for d in 1:dim
_evaluate_1d!(c,x,orders[d],d)
_gradient_1d!(g,x,orders[d],d)
_hessian_1d!(h,x,orders[d],d)
end
z = zero(mutable(TensorValue{D,D,T}))
o = one(T)
k = 1
for ci in terms
s = z
for i in eachindex(s)
@inbounds s[i] = o
end
for r in 1:dim
for q in 1:dim
for d in 1:dim
if d != q && d != r
@inbounds s[r,q] *= c[d,ci[d]]
elseif d == q && d ==r
@inbounds s[r,q] *= h[d,ci[d]]
else
@inbounds s[r,q] *= g[d,ci[d]]
end
end
end
end
k = _set_gradient!(v,s,k,V)
end
end
_maximum(orders::Tuple{}) = 0
_maximum(orders) = maximum(orders)