changeset 354:db5f7b51038c refactor/remove_dynamic_size_tensormapping

Merge in default
author Jonatan Werpers <jonatan@werpers.com>
date Sun, 27 Sep 2020 13:48:13 +0200
parents 7fe43d902a27 (diff) 28e71a861531 (current diff)
children 5c9212a8ee4f
files
diffstat 6 files changed, 65 insertions(+), 62 deletions(-) [+]
line wrap: on
line diff
--- a/src/LazyTensors/lazy_tensor_operations.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/src/LazyTensors/lazy_tensor_operations.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -11,12 +11,13 @@
     t::TensorMapping{T,R,D}
     o::AbstractArray{T,D}
 end
+# TODO: Do boundschecking on creation!
 export LazyTensorMappingApplication
 
 Base.:*(tm::TensorMapping{T,R,D}, o::AbstractArray{T,D}) where {T,R,D} = LazyTensorMappingApplication(tm,o)
 Base.getindex(ta::LazyTensorMappingApplication{T,R,D}, I::Vararg{Index,R}) where {T,R,D} = apply(ta.t, ta.o, I...)
 Base.getindex(ta::LazyTensorMappingApplication{T,R,D}, I::Vararg{Int,R}) where {T,R,D} = apply(ta.t, ta.o, Index{Unknown}.(I)...)
-Base.size(ta::LazyTensorMappingApplication{T,R,D}) where {T,R,D} = range_size(ta.t,size(ta.o))
+Base.size(ta::LazyTensorMappingApplication{T,R,D}) where {T,R,D} = range_size(ta.t)
 # TODO: What else is needed to implement the AbstractArray interface?
 
 # # We need the associativity to be a→b→c = a→(b→c), which is the case for '→'
@@ -41,16 +42,15 @@
 export LazyTensorMappingTranspose
 
 # # TBD: Should this be implemented on a type by type basis or through a trait to provide earlier errors?
+# Jonatan 2020-09-25: Is the problem that you can take the transpose of any TensorMapping even if it doesn't implement `apply_transpose`?
 Base.adjoint(tm::TensorMapping) = LazyTensorMappingTranspose(tm)
 Base.adjoint(tmt::LazyTensorMappingTranspose) = tmt.tm
 
 apply(tmt::LazyTensorMappingTranspose{T,R,D}, v::AbstractArray{T,R}, I::Vararg{Index,D}) where {T,R,D} = apply_transpose(tmt.tm, v, I...)
 apply_transpose(tmt::LazyTensorMappingTranspose{T,R,D}, v::AbstractArray{T,D}, I::Vararg{Index,R}) where {T,R,D} = apply(tmt.tm, v, I...)
 
-range_size(tmt::LazyTensorMappingTranspose{T,R,D}, d_size::NTuple{R,Integer}) where {T,R,D} = domain_size(tmt.tm, d_size)
-domain_size(tmt::LazyTensorMappingTranspose{T,R,D}, r_size::NTuple{D,Integer}) where {T,R,D} = range_size(tmt.tm, r_size)
-
-
+range_size(tmt::LazyTensorMappingTranspose{T,R,D}) where {T,R,D} = domain_size(tmt.tm)
+domain_size(tmt::LazyTensorMappingTranspose{T,R,D}) where {T,R,D} = range_size(tmt.tm)
 
 
 struct LazyTensorMappingBinaryOperation{Op,T,R,D,T1<:TensorMapping{T,R,D},T2<:TensorMapping{T,R,D}} <: TensorMapping{T,D,R}
@@ -61,12 +61,13 @@
         return new{Op,T,R,D,T1,T2}(tm1,tm2)
     end
 end
+# TODO: Boundschecking in constructor.
 
 apply(tmBinOp::LazyTensorMappingBinaryOperation{:+,T,R,D}, v::AbstractArray{T,D}, I::Vararg{Index,R}) where {T,R,D} = apply(tmBinOp.tm1, v, I...) + apply(tmBinOp.tm2, v, I...)
 apply(tmBinOp::LazyTensorMappingBinaryOperation{:-,T,R,D}, v::AbstractArray{T,D}, I::Vararg{Index,R}) where {T,R,D} = apply(tmBinOp.tm1, v, I...) - apply(tmBinOp.tm2, v, I...)
 
-range_size(tmBinOp::LazyTensorMappingBinaryOperation{Op,T,R,D}, domain_size::NTuple{D,Integer}) where {Op,T,R,D} = range_size(tmBinOp.tm1, domain_size)
-domain_size(tmBinOp::LazyTensorMappingBinaryOperation{Op,T,R,D}, range_size::NTuple{R,Integer}) where {Op,T,R,D} = domain_size(tmBinOp.tm2, range_size)
+range_size(tmBinOp::LazyTensorMappingBinaryOperation{Op,T,R,D}) where {Op,T,R,D} = range_size(tmBinOp.tm1)
+domain_size(tmBinOp::LazyTensorMappingBinaryOperation{Op,T,R,D}) where {Op,T,R,D} = domain_size(tmBinOp.tm1)
 
 Base.:+(tm1::TensorMapping{T,R,D}, tm2::TensorMapping{T,R,D}) where {T,R,D} = LazyTensorMappingBinaryOperation{:+,T,R,D}(tm1,tm2)
 Base.:-(tm1::TensorMapping{T,R,D}, tm2::TensorMapping{T,R,D}) where {T,R,D} = LazyTensorMappingBinaryOperation{:-,T,R,D}(tm1,tm2)
--- a/src/LazyTensors/tensor_mapping.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/src/LazyTensors/tensor_mapping.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -6,11 +6,11 @@
 
     apply(t::TensorMapping{T,R,D}, v::AbstractArray{T,D}, I::Vararg) where {R,D,T}
 
-The size of range tensor should be dependent on the size of the domain tensor
-and the type should implement the methods
+The size of the range and domain that the operator works with should be returned by
+the functions
 
-    range_size(::TensorMapping{T,R,D}, domain_size::NTuple{D,Integer}) where {T,R,D}
-    domain_size(::TensorMapping{T,R,D}, range_size::NTuple{R,Integer}) where {T,R,D}
+    range_size(::TensorMapping)
+    domain_size(::TensorMapping)
 
 to allow querying for one or the other.
 
@@ -49,35 +49,19 @@
 export range_dim, domain_dim
 
 """
-    range_size(M::TensorMapping, domain_size)
+    range_size(M::TensorMapping)
 
-Return the resulting range size for the mapping applied to a given domain_size
+Return the range size for the mapping.
 """
 function range_size end
 
 """
-    domain_size(M::TensorMapping, range_size)
+    domain_size(M::TensorMapping)
 
-Return the resulting domain size for the mapping applied to a given range_size
+Return the domain size for the mapping.
 """
 function domain_size end
 
-"""
-    Dummy type for representing dimensions of tensormappings when domain_size is unknown
-"""
-struct UnknownDim end
-export range_size, domain_size, TensorMappingDim, UnknownDim
+export range_size, domain_size
 
 # TODO: Think about boundschecking!
-
-
-"""
-    TensorOperator{T,D}
-
-A `TensorMapping{T,D,D}` where the range and domain tensor have the same number of
-dimensions and the same size.
-"""
-abstract type TensorOperator{T,D} <: TensorMapping{T,D,D} end
-export TensorOperator
-domain_size(::TensorOperator{T,D}, range_size::NTuple{D,Integer}) where {T,D} = range_size
-range_size(::TensorOperator{T,D}, domain_size::NTuple{D,Integer}) where {T,D} = domain_size
--- a/src/SbpOperators/SbpOperators.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/src/SbpOperators/SbpOperators.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -13,4 +13,21 @@
 include("quadrature/quadrature.jl")
 include("quadrature/inverse_diagonal_inner_product.jl")
 include("quadrature/inverse_quadrature.jl")
+
+abstract type SbpOperator{T,R,D} <: TensorMapping{T,R,D} end
+
+"""
+    grid(::ColocationOperator)
+
+Return the the grid which the sbp-operator lives on
+"""
+function grid end
+
+abstract type ColocationOperator{T,R,D} <: SbpOperator{T,R,D} end
+
+LazyTensors.range_size(co::ColocationOperator) = size(grid(co))
+LazyTensors.domain_size(co::ColocationOperator) = size(grid(co))
+
+# Allt ovan kanske är overkill.. Eventuellt bara lättare och tydligare att alla typer definerar sina range och domain size hur dom vill. (I praktiken typ alltid genom att ha gridden som ett fält i structen.)
+
 end # module
--- a/src/SbpOperators/laplace/secondderivative.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/src/SbpOperators/laplace/secondderivative.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -4,7 +4,7 @@
 in 1D dimension
 """
 
-struct SecondDerivative{T,N,M,K} <: TensorOperator{T,1}
+struct SecondDerivative{T,N,M,K} <: TensorMapping{T,1,1}
     h_inv::T # The grid spacing could be included in the stencil already. Preferable?
     innerStencil::Stencil{T,N}
     closureStencils::NTuple{M,Stencil{T,K}}
--- a/src/Sbplib.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/src/Sbplib.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -3,8 +3,8 @@
 include("RegionIndices/RegionIndices.jl")
 include("LazyTensors/LazyTensors.jl")
 include("Grids/Grids.jl")
-include("SbpOperators/SbpOperators.jl")
-include("DiffOps/DiffOps.jl")
+# include("SbpOperators/SbpOperators.jl")
+# include("DiffOps/DiffOps.jl")
 
 export RegionIndices
 export LazyTensors
--- a/test/testLazyTensors.jl	Sun Sep 27 13:47:40 2020 +0200
+++ b/test/testLazyTensors.jl	Sun Sep 27 13:48:13 2020 +0200
@@ -12,20 +12,14 @@
     @test apply(DummyMapping{Int,2,3}(), zeros(Int, (0,0,0)),(Index{Unknown}(0),Index{Unknown}(0))) == :apply
 end
 
-@testset "Generic Operator methods" begin
-    struct DummyOperator{T,D} <: TensorOperator{T,D} end
-    @test range_size(DummyOperator{Int,2}(), (3,5)) == (3,5)
-    @test domain_size(DummyOperator{Float64, 3}(), (3,3,1)) == (3,3,1)
-end
-
 @testset "Mapping transpose" begin
     struct DummyMapping{T,R,D} <: TensorMapping{T,R,D} end
 
     LazyTensors.apply(m::DummyMapping{T,R,D}, v, I::Vararg{Index{<:Region},R}) where {T,R,D} = :apply
     LazyTensors.apply_transpose(m::DummyMapping{T,R,D}, v, I::Vararg{Index{<:Region},D}) where {T,R,D} = :apply_transpose
 
-    LazyTensors.range_size(m::DummyMapping{T,R,D}, domain_size::NTuple{D,Integer}) where {T,R,D} = :range_size
-    LazyTensors.domain_size(m::DummyMapping{T,R,D}, range_size::NTuple{R,Integer}) where {T,R,D} = :domain_size
+    LazyTensors.range_size(m::DummyMapping{T,R,D}) where {T,R,D} = :range_size
+    LazyTensors.domain_size(m::DummyMapping{T,R,D}) where {T,R,D} = :domain_size
 
     m = DummyMapping{Float64,2,3}()
     I = Index{Unknown}(0)
@@ -35,19 +29,21 @@
     @test apply(m'',zeros(Float64,(0,0,0)), I, I) == :apply
     @test apply_transpose(m', zeros(Float64,(0,0,0)), I, I) == :apply
 
-    @test range_size(m', (0,0)) == :domain_size
-    @test domain_size(m', (0,0,0)) == :range_size
+    @test range_size(m') == :domain_size
+    @test domain_size(m') == :range_size
 end
 
 @testset "TensorApplication" begin
-    struct DummyMapping{T,R,D} <: TensorMapping{T,R,D} end
+    struct SizeDoublingMapping{T,R,D} <: TensorMapping{T,R,D}
+        domain_size::NTuple{D,Int}
+    end
 
-    LazyTensors.apply(m::DummyMapping{T,R,D}, v, i::Vararg{Index{<:Region},R}) where {T,R,D} = (:apply,v,i)
-    LazyTensors.range_size(m::DummyMapping{T,R,D}, domain_size::NTuple{D,Integer}) where {T,R,D} = 2 .* domain_size
-    LazyTensors.domain_size(m::DummyMapping{T,R,D}, range_size::NTuple{R,Integer}) where {T,R,D} = range_size.÷2
+    LazyTensors.apply(m::SizeDoublingMapping{T,R,D}, v, i::Vararg{Index{<:Region},R}) where {T,R,D} = (:apply,v,i)
+    LazyTensors.range_size(m::SizeDoublingMapping) = 2 .* m.domain_size
+    LazyTensors.domain_size(m::SizeDoublingMapping) = m.domain_size
 
 
-    m = DummyMapping{Int, 1, 1}()
+    m = SizeDoublingMapping{Int, 1, 1}((3,))
     v = [0,1,2]
     @test m*v isa AbstractVector{Int}
     @test size(m*v) == 2 .*size(v)
@@ -63,28 +59,31 @@
     @test_broken BoundsError == (m*m*v)[0]
     @test_broken BoundsError == (m*m*v)[7]
 
-    m = DummyMapping{Int, 2, 1}()
+    m = SizeDoublingMapping{Int, 2, 1}((3,))
     @test_throws MethodError m*ones(Int,2,2)
     @test_throws MethodError m*m*v
 
-    m = DummyMapping{Float64, 2, 2}()
+    m = SizeDoublingMapping{Float64, 2, 2}((3,3))
     v = ones(3,3)
     I = (Index{Lower}(1),Index{Interior}(2));
     @test size(m*v) == 2 .*size(v)
     @test (m*v)[I] == (:apply,v,I)
 
-    struct ScalingOperator{T,D} <: TensorOperator{T,D}
+    struct ScalingOperator{T,D} <: TensorMapping{T,D,D}
         λ::T
+        size::NTuple{D,Int}
     end
 
     LazyTensors.apply(m::ScalingOperator{T,D}, v, I::Vararg{Index,D}) where {T,D} = m.λ*v[I]
+    LazyTensors.range_size(m::ScalingOperator) = m.size
+    LazyTensors.domain_size(m::ScalingOperator) = m.size
 
-    m = ScalingOperator{Int,1}(2)
+    m = ScalingOperator{Int,1}(2,(3,))
     v = [1,2,3]
     @test m*v isa AbstractVector
     @test m*v == [2,4,6]
 
-    m = ScalingOperator{Int,2}(2)
+    m = ScalingOperator{Int,2}(2,(2,2))
     v = [[1 2];[3 4]]
     @test m*v == [[2 4];[6 8]]
     I = (Index{Upper}(2),Index{Lower}(1))
@@ -94,14 +93,16 @@
 @testset "TensorMapping binary operations" begin
     struct ScalarMapping{T,R,D} <: TensorMapping{T,R,D}
         λ::T
+        range_size::NTuple{R,Int}
+        domain_size::NTuple{D,Int}
     end
 
     LazyTensors.apply(m::ScalarMapping{T,R,D}, v, I::Vararg{Index{<:Region}}) where {T,R,D} = m.λ*v[I...]
-    LazyTensors.range_size(m::ScalarMapping, domain_size) = domain_size
-    LazyTensors.domain_size(m::ScalarMapping, range_sizes) = range_sizes
+    LazyTensors.range_size(m::ScalarMapping) = m.domain_size
+    LazyTensors.domain_size(m::ScalarMapping) = m.range_size
 
-    A = ScalarMapping{Float64,1,1}(2.0)
-    B = ScalarMapping{Float64,1,1}(3.0)
+    A = ScalarMapping{Float64,1,1}(2.0, (3,), (3,))
+    B = ScalarMapping{Float64,1,1}(3.0, (3,), (3,))
 
     v = [1.1,1.2,1.3]
     for i ∈ eachindex(v)
@@ -112,8 +113,8 @@
         @test ((A-B)*v)[i] == 2*v[i] - 3*v[i]
     end
 
-    @test range_size(A+B, (3,)) == range_size(A, (3,)) == range_size(B,(3,))
-    @test domain_size(A+B, (3,)) == domain_size(A, (3,)) == domain_size(B,(3,))
+    @test range_size(A+B) == range_size(A) == range_size(B)
+    @test domain_size(A+B) == domain_size(A) == domain_size(B)
 end
 
 @testset "LazyArray" begin
@@ -192,4 +193,4 @@
     @test_throws DimensionMismatch v1 + v2
 end
 
-end
\ No newline at end of file
+end