-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathfieldvector.jl
571 lines (486 loc) · 17.2 KB
/
fieldvector.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
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
import BlockArrays
"""
FieldVector
A `FieldVector` is a wrapper around one or more `Field`s that acts like vector
of the underlying arrays.
It is similar in spirit to [`ArrayPartition` from
RecursiveArrayTools.jl](https://github.com/SciML/RecursiveArrayTools.jl#arraypartition),
but allows referring to fields by name.
# Constructors
FieldVector(;name1=field1, name2=field2, ...)
Construct a `FieldVector`, wrapping `field1, field2, ...` using the names
`name1, name2, ...`.
"""
struct FieldVector{T, M} <: BlockArrays.AbstractBlockVector{T}
values::M
end
FieldVector{T}(values::M) where {T, M} = FieldVector{T, M}(values)
function Adapt.adapt_structure(to, fv::FieldVector)
pn = propertynames(fv)
vals = map(key -> Adapt.adapt(to, getproperty(fv, key)), pn)
return FieldVector(; NamedTuple{pn}(vals)...)
end
"""
Fields.ScalarWrapper(val) <: AbstractArray{T,0}
This is a wrapper around scalar values that allows them to be mutated as part of
a FieldVector. A call `getproperty` on a `FieldVector` with this component will
return a scalar, instead of the boxed object.
"""
mutable struct ScalarWrapper{T} <: AbstractArray{T, 0}
val::T
end
Base.size(::ScalarWrapper) = ()
Base.getindex(s::ScalarWrapper) = s.val
Base.setindex!(s::ScalarWrapper, value) = s.val = value
Base.similar(s::ScalarWrapper) = ScalarWrapper(s.val)
"""
Fields.wrap(x)
Construct a mutable wrapper around `x`. This can be extended for new types
(especially immutable ones).
"""
wrap(x) = x
wrap(x::Real) = ScalarWrapper(x)
wrap(x::NamedTuple) = FieldVector(; pairs(x)...)
"""
Fields.unwrap(x::T)
This is called when calling `getproperty` on a `FieldVector` property of element
type `T`.
"""
unwrap(x) = x
unwrap(x::ScalarWrapper) = x[]
function FieldVector(; kwargs...)
values = map(wrap, NamedTuple(kwargs))
T = promote_type(
map(RecursiveArrayTools.recursive_bottom_eltype, values)...,
)
return FieldVector{T}(values)
end
_values(fv::FieldVector) = getfield(fv, :values)
"""
backing_array(x)
The `AbstractArray` that is backs an object `x`, allowing it to be treated as a
component of a `FieldVector`.
"""
backing_array(x) = x
backing_array(x::Field) = parent(x)
Base.propertynames(fv::FieldVector) = propertynames(_values(fv))
@inline function Base.getproperty(fv::FieldVector, name::Symbol)
unwrap(getfield(_values(fv), name))
end
@inline function Base.setproperty!(fv::FieldVector, name::Symbol, value)
x = getfield(_values(fv), name)
x .= value
end
BlockArrays.blockaxes(fv::FieldVector) =
(BlockArrays.BlockRange(1:length(_values(fv))),)
Base.axes(fv::FieldVector) =
(BlockArrays.blockedrange(map(length ∘ backing_array, Tuple(_values(fv)))),)
Base.@propagate_inbounds Base.getindex(
fv::FieldVector,
block::BlockArrays.Block{1},
) = backing_array(_values(fv)[block.n...])
Base.@propagate_inbounds function Base.getindex(
fv::FieldVector,
bidx::BlockArrays.BlockIndex{1},
)
X = fv[BlockArrays.block(bidx)]
X[bidx.α...]
end
# TODO: drop support for this
Base.@propagate_inbounds Base.getindex(fv::FieldVector, i::Integer) =
getindex(fv, BlockArrays.findblockindex(axes(fv, 1), i))
Base.@propagate_inbounds function Base.setindex!(
fv::FieldVector,
val,
bidx::BlockArrays.BlockIndex{1},
)
X = fv[BlockArrays.block(bidx)]
X[bidx.α...] = val
end
# TODO: drop support for this
Base.@propagate_inbounds Base.setindex!(fv::FieldVector, val, i::Integer) =
setindex!(fv, val, BlockArrays.findblockindex(axes(fv, 1), i))
Base.similar(fv::FieldVector{T}) where {T} =
FieldVector{T}(map(similar, _values(fv)))
Base.similar(fv::FieldVector{T}, ::Type{T}) where {T} =
FieldVector{T}(map(similar, _values(fv)))
_similar(x, ::Type{T}) where {T} = similar(x, T)
_similar(x::Field, ::Type{T}) where {T} =
Field(DataLayouts.replace_basetype(field_values(x), T), axes(x))
Base.similar(fv::FieldVector{T}, ::Type{T′}) where {T, T′} =
FieldVector{T′}(map(x -> _similar(x, T′), _values(fv)))
Base.copy(fv::FieldVector{T}) where {T} = FieldVector{T}(map(copy, _values(fv)))
Base.zero(fv::FieldVector{T}) where {T} = FieldVector{T}(map(zero, _values(fv)))
Base.@propagate_inbounds slab(fv::FieldVector{T}, inds...) where {T} =
FieldVector{T}(slab_args(_values(fv), inds...))
Base.@propagate_inbounds column(fv::FieldVector{T}, inds...) where {T} =
FieldVector{T}(column_args(_values(fv), inds...))
struct FieldVectorStyle <: Base.Broadcast.AbstractArrayStyle{1} end
Base.Broadcast.BroadcastStyle(::Type{<:FieldVector}) = FieldVectorStyle()
Base.Broadcast.BroadcastStyle(
fs::FieldVectorStyle,
as::Base.Broadcast.DefaultArrayStyle{0},
) = fs
Base.Broadcast.BroadcastStyle(
fs::FieldVectorStyle,
as::Base.Broadcast.AbstractArrayStyle{0},
) = fs
Base.Broadcast.BroadcastStyle(
fs::FieldVectorStyle,
as::Base.Broadcast.DefaultArrayStyle,
) = as
Base.Broadcast.BroadcastStyle(
fs::FieldVectorStyle,
as::Base.Broadcast.AbstractArrayStyle,
) = as
function Base.similar(
bc::Base.Broadcast.Broadcasted{FieldVectorStyle},
::Type{T},
) where {T}
for arg in bc.args
if arg isa FieldVector ||
arg isa Base.Broadcast.Broadcasted{FieldVectorStyle}
return similar(arg, T)
end
end
error("Cannot construct FieldVector")
end
"""
Spaces.create_dss_buffer(fv::FieldVector)
Create a NamedTuple of buffers for communicating neighbour information of
each Field in `fv`. In this NamedTuple, the name of each field is mapped
to the buffer.
"""
function Spaces.create_dss_buffer(fv::FieldVector)
NamedTuple{propertynames(fv)}(
map(
key -> Spaces.create_dss_buffer(getproperty(fv, key)),
propertynames(fv),
),
)
end
"""
Spaces.weighted_dss!(fv::FieldVector, dss_buffer = Spaces.create_dss_buffer(fv))
Apply weighted direct stiffness summation (DSS) to each field in `fv`.
If a `dss_buffer` object is not provided, a buffer will be created for each
field in `fv`.
Note that using the `Pair` interface here parallelizes the `weighted_dss!` calls.
"""
function Spaces.weighted_dss!(
fv::FieldVector,
dss_buffer = Spaces.create_dss_buffer(fv),
)
pairs = map(propertynames(fv)) do key
Pair(getproperty(fv, key), getproperty(dss_buffer, key))
end
Spaces.weighted_dss!(pairs...)
end
@inline function first_fieldvector_in_bc(args::Tuple, rargs...)
x1 = first_fieldvector_in_bc(args[1], rargs...)
x1 isa FieldVector && return x1
return first_fieldvector_in_bc(Base.tail(args), rargs...)
end
@inline first_fieldvector_in_bc(args::Tuple{Any}, rargs...) =
first_fieldvector_in_bc(args[1], rargs...)
@inline first_fieldvector_in_bc(args::Tuple{}, rargs...) = nothing
@inline first_fieldvector_in_bc(x) = nothing
@inline first_fieldvector_in_bc(x::FieldVector) = x
@inline first_fieldvector_in_bc(
bc::Base.Broadcast.Broadcasted{FieldVectorStyle},
) = first_fieldvector_in_bc(bc.args)
@inline _is_diagonal_bc_args(
truesofar,
::Type{TStart},
args::Tuple,
rargs...,
) where {TStart} =
truesofar &&
_is_diagonal_bc(truesofar, TStart, args[1], rargs...) &&
_is_diagonal_bc_args(truesofar, TStart, Base.tail(args), rargs...)
@inline _is_diagonal_bc_args(
truesofar,
::Type{TStart},
args::Tuple{Any},
rargs...,
) where {TStart} =
truesofar && _is_diagonal_bc(truesofar, TStart, args[1], rargs...)
@inline _is_diagonal_bc_args(
truesofar,
::Type{TStart},
args::Tuple{},
rargs...,
) where {TStart} = truesofar
@inline function _is_diagonal_bc(
truesofar,
::Type{TStart},
bc::Base.Broadcast.Broadcasted{FieldVectorStyle},
) where {TStart}
return truesofar && _is_diagonal_bc_args(truesofar, TStart, bc.args)
end
@inline _is_diagonal_bc(
truesofar,
::Type{TStart},
::TStart,
) where {TStart <: FieldVector} = true
@inline _is_diagonal_bc(
truesofar,
::Type{TStart},
x::FieldVector,
) where {TStart} = false
@inline _is_diagonal_bc(truesofar, ::Type{TStart}, x) where {TStart} = truesofar
# Find the first fieldvector in the broadcast expression (BCE),
# and compare against every other fieldvector in the BCE
@inline is_diagonal_bc(bc::Base.Broadcast.Broadcasted{FieldVectorStyle}) =
_is_diagonal_bc_args(true, typeof(first_fieldvector_in_bc(bc)), bc.args)
# Specialize on FieldVectorStyle to avoid inference failure
# in fieldvector broadcast expressions:
# https://github.com/JuliaArrays/BlockArrays.jl/issues/310
function Base.Broadcast.instantiate(
bc::Base.Broadcast.Broadcasted{FieldVectorStyle},
)
if bc.axes isa Nothing # Not done via dispatch to make it easier to extend instantiate(::Broadcasted{Style})
axes = Base.Broadcast.combine_axes(bc.args...)
else
axes = bc.axes
# Base.Broadcast.check_broadcast_axes is type-unstable
# for broadcast expressions with multiple fieldvectors.
# So, let's statically elide this when we have "diagonal"
# broadcast expressions:
if !is_diagonal_bc(bc)
Base.Broadcast.check_broadcast_axes(axes, bc.args...)
end
end
return Base.Broadcast.Broadcasted{FieldVectorStyle}(bc.f, bc.args, axes)
end
# Recursively call transform_bc_args() on broadcast arguments in a way that is statically reducible by the optimizer
# see Base.Broadcast.preprocess_args
@inline transform_bc_args(args::Tuple, inds...) = (
transform_broadcasted(args[1], inds...),
transform_bc_args(Base.tail(args), inds...)...,
)
@inline transform_bc_args(args::Tuple{Any}, inds...) =
(transform_broadcasted(args[1], inds...),)
@inline transform_bc_args(args::Tuple{}, inds...) = ()
@inline function transform_broadcasted(
bc::Base.Broadcast.Broadcasted{FieldVectorStyle},
symb,
axes,
)
Base.Broadcast.Broadcasted(
bc.f,
transform_bc_args(bc.args, symb, axes),
axes,
)
end
@inline transform_broadcasted(fv::FieldVector, symb, axes) =
parent(getfield(_values(fv), symb))
@inline transform_broadcasted(x, symb, axes) = x
@inline function Base.copyto!(
dest::FieldVector,
bc::Union{FieldVector, Base.Broadcast.Broadcasted{FieldVectorStyle}},
)
copyto_per_field!(dest, bc)
call_post_op_callback() && post_op_callback(dest, dest, bc)
return dest
end
@inline function copyto_per_field!(
dest::FieldVector,
bc::Union{FieldVector, Base.Broadcast.Broadcasted{FieldVectorStyle}},
)
map(propertynames(dest)) do symb
Base.@_inline_meta
array = parent(getfield(_values(dest), symb))
bct = transform_broadcasted(bc, symb, axes(array))
if array isa FieldVector # recurse
copyto_per_field!(array, bct)
else
copyto_per_field!(
array,
Base.Broadcast.instantiate(bct),
DataLayouts.device_dispatch(array),
)
end
end
return dest
end
@inline function Base.copyto!(
dest::FieldVector,
bc::Base.Broadcast.Broadcasted{<:Base.Broadcast.Style{Tuple}},
)
copyto_per_field_scalar!(dest, bc)
call_post_op_callback() && post_op_callback(dest, dest, bc)
return dest
end
@inline function Base.copyto!(
dest::FieldVector,
bc::Base.Broadcast.Broadcasted{<:Base.Broadcast.AbstractArrayStyle{0}},
)
copyto_per_field_scalar!(dest, bc)
call_post_op_callback() && post_op_callback(dest, dest, bc)
return dest
end
@inline function Base.copyto!(dest::FieldVector, bc::Real)
copyto_per_field_scalar!(dest, bc)
call_post_op_callback() && post_op_callback(dest, dest, bc)
return dest
end
@inline function copyto_per_field_scalar!(dest::FieldVector, bc)
map(propertynames(dest)) do symb
Base.@_inline_meta
array = parent((getfield(_values(dest), symb)))
if array isa FieldVector # recurse
copyto_per_field_scalar!(array, bc)
else
copyto_per_field_scalar!(
array,
Base.Broadcast.instantiate(bc),
DataLayouts.device_dispatch(array),
)
end
nothing
end
return dest
end
Base.fill!(dest::FieldVector, value) = dest .= value
Base.mapreduce(f, op, fv::FieldVector) =
mapreduce(x -> mapreduce(f, op, backing_array(x)), op, _values(fv))
Base.any(f, fv::FieldVector) = any(x -> any(f, backing_array(x)), _values(fv))
Base.any(f::Function, fv::FieldVector) = # avoid ambiguities
any(x -> any(f, backing_array(x)), _values(fv))
Base.any(fv::FieldVector) = any(identity, A)
Base.all(f, fv::FieldVector) = all(x -> all(f, backing_array(x)), _values(fv))
Base.all(f::Function, fv::FieldVector) =
all(x -> all(f, backing_array(x)), _values(fv))
Base.all(fv::FieldVector) = all(identity, fv)
# TODO: figure out a better way to handle these
# https://github.com/JuliaArrays/BlockArrays.jl/issues/185
LinearAlgebra.ldiv!(
x::FieldVector,
A::LinearAlgebra.QRCompactWY,
b::FieldVector,
) = x .= LinearAlgebra.ldiv!(A, Vector(b))
LinearAlgebra.ldiv!(A::LinearAlgebra.QRCompactWY, x::FieldVector) =
x .= LinearAlgebra.ldiv!(A, Vector(x))
LinearAlgebra.ldiv!(x::FieldVector, A::LinearAlgebra.LU, b::FieldVector) =
x .= LinearAlgebra.ldiv!(A, Vector(b))
LinearAlgebra.ldiv!(A::LinearAlgebra.LU, x::FieldVector) =
x .= LinearAlgebra.ldiv!(A, Vector(x))
function LinearAlgebra.norm_sqr(x::FieldVector)
value_norm_sqrs = unrolled_map(_values(x)) do value
LinearAlgebra.norm_sqr(backing_array(value))
end
return sum(value_norm_sqrs; init = zero(eltype(x)))
end
function LinearAlgebra.norm(x::FieldVector)
sqrt(LinearAlgebra.norm_sqr(x))
end
import ClimaComms
ClimaComms.array_type(x::FieldVector) =
promote_type(unrolled_map(ClimaComms.array_type, _values(x))...)
function __rprint_diff(
io::IO,
x::T,
y::T;
pc,
xname,
yname,
) where {T <: Union{FieldVector, Field, DataLayouts.AbstractData, NamedTuple}}
for pn in propertynames(x)
pc_full = (pc..., ".", pn)
xi = getproperty(x, pn)
yi = getproperty(y, pn)
__rprint_diff(io, xi, yi; pc = pc_full, xname, yname)
end
end;
function __rprint_diff(io::IO, xi, yi; pc, xname, yname) # assume we can compute difference here
if !(xi == yi)
xs = xname * string(join(pc))
ys = yname * string(join(pc))
println(io, "==================== Difference found:")
println(io, "$xs: ", xi)
println(io, "$ys: ", yi)
println(io, "($xs .- $ys): ", (xi .- yi))
end
return nothing
end
"""
rprint_diff(io::IO, ::T, ::T) where {T <: Union{FieldVector, NamedTuple}}
rprint_diff(::T, ::T) where {T <: Union{FieldVector, NamedTuple}}
Recursively print differences in given `Union{FieldVector, NamedTuple}`.
"""
_rprint_diff(
io::IO,
x::T,
y::T,
xname,
yname,
) where {T <: Union{FieldVector, NamedTuple}} =
__rprint_diff(io, x, y; pc = (), xname, yname)
_rprint_diff(
x::T,
y::T,
xname,
yname,
) where {T <: Union{FieldVector, NamedTuple}} =
_rprint_diff(stdout, x, y, xname, yname)
"""
@rprint_diff(::T, ::T) where {T <: Union{FieldVector, NamedTuple}}
Recursively print differences in given `Union{FieldVector, NamedTuple}`.
"""
macro rprint_diff(x, y)
return :(_rprint_diff(
stdout,
$(esc(x)),
$(esc(y)),
$(string(x)),
$(string(y)),
))
end
# Recursively compare contents of similar fieldvectors
_rcompare(pass, x::T, y::T; strict) where {T <: Field} =
pass && _rcompare(pass, field_values(x), field_values(y); strict)
_rcompare(pass, x::T, y::T; strict) where {T <: DataLayouts.AbstractData} =
pass && (parent(x) == parent(y))
_rcompare(pass, x::T, y::T; strict) where {T} = pass && (x == y)
_rcompare(pass, x::NamedTuple, y::NamedTuple; strict) =
_rcompare_nt(pass, x, y; strict)
_rcompare(pass, x::FieldVector, y::FieldVector; strict) =
_rcompare_nt(pass, x, y; strict)
function _rcompare_nt(pass, x, y; strict)
length(propertynames(x)) ≠ length(propertynames(y)) && return false
if strict
typeof(x) == typeof(y) || return false
end
for pn in propertynames(x)
pass &= _rcompare(pass, getproperty(x, pn), getproperty(y, pn); strict)
end
return pass
end
"""
rcompare(x::T, y::T; strict = true) where {T <: Union{FieldVector, NamedTuple}}
Recursively compare given fieldvectors via `==`.
Returns `true` if `x == y` recursively.
The keyword `strict = true` allows users to additionally
check that the types match. If `strict = false`, then
`rcompare` will return `true` for `FieldVector`s and
`NamedTuple`s with the same properties but permuted order.
For example:
- `rcompare((;a=1,b=2), (;b=2,a=1); strict = true)` will return `false` and
- `rcompare((;a=1,b=2), (;b=2,a=1); strict = false)` will return `true`
"""
rcompare(
x::T,
y::T;
strict = true,
) where {T <: Union{FieldVector, NamedTuple}} = _rcompare(true, x, y; strict)
rcompare(x::T, y::T; strict = true) where {T <: FieldVector} =
_rcompare(true, x, y; strict)
rcompare(x::T, y::T; strict = true) where {T <: NamedTuple} =
_rcompare(true, x, y; strict)
# FieldVectors with different types are always different
rcompare(x::FieldVector, y::FieldVector; strict::Bool = true) =
strict ? false : _rcompare(true, x, y; strict)
rcompare(x::NamedTuple, y::NamedTuple; strict::Bool = true) =
strict ? false : _rcompare(true, x, y; strict)
# Define == to call rcompare for two fieldvectors
Base.:(==)(x::FieldVector, y::FieldVector) = rcompare(x, y; strict = true)