changeset 1329:e94ddef5e72f refactor/grids

Clean up documentation for changed types in SbpOperatorClean up documentation for changed types in SbpOperatorss
author Jonatan Werpers <jonatan@werpers.com>
date Tue, 02 May 2023 22:09:33 +0200
parents f00a205ae347
children 5f05a708d730
files grid_refactor.md src/SbpOperators/boundaryops/boundary_operator.jl src/SbpOperators/boundaryops/boundary_restriction.jl src/SbpOperators/boundaryops/normal_derivative.jl src/SbpOperators/volumeops/constant_interior_scaling_operator.jl src/SbpOperators/volumeops/derivatives/dissipation.jl src/SbpOperators/volumeops/derivatives/first_derivative.jl src/SbpOperators/volumeops/derivatives/second_derivative.jl src/SbpOperators/volumeops/inner_products/inner_product.jl src/SbpOperators/volumeops/inner_products/inverse_inner_product.jl src/SbpOperators/volumeops/laplace/laplace.jl src/SbpOperators/volumeops/stencil_operator_distinct_closures.jl src/SbpOperators/volumeops/volume_operator.jl
diffstat 13 files changed, 129 insertions(+), 63 deletions(-) [+]
line wrap: on
line diff
--- a/grid_refactor.md	Tue May 02 20:14:39 2023 +0200
+++ b/grid_refactor.md	Tue May 02 22:09:33 2023 +0200
@@ -21,7 +21,6 @@
 * Write about the design choices in the docs.
 * Merge and run benchmarks
 
-* Check all the docstring of all types that have been changed
 * Clean out Notes.md of any solved issues
 * Delete this document, move remaining notes to Notes.md
 
--- a/src/SbpOperators/boundaryops/boundary_operator.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/boundaryops/boundary_operator.jl	Tue May 02 22:09:33 2023 +0200
@@ -3,9 +3,10 @@
 
 Implements the boundary operator `op` for 1D as a `LazyTensor`
 
-`op` is the restriction of a grid function to the boundary using some closure `Stencil{T,N}`.
-The boundary to restrict to is determined by `R`.
-`op'` is the prolongation of a zero dimensional array to the whole grid using the same closure stencil.
+`op` is the restriction of a grid function to the boundary using some closure
+`Stencil{T,N}`. The boundary to restrict to is determined by `R`. `op'` is the
+prolongation of a zero dimensional array to the whole grid using the same
+closure stencil.
 """
 struct BoundaryOperator{T,R<:Region,N} <: LazyTensor{T,0,1}
     stencil::Stencil{T,N}
@@ -15,8 +16,9 @@
 """
     BoundaryOperator(grid::EquidistantGrid, closure_stencil, region)
 
-Constructs the BoundaryOperator with stencil `closure_stencil` for a `EquidistantGrid` `grid`, restricting to
-to the boundary specified by `region`.
+Constructs the BoundaryOperator with stencil `closure_stencil` for a
+`EquidistantGrid` `grid`, restricting to to the boundary specified by
+`region`.
 """
 function BoundaryOperator(grid::EquidistantGrid, closure_stencil::Stencil{T,N}, region::Region) where {T,N}
     return BoundaryOperator{T,typeof(region),N}(closure_stencil,size(grid)[1])
--- a/src/SbpOperators/boundaryops/boundary_restriction.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/boundaryops/boundary_restriction.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,16 +1,20 @@
 """
-    boundary_restriction(g, closure_stencil::Stencil, boundary)
+    boundary_restriction(g, stencil_set::StencilSet, boundary)
+    boundary_restriction(g::TensorGrid, stencil_set::StencilSet, boundary::TensorGridBoundary)
+    boundary_restriction(g::EquidistantGrid, stencil_set::StencilSet, boundary)
 
 Creates boundary restriction operators `e` as `LazyTensor`s on `boundary`
 
-`e` is the restriction of a grid function to `boundary` using a `Stencil` `closure_stencil`.
-`e'` is the prolongation of a grid function on `boundary` to the whole grid using the same `closure_stencil`.
-On a one-dimensional grid, `e` is a `BoundaryOperator`. On a multi-dimensional grid, `e` is the inflation of
-a `BoundaryOperator`.
+`e` is the restriction of a grid function to `boundary` using the 'e' stencil
+in the guven stencil set. `e'` is the prolongation of a grid function on
+`boundary` to the whole grid using the same stencil. On a one-dimensional
+grid, `e` is a `BoundaryOperator`. On a multi-dimensional grid, `e` is the
+inflation of a `BoundaryOperator`.
 
 See also: [`BoundaryOperator`](@ref), [`LazyTensors.inflate`](@ref).
 """
-#TODO: Check docstring
+function boundary_restriction end
+
 function boundary_restriction(g::TensorGrid, stencil_set::StencilSet, boundary::TensorGridBoundary)
     op = boundary_restriction(g.grids[grid_id(boundary)], stencil_set, boundary_id(boundary))
     return LazyTensors.inflate(op, size(g), grid_id(boundary))
--- a/src/SbpOperators/boundaryops/normal_derivative.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/boundaryops/normal_derivative.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,16 +1,21 @@
 """
-    normal_derivative(g, closure_stencil::Stencil, boundary)
+    normal_derivative(g, stencil_set::StencilSet, boundary)
+    normal_derivative(g::TensorGrid, stencil_set::StencilSet, boundary::TensorGridBoundary)
+    normal_derivative(g::EquidistantGrid, stencil_set::StencilSet, boundary)
 
 Creates the normal derivative boundary operator `d` as a `LazyTensor`
 
-`d` computes the normal derivative of a grid function  on `boundary` a `Stencil` `closure_stencil`.
-`d'` is the prolongation of the normal derivative of a grid function to the whole grid using the same `closure_stencil`.
-On a one-dimensional `grid`, `d` is a `BoundaryOperator`. On a multi-dimensional `grid`, `d` is the inflation of
-a `BoundaryOperator`.
+`d` computes the normal derivative of a grid function  on `boundary` using the
+'d1' stencil in the given stencil_set. `d'` is the prolongation of the normal
+derivative of a grid function to the whole grid using the same stencil. On a
+one-dimensional `grid`, `d` is a `BoundaryOperator`. On a multi-dimensional
+`grid`, `d` is the inflation of a `BoundaryOperator`.
 
 See also: [`BoundaryOperator`](@ref), [`LazyTensors.inflate`](@ref).
 """
-#TODO: Check docstring
+function normal_derivative end
+
+
 function normal_derivative(g::TensorGrid, stencil_set::StencilSet, boundary::TensorGridBoundary)
     op = normal_derivative(g.grids[grid_id(boundary)], stencil_set, boundary_id(boundary))
     return LazyTensors.inflate(op, size(g), grid_id(boundary))
--- a/src/SbpOperators/volumeops/constant_interior_scaling_operator.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/constant_interior_scaling_operator.jl	Tue May 02 22:09:33 2023 +0200
@@ -2,7 +2,8 @@
     ConstantInteriorScalingOperator{T,N} <: LazyTensor{T,1,1}
 
 A one-dimensional operator scaling a vector. The first and last `N` points are
-scaled with individual weights while all interior points are scaled the same.
+scaled with individual weights while all interior points are scaled using the
+same factor.
 """
 struct ConstantInteriorScalingOperator{T,N} <: LazyTensor{T,1,1}
     interior_weight::T
--- a/src/SbpOperators/volumeops/derivatives/dissipation.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/derivatives/dissipation.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,7 +1,8 @@
 """
-    undivided_skewed04(g::EquidistantGrid, p, direction)
+    undivided_skewed04(g::TensorGrid, p, direction)
+    undivided_skewed04(g::EquidistantGrid, p)
 
-Create undivided difference operators approximating the `p`th derivative. The
+Undivided difference operators approximating the `p`th derivative. The
 operators do not satisfy any SBP-property and are meant to be used for
 building artificial dissipation terms.
 
@@ -10,6 +11,8 @@
 Artificial Dissipation,” Journal of Scientific Computing, vol. 21, no. 1, pp.
 57–79, Aug. 2004"
 """
+function undivided_skewed04 end
+
 function undivided_skewed04(g::TensorGrid, p, direction)
     op = undivided_skewed04(g.grids[direction], p)
     return LazyTensors.inflate(op, size(g), direction)
--- a/src/SbpOperators/volumeops/derivatives/first_derivative.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/derivatives/first_derivative.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,12 +1,16 @@
+"""
+    first_derivative(g, ..., [direction])
+
+The first-derivative operator `D1` as a `LazyTensor` on the given grid.
+
+`D1` approximates the first-derivative d/dξ on `g` along the coordinate
+dimension specified by `direction`.
+"""
+function first_derivative end
+
 """
     first_derivative(g::TensorGrid, stencil_set, direction)
 
-Creates the first-derivative operator `D1` as a `LazyTensor`
-
-`D1` approximates the first-derivative d/dξ on `g` along the coordinate dimension specified by
-`direction`, using the stencil `inner_stencil` in the interior and a set of stencils `closure_stencils`
-for the points in the closure regions.
-
 See also: [`VolumeOperator`](@ref), [`LazyTensors.inflate`](@ref).
 """
 function first_derivative(g::TensorGrid, stencil_set, direction)
@@ -17,7 +21,8 @@
 """
     first_derivative(g::EquidistantGrid, stencil_set)
 
-Creates a `first_derivative` operator on a `EquidistantGrid` given a `StencilSet`.
+The first derivative operator on a `EquidistantGrid` given a
+`StencilSet`. Uses the `D1` stencil in the stencil set.
 """
 function first_derivative(g::EquidistantGrid, stencil_set::StencilSet)
     inner_stencil = parse_stencil(stencil_set["D1"]["inner_stencil"])
@@ -28,7 +33,8 @@
 """
     first_derivative(g::EquidistantGrid, inner_stencil, closure_stencils)
 
-Creates a `first_derivative` operator on a `EquidistantGrid` given an `inner_stencil` and a`closure_stencils`.
+The first derivative operator on a `EquidistantGrid` given an
+`inner_stencil` and a`closure_stencils`.
 """
 function first_derivative(g::EquidistantGrid, inner_stencil::Stencil, closure_stencils)
     h⁻¹ = inverse_spacing(g)
--- a/src/SbpOperators/volumeops/derivatives/second_derivative.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/derivatives/second_derivative.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,11 +1,10 @@
 """
-    second_derivative(g::EquidistantGrid, inner_stencil, closure_stencils, direction)
+    second_derivative(g::EquidistantGrid, stencil_set, direction)
 
 Creates the second-derivative operator `D2` as a `LazyTensor`
 
-`D2` approximates the second-derivative d²/dξ² on `g` along the coordinate dimension specified by
-`direction`, using the stencil `inner_stencil` in the interior and a set of stencils `closure_stencils`
-for the points in the closure regions.
+`D2` approximates the second-derivative d²/dξ² on `g` along the coordinate
+dimension specified by `direction`.
 
 See also: [`VolumeOperator`](@ref), [`LazyTensors.inflate`](@ref).
 """
@@ -17,7 +16,8 @@
 """
     second_derivative(g, stencil_set)
 
-Creates a `second_derivative` operator on a 1D `g` given a `stencil_set`.
+Creates a `second_derivative` operator on a 1D `g` given a `stencil_set`. Uses
+the `D2` stencil in the stencil set.
 """
 function second_derivative(g::EquidistantGrid, stencil_set::StencilSet)
     inner_stencil = parse_stencil(stencil_set["D2"]["inner_stencil"])
@@ -28,7 +28,8 @@
 """
     second_derivative(g, inner_stencil, closure_stencils)
 
-Creates a `second_derivative` operator on a 1D `g` given `inner_stencil` and `closure_stencils`.
+Creates a `second_derivative` operator on a 1D `g` given `inner_stencil` and
+`closure_stencils`.
 """
 function second_derivative(g::EquidistantGrid, inner_stencil::Stencil, closure_stencils)
     h⁻¹ = inverse_spacing(g)
--- a/src/SbpOperators/volumeops/inner_products/inner_product.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/inner_products/inner_product.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,31 +1,52 @@
 """
-    inner_product(grid::EquidistantGrid, interior_weight, closure_weights)
-
-Creates the discrete inner product operator `H` as a `LazyTensor` on an
-equidistant grid, defined as `(u,v)  = u'Hv` for grid functions `u,v`.
+    inner_product(grid, ...)
 
-`inner_product` creates `H` on `grid` using the `interior_weight` for the
-interior points and the `closure_weights` for the points close to the
-boundary.
+The inner product on a given grid with weights from a stencils set or given
+explicitly.
+"""
+function inner_product end
 
-On a 1-dimensional grid, `H` is a `ConstantInteriorScalingOperator`. On a
-N-dimensional grid, `H` is the outer product of the 1-dimensional inner
-product operators for each coordinate direction. On a 0-dimensional grid,
-`H` is a 0-dimensional `IdentityTensor`.
+"""
+    inner_product(tg::TensorGrid, stencil_set::StencilSet)
 
-See also: [`ConstantInteriorScalingOperator`](@ref).
+The inner product on the given tensor grid which is the tensor product of the
+individual grids' inner products.
 """
 function inner_product(tg::TensorGrid, stencil_set::StencilSet)
     return mapreduce(g->inner_product(g,stencil_set), ⊗, tg.grids)
 end
 
+"""
+    inner_product(g::EquidistantGrid, stencil_set::StencilSet)
+
+The inner product using weights `H` from `stencil_set`.
+
+See also: [`ConstantInteriorScalingOperator`](@ref).
+"""
 function inner_product(g::EquidistantGrid, stencil_set::StencilSet)
     interior_weight = parse_scalar(stencil_set["H"]["inner"])
     closure_weights = parse_tuple(stencil_set["H"]["closure"])
+    return inner_product(g, interior_weight, closure_weights)
+end
 
+"""
+    inner_product(g::EquidistantGrid, interior_weight, closure_weights)
+
+The inner product on `g` with explicit weights.
+
+See also: [`ConstantInteriorScalingOperator`](@ref).
+"""
+function inner_product(g::EquidistantGrid, interior_weight, closure_weights)
     h = spacing(g)
     return SbpOperators.ConstantInteriorScalingOperator(g, h*interior_weight, h.*closure_weights)
 end
 
+"""
+    inner_product(g::ZeroDimGrid, stencil_set::StencilSet)
+
+The identity tensor with the correct type parameters.
+
+Implemented to simplify 1D code for sbp-operators.
+"""
 inner_product(g::ZeroDimGrid, stencil_set::StencilSet) = IdentityTensor{component_type(g)}()
 
--- a/src/SbpOperators/volumeops/inner_products/inverse_inner_product.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/inner_products/inverse_inner_product.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,27 +1,51 @@
 """
-    inverse_inner_product(grid::EquidistantGrid, interior_weight, closure_weights)
-
-Constructs the inverse inner product operator `H⁻¹` as a `LazyTensor` using
-the weights of `H`, `interior_weight`, `closure_weights`. `H⁻¹` is inverse of
-the inner product operator `H`.
+    inverse_inner_product(grid, ...)
 
-On a 1-dimensional grid, `H⁻¹` is a `ConstantInteriorScalingOperator`. On an
-N-dimensional grid, `H⁻¹` is the outer product of the 1-dimensional inverse
-inner product operators for each coordinate direction. On a 0-dimensional
-`grid`, `H⁻¹` is a 0-dimensional `IdentityTensor`.
+The inverse inner product on a given grid with weights from a stencils set or given
+explicitly.
+"""
+function inverse_inner_product end
 
-See also: [`ConstantInteriorScalingOperator`](@ref).
+"""
+    inverse_inner_product(tg::TensorGrid, stencil_set::StencilSet)
+
+The inverse inner product on the given tensor grid which is the tensor product of the
+individual grids' inner products.
 """
 function inverse_inner_product(tg::TensorGrid, stencil_set::StencilSet)
     return mapreduce(g->inverse_inner_product(g,stencil_set), ⊗, tg.grids)
 end
 
+"""
+    inverse_inner_product(g::EquidistantGrid, stencil_set::StencilSet)
+
+The inverse inner product using weights `H` from `stencil_set`.
+
+See also: [`ConstantInteriorScalingOperator`](@ref).
+"""
 function inverse_inner_product(g::EquidistantGrid, stencil_set::StencilSet)
     interior_weight = parse_scalar(stencil_set["H"]["inner"])
     closure_weights = parse_tuple(stencil_set["H"]["closure"])
+    return inverse_inner_product(g, interior_weight, closure_weights)
+end
 
+"""
+    inverse_inner_product(g::EquidistantGrid, interior_weight, closure_weights)
+
+The inverse inner product on `g` with explicit weights.
+
+See also: [`ConstantInteriorScalingOperator`](@ref).
+"""
+function inverse_inner_product(g::EquidistantGrid, interior_weight, closure_weights)
     h⁻¹ = inverse_spacing(g)
     return SbpOperators.ConstantInteriorScalingOperator(g, h⁻¹*1/interior_weight, h⁻¹./closure_weights)
 end
 
+"""
+    inverse_inner_product(g::ZeroDimGrid, stencil_set::StencilSet)
+
+The identity tensor with the correct type parameters.
+
+Implemented to simplify 1D code for sbp-operators.
+"""
 inverse_inner_product(g::ZeroDimGrid, stencil_set::StencilSet) = IdentityTensor{component_type(g)}()
--- a/src/SbpOperators/volumeops/laplace/laplace.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/laplace/laplace.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,9 +1,8 @@
 """
     Laplace{T, Dim, TM} <: LazyTensor{T, Dim, Dim}
 
-Implements the Laplace operator, approximating ∑d²/xᵢ² , i = 1,...,`Dim` as a
-`LazyTensor`. Additionally `Laplace` stores the `StencilSet`
-used to construct the `LazyTensor`.
+The Laplace operator, approximating ∑d²/xᵢ² , i = 1,...,`Dim` as a
+`LazyTensor`.
 """
 struct Laplace{T, Dim, TM<:LazyTensor{T, Dim, Dim}} <: LazyTensor{T, Dim, Dim}
     D::TM       # Difference operator
@@ -32,7 +31,7 @@
 """
     laplace(g::Grid, stencil_set)
 
-Creates the Laplace operator operator `Δ` as a `LazyTensor` on the given grid
+Creates the Laplace operator operator `Δ` as a `LazyTensor` on the given grid.
 
 `Δ` approximates the Laplace operator ∑d²/xᵢ² , i = 1,...,`Dim` on `g`. The
 approximation depends on the type of grid and the stencil set.
--- a/src/SbpOperators/volumeops/stencil_operator_distinct_closures.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/stencil_operator_distinct_closures.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,7 +1,8 @@
 """
     StencilOperatorDistinctClosures{T,K,N,M,L} <: LazyTensor{T,1}
 
-A one dimensional stencil operator with separate closures for the two boundaries.
+A one dimensional stencil operator with separate closures for the two
+boundaries.
 
 `StencilOperatorDistinctClosures` can be contrasted to `VolumeOperator` in
 that it has different closure stencils for the upper and lower boundary.
--- a/src/SbpOperators/volumeops/volume_operator.jl	Tue May 02 20:14:39 2023 +0200
+++ b/src/SbpOperators/volumeops/volume_operator.jl	Tue May 02 22:09:33 2023 +0200
@@ -1,7 +1,7 @@
 """
     VolumeOperator{T,N,M,K} <: LazyTensor{T,1,1}
 
-Implements a one-dimensional constant coefficients volume operator
+A one-dimensional constant coefficients stencil operator.
 """
 struct VolumeOperator{T,N,M,K} <: LazyTensor{T,1,1}
     inner_stencil::Stencil{T,N}