|
| 1 | +module OceananigansXESMFExt |
| 2 | + |
| 3 | +using XESMF |
| 4 | +using Oceananigans |
| 5 | +using Oceananigans.Architectures: architecture, on_architecture |
| 6 | +using Oceananigans.Fields: AbstractField, topology, location |
| 7 | +using Oceananigans.Grids: λnodes, φnodes, Center, Face, total_length |
| 8 | + |
| 9 | +import Oceananigans.Fields: regrid! |
| 10 | +import XESMF: Regridder, extract_xesmf_coordinates_structure |
| 11 | + |
| 12 | +function x_node_array(x::AbstractVector, Nx, Ny) |
| 13 | + return Array(repeat(view(x, 1:Nx), 1, Ny))' |
| 14 | +end |
| 15 | +function y_node_array(x::AbstractVector, Nx, Ny) |
| 16 | + return Array(repeat(view(x, 1:Ny)', Nx, 1))' |
| 17 | +end |
| 18 | +x_node_array(x::AbstractMatrix, Nx, Ny) = Array(view(x, 1:Nx, 1:Ny))' |
| 19 | + |
| 20 | +function x_vertex_array(x::AbstractVector, Nx, Ny) |
| 21 | + return Array(repeat(view(x, 1:Nx+1), 1, Ny+1))' |
| 22 | +end |
| 23 | +function y_vertex_array(x::AbstractVector, Nx, Ny) |
| 24 | + return Array(repeat(view(x, 1:Ny+1)', Nx+1, 1))' |
| 25 | +end |
| 26 | +x_vertex_array(x::AbstractMatrix, Nx, Ny) = Array(view(x, 1:Nx+1, 1:Ny+1))' |
| 27 | + |
| 28 | +y_node_array(x::AbstractMatrix, Nx, Ny) = x_node_array(x, Nx, Ny) |
| 29 | +y_vertex_array(x::AbstractMatrix, Nx, Ny) = x_vertex_array(x, Nx, Ny) |
| 30 | + |
| 31 | +function extract_xesmf_coordinates_structure(dst_field::AbstractField, src_field::AbstractField) |
| 32 | + |
| 33 | + ℓx, ℓy, ℓz = Oceananigans.Fields.instantiated_location(src_field) |
| 34 | + |
| 35 | + dst_grid = dst_field.grid |
| 36 | + src_grid = src_field.grid |
| 37 | + |
| 38 | + # Extract center coordinates from both fields |
| 39 | + λᵈ = λnodes(dst_grid, Center(), Center(), ℓz, with_halos=true) |
| 40 | + φᵈ = φnodes(dst_grid, Center(), Center(), ℓz, with_halos=true) |
| 41 | + λˢ = λnodes(src_grid, Center(), Center(), ℓz, with_halos=true) |
| 42 | + φˢ = φnodes(src_grid, Center(), Center(), ℓz, with_halos=true) |
| 43 | + |
| 44 | + # Extract cell vertices |
| 45 | + λvᵈ = λnodes(dst_grid, Face(), Face(), ℓz, with_halos=true) |
| 46 | + φvᵈ = φnodes(dst_grid, Face(), Face(), ℓz, with_halos=true) |
| 47 | + λvˢ = λnodes(src_grid, Face(), Face(), ℓz, with_halos=true) |
| 48 | + φvˢ = φnodes(src_grid, Face(), Face(), ℓz, with_halos=true) |
| 49 | + |
| 50 | + # Build data structures expected by xESMF |
| 51 | + Nˢx, Nˢy, Nˢz = size(src_field) |
| 52 | + Nᵈx, Nᵈy, Nᵈz = size(dst_field) |
| 53 | + |
| 54 | + λᵈ = x_node_array(λᵈ, Nᵈx, Nᵈy) |
| 55 | + φᵈ = y_node_array(φᵈ, Nᵈx, Nᵈy) |
| 56 | + λˢ = x_node_array(λˢ, Nˢx, Nˢy) |
| 57 | + φˢ = y_node_array(φˢ, Nˢx, Nˢy) |
| 58 | + |
| 59 | + λvᵈ = x_vertex_array(λvᵈ, Nᵈx, Nᵈy) |
| 60 | + φvᵈ = y_vertex_array(φvᵈ, Nᵈx, Nᵈy) |
| 61 | + λvˢ = x_vertex_array(λvˢ, Nˢx, Nˢy) |
| 62 | + φvˢ = y_vertex_array(φvˢ, Nˢx, Nˢy) |
| 63 | + |
| 64 | + dst_coordinates = Dict("lat" => φᵈ, # φ is latitude |
| 65 | + "lon" => λᵈ, # λ is longitude |
| 66 | + "lat_b" => φvᵈ, |
| 67 | + "lon_b" => λvᵈ) |
| 68 | + |
| 69 | + src_coordinates = Dict("lat" => φˢ, # φ is latitude |
| 70 | + "lon" => λˢ, # λ is longitude |
| 71 | + "lat_b" => φvˢ, |
| 72 | + "lon_b" => λvˢ) |
| 73 | + |
| 74 | + return dst_coordinates, src_coordinates |
| 75 | +end |
| 76 | + |
| 77 | +""" |
| 78 | + Regridder(dst_field::AbstractField, src_field::AbstractField; method="conservative") |
| 79 | +
|
| 80 | +Return a regridder from `src_field` to `dst_field` using the specified `method`. |
| 81 | +The regridder contains a sparse matrix with the regridding weights. |
| 82 | +The regridding weights are obtained via xESMF Python package. |
| 83 | +xESMF exposes five different regridding algorithms from the ESMF library, |
| 84 | +specified with the `method` keyword argument: |
| 85 | +
|
| 86 | +* `"bilinear"`: `ESMF.RegridMethod.BILINEAR` |
| 87 | +* `"conservative"`: `ESMF.RegridMethod.CONSERVE` |
| 88 | +* `"conservative_normed"`: `ESMF.RegridMethod.CONSERVE` |
| 89 | +* `"patch"`: `ESMF.RegridMethod.PATCH` |
| 90 | +* `"nearest_s2d"`: `ESMF.RegridMethod.NEAREST_STOD` |
| 91 | +* `"nearest_d2s"`: `ESMF.RegridMethod.NEAREST_DTOS` |
| 92 | +
|
| 93 | +where `conservative_normed` is just the conservative method with the normalization set to |
| 94 | +`ESMF.NormType.FRACAREA` instead of the default `norm_type = ESMF.NormType.DSTAREA`. |
| 95 | +
|
| 96 | +For more information, see the Python xESMF documentation at: |
| 97 | +
|
| 98 | +> https://xesmf.readthedocs.io/en/latest/notebooks/Compare_algorithms.html |
| 99 | +
|
| 100 | +Example |
| 101 | +======= |
| 102 | +
|
| 103 | +```@example |
| 104 | +using Oceananigans |
| 105 | +using XESMF |
| 106 | +
|
| 107 | +z = (-1, 0) |
| 108 | +tg = TripolarGrid(; size=(360, 170, 1), z, southernmost_latitude = -80) |
| 109 | +llg = LatitudeLongitudeGrid(; size=(360, 180, 1), z, |
| 110 | + longitude=(0, 360), latitude=(-82, 90)) |
| 111 | +
|
| 112 | +src_field = CenterField(tg) |
| 113 | +dst_field = CenterField(llg) |
| 114 | +
|
| 115 | +regridder = Oceananigans.Fields.Regridder(dst_field, src_field, method="conservative") |
| 116 | +``` |
| 117 | +""" |
| 118 | +function Regridder(dst_field::AbstractField, src_field::AbstractField; method="conservative") |
| 119 | + |
| 120 | + ℓx, ℓy, ℓz = Oceananigans.Fields.instantiated_location(src_field) |
| 121 | + |
| 122 | + # We only support regridding between centered fields |
| 123 | + @assert ℓx isa Center |
| 124 | + @assert ℓy isa Center |
| 125 | + @assert (ℓx, ℓy, ℓz) == Oceananigans.Fields.instantiated_location(dst_field) |
| 126 | + |
| 127 | + src_Nz = size(src_field)[3] |
| 128 | + dst_Nz = size(dst_field)[3] |
| 129 | + @assert src_field.grid.z.cᵃᵃᶠ[1:src_Nz+1] == dst_field.grid.z.cᵃᵃᶠ[1:dst_Nz+1] |
| 130 | + |
| 131 | + dst_coordinates, src_coordinates = extract_xesmf_coordinates_structure(dst_field, src_field) |
| 132 | + periodic = Oceananigans.Grids.topology(src_field.grid, 1) === Periodic ? true : false |
| 133 | + |
| 134 | + regridder = XESMF.Regridder(src_coordinates, dst_coordinates; method, periodic) |
| 135 | + weights = regridder.weights |
| 136 | + |
| 137 | + arch = architecture(src_field) |
| 138 | + |
| 139 | + weights = on_architecture(arch, weights) |
| 140 | + |
| 141 | + temp_src = on_architecture(architecture(src_field), regridder.src_temp) |
| 142 | + temp_dst = on_architecture(architecture(dst_field), regridder.dst_temp) |
| 143 | + |
| 144 | + return XESMF.Regridder(method, weights, temp_src, temp_dst) |
| 145 | +end |
| 146 | + |
| 147 | +""" |
| 148 | + regrid!(dst_field, regrider::XESMF.Regridder, src_field) |
| 149 | +
|
| 150 | +Regrid `src_field` onto the grid of field `dst_field` using the regrider `r`. |
| 151 | +
|
| 152 | +Example |
| 153 | +======= |
| 154 | +
|
| 155 | +```@example |
| 156 | +using Oceananigans |
| 157 | +using XESMF |
| 158 | +
|
| 159 | +z = (-1, 0) |
| 160 | +
|
| 161 | +tg = TripolarGrid(; size=(360, 170, 1), z, southernmost_latitude = -80) |
| 162 | +
|
| 163 | +llg = LatitudeLongitudeGrid(; size=(360, 180, 1), z, |
| 164 | + longitude=(0, 360), latitude=(-82, 90)) |
| 165 | +
|
| 166 | +src_field = CenterField(tg) |
| 167 | +dst_field = CenterField(llg) |
| 168 | +
|
| 169 | +λ₀, φ₀ = 150, 30. # degrees |
| 170 | +width = 12 # degrees |
| 171 | +set!(src_field, (λ, φ, z) -> exp(-((λ - λ₀)^2 + (φ - φ₀)^2) / 2width^2)) |
| 172 | +
|
| 173 | +regridder = XESMF.Regridder(dst_field, src_field, method="conservative") |
| 174 | +
|
| 175 | +regrid!(dst_field, regridder, src_field) |
| 176 | +
|
| 177 | +first(Field(Integral(dst_field, dims=(1, 2)))) |
| 178 | +``` |
| 179 | +""" |
| 180 | +function regrid!(dst_field, regridder::XESMF.Regridder, src_field) |
| 181 | + Nz = size(src_field.grid)[3] |
| 182 | + topo_z = topology(src_field)[3]() |
| 183 | + ℓz = location(src_field)[3]() |
| 184 | + |
| 185 | + for k in 1:total_length(ℓz, topo_z, Nz) |
| 186 | + src = vec(interior(src_field, :, :, k)) |
| 187 | + dst = vec(interior(dst_field, :, :, k)) |
| 188 | + regridder(dst, src) |
| 189 | + end |
| 190 | + |
| 191 | + return dst_field |
| 192 | +end |
| 193 | + |
| 194 | +end # module |
0 commit comments