comparison src/Grids/EquidistantGrid.jl @ 912:f800f97b3a45 feature/variable_derivatives

Merge default
author Jonatan Werpers <jonatan@werpers.com>
date Tue, 15 Feb 2022 15:15:11 +0100
parents b900fea1c057
children 9b40aeac4269 5b3d4a8ec3ab
comparison
equal deleted inserted replaced
906:33c7e266e1a9 912:f800f97b3a45
3 export inverse_spacing 3 export inverse_spacing
4 export restrict 4 export restrict
5 export boundary_identifiers 5 export boundary_identifiers
6 export boundary_grid 6 export boundary_grid
7 export refine 7 export refine
8 export coarsen
8 9
9 """ 10 """
10 EquidistantGrid(size::NTuple{Dim, Int}, limit_lower::NTuple{Dim, T}, limit_upper::NTuple{Dim, T}) 11 EquidistantGrid{Dim,T<:Real} <: AbstractGrid
11 EquidistantGrid{T}()
12 12
13 `EquidistantGrid` is a grid with equidistant grid spacing per coordinat direction. 13 `Dim`-dimensional equidistant grid with coordinates of type `T`.
14
15 `EquidistantGrid(size, limit_lower, limit_upper)` construct the grid with the
16 domain defined by the two points P1, and P2 given by `limit_lower` and
17 `limit_upper`. The length of the domain sides are given by the components of
18 (P2-P1). E.g for a 2D grid with P1=(-1,0) and P2=(1,2) the domain is defined
19 as (-1,1)x(0,2). The side lengths of the grid are not allowed to be negative.
20 The number of equidistantly spaced points in each coordinate direction are given
21 by `size`.
22
23 `EquidistantGrid{T}()` constructs a 0-dimensional grid.
24
25 """ 14 """
26 struct EquidistantGrid{Dim,T<:Real} <: AbstractGrid 15 struct EquidistantGrid{Dim,T<:Real} <: AbstractGrid
27 size::NTuple{Dim, Int} 16 size::NTuple{Dim, Int}
28 limit_lower::NTuple{Dim, T} 17 limit_lower::NTuple{Dim, T}
29 limit_upper::NTuple{Dim, T} 18 limit_upper::NTuple{Dim, T}
30 19
31 # General constructor
32 function EquidistantGrid{Dim,T}(size::NTuple{Dim, Int}, limit_lower::NTuple{Dim, T}, limit_upper::NTuple{Dim, T}) where {Dim,T} 20 function EquidistantGrid{Dim,T}(size::NTuple{Dim, Int}, limit_lower::NTuple{Dim, T}, limit_upper::NTuple{Dim, T}) where {Dim,T}
33 if any(size .<= 0) 21 if any(size .<= 0)
34 throw(DomainError("all components of size must be postive")) 22 throw(DomainError("all components of size must be postive"))
35 end 23 end
36 if any(limit_upper.-limit_lower .<= 0) 24 if any(limit_upper.-limit_lower .<= 0)
38 end 26 end
39 return new{Dim,T}(size, limit_lower, limit_upper) 27 return new{Dim,T}(size, limit_lower, limit_upper)
40 end 28 end
41 end 29 end
42 30
31 """
32 EquidistantGrid(size, limit_lower, limit_upper)
33
34 Construct an equidistant grid with corners at the coordinates `limit_lower` and
35 `limit_upper`.
36
37 The length of the domain sides are given by the components of
38 `limit_upper-limit_lower`. E.g for a 2D grid with `limit_lower=(-1,0)` and `limit_upper=(1,2)` the domain is defined
39 as `(-1,1)x(0,2)`. The side lengths of the grid are not allowed to be negative.
40
41 The number of equidistantly spaced points in each coordinate direction are given
42 by the tuple `size`.
43 """
43 function EquidistantGrid(size, limit_lower, limit_upper) 44 function EquidistantGrid(size, limit_lower, limit_upper)
44 return EquidistantGrid{length(size), eltype(limit_lower)}(size, limit_lower, limit_upper) 45 return EquidistantGrid{length(size), eltype(limit_lower)}(size, limit_lower, limit_upper)
45 end 46 end
47
48 """
49 EquidistantGrid{T}()
50
51 Constructs a 0-dimensional grid.
52 """
46 EquidistantGrid{T}() where T = EquidistantGrid{0,T}((),(),()) # Convenience constructor for 0-dim grid 53 EquidistantGrid{T}() where T = EquidistantGrid{0,T}((),(),()) # Convenience constructor for 0-dim grid
47 54
48 """ 55 """
49 EquidistantGrid(size::Int, limit_lower::T, limit_upper::T) 56 EquidistantGrid(size::Int, limit_lower::T, limit_upper::T)
50 57
68 dimension(grid::EquidistantGrid{Dim}) where Dim = Dim 75 dimension(grid::EquidistantGrid{Dim}) where Dim = Dim
69 76
70 """ 77 """
71 spacing(grid::EquidistantGrid) 78 spacing(grid::EquidistantGrid)
72 79
73 The spacing between the grid points of the grid. 80 The spacing between grid points.
74 """ 81 """
75 spacing(grid::EquidistantGrid) = (grid.limit_upper.-grid.limit_lower)./(grid.size.-1) 82 spacing(grid::EquidistantGrid) = (grid.limit_upper.-grid.limit_lower)./(grid.size.-1)
76 83
77 """ 84 """
78 inverse_spacing(grid::EquidistantGrid) 85 inverse_spacing(grid::EquidistantGrid)
79 86
80 The reciprocal of the spacing between the grid points of the grid. 87 The reciprocal of the spacing between grid points.
81 """ 88 """
82 inverse_spacing(grid::EquidistantGrid) = 1 ./ spacing(grid) 89 inverse_spacing(grid::EquidistantGrid) = 1 ./ spacing(grid)
83 90
84 """ 91 """
85 points(grid::EquidistantGrid) 92 points(grid::EquidistantGrid)
97 end 104 end
98 105
99 """ 106 """
100 restrict(::EquidistantGrid, dim) 107 restrict(::EquidistantGrid, dim)
101 108
102 Pick out given dimensions from the grid and return a grid for them 109 Pick out given dimensions from the grid and return a grid for them.
103 """ 110 """
104 function restrict(grid::EquidistantGrid, dim) 111 function restrict(grid::EquidistantGrid, dim)
105 size = grid.size[dim] 112 size = grid.size[dim]
106 limit_lower = grid.limit_lower[dim] 113 limit_lower = grid.limit_lower[dim]
107 limit_upper = grid.limit_upper[dim] 114 limit_upper = grid.limit_upper[dim]
143 """ 150 """
144 refine(grid::EquidistantGrid, r::Int) 151 refine(grid::EquidistantGrid, r::Int)
145 152
146 Refines `grid` by a factor `r`. The factor is applied to the number of 153 Refines `grid` by a factor `r`. The factor is applied to the number of
147 intervals which is 1 less than the size of the grid. 154 intervals which is 1 less than the size of the grid.
155
156 See also: [`coarsen`](@ref)
148 """ 157 """
149 function refine(grid::EquidistantGrid, r::Int) 158 function refine(grid::EquidistantGrid, r::Int)
150 sz = size(grid) 159 sz = size(grid)
151 new_sz = (sz .- 1).*r .+ 1 160 new_sz = (sz .- 1).*r .+ 1
152 return EquidistantGrid{dimension(grid), eltype(grid)}(new_sz, grid.limit_lower, grid.limit_upper) 161 return EquidistantGrid{dimension(grid), eltype(grid)}(new_sz, grid.limit_lower, grid.limit_upper)
153 end 162 end
163
164 """
165 coarsen(grid::EquidistantGrid, r::Int)
166
167 Coarsens `grid` by a factor `r`. The factor is applied to the number of
168 intervals which is 1 less than the size of the grid. If the number of
169 intervals are not divisible by `r` an error is raised.
170
171 See also: [`refine`](@ref)
172 """
173 function coarsen(grid::EquidistantGrid, r::Int)
174 sz = size(grid)
175
176 if !all(n -> (n % r == 0), sz.-1)
177 throw(DomainError(r, "Size minus 1 must be divisible by the ratio."))
178 end
179
180 new_sz = (sz .- 1).÷r .+ 1
181
182 return EquidistantGrid{dimension(grid), eltype(grid)}(new_sz, grid.limit_lower, grid.limit_upper)
183 end