MOM_transcribe_grid.F90
1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> Module with routines for copying information from a shared dynamic horizontal
6!! grid to an ocean-specific horizontal grid and the reverse.
7module mom_transcribe_grid
8
9use mom_array_transform, only : rotate_array, rotate_array_pair
10use mom_domains, only : pass_var, pass_vector
11use mom_domains, only : to_all, scalar_pair, cgrid_ne, agrid, bgrid_ne, corner
12use mom_dyn_horgrid, only : dyn_horgrid_type, set_derived_dyn_horgrid
13use mom_dyn_horgrid, only : rotate_dyngrid=>rotate_dyn_horgrid
14use mom_error_handler, only : mom_error, mom_mesg, fatal, warning
15use mom_grid, only : ocean_grid_type, set_derived_metrics
16use mom_unit_scaling, only : unit_scale_type
17
18implicit none ; private
19
20public copy_dyngrid_to_mom_grid, copy_mom_grid_to_dyngrid, rotate_dyngrid
21
22contains
23
24!> Copies information from a dynamic (shared) horizontal grid type into an
25!! ocean_grid_type. There may also be a change in the reference
26!! height for topography between the two grids.
27subroutine copy_dyngrid_to_mom_grid(dG, oG, US)
28 type(dyn_horgrid_type), intent(in) :: dg !< Common horizontal grid type
29 type(ocean_grid_type), intent(inout) :: og !< Ocean grid type
30 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
31
32 integer :: isd, ied, jsd, jed ! Common data domains.
33 integer :: isdb, iedb, jsdb, jedb ! Common data domains.
34 integer :: ido, jdo, ido2, jdo2 ! Indexing offsets between the grids.
35 integer :: igst, jgst ! Global starting indices.
36 integer :: i, j
37
38 ! MOM_grid_init and create_dyn_horgrid are called outside of this routine.
39 ! This routine copies over the fields that were set by MOM_initialized_fixed.
40
41 ! Determine the indexing offsets between the grids.
42 ido = dg%idg_offset - og%idg_offset
43 jdo = dg%jdg_offset - og%jdg_offset
44
45 isd = max(og%isd, dg%isd+ido) ; jsd = max(og%jsd, dg%jsd+jdo)
46 ied = min(og%ied, dg%ied+ido) ; jed = min(og%jed, dg%jed+jdo)
47 isdb = max(og%IsdB, dg%IsdB+ido) ; jsdb = max(og%JsdB, dg%JsdB+jdo)
48 iedb = min(og%IedB, dg%IedB+ido) ; jedb = min(og%JedB, dg%JedB+jdo)
49
50 ! Check that the grids conform.
51 if ((isd > og%isc) .or. (ied < og%ied) .or. (jsd > og%jsc) .or. (jed > og%jed)) &
52 call mom_error(fatal, "copy_dyngrid_to_MOM_grid called with incompatible grids.")
53
54 do i=isd,ied ; do j=jsd,jed
55 og%geoLonT(i,j) = dg%geoLonT(i+ido,j+jdo)
56 og%geoLatT(i,j) = dg%geoLatT(i+ido,j+jdo)
57 og%dxT(i,j) = dg%dxT(i+ido,j+jdo)
58 og%dyT(i,j) = dg%dyT(i+ido,j+jdo)
59 og%areaT(i,j) = dg%areaT(i+ido,j+jdo)
60 og%bathyT(i,j) = dg%bathyT(i+ido,j+jdo) - og%Z_ref
61 og%meanSL(i,j) = dg%meanSL(i+ido,j+jdo) + og%Z_ref
62
63 og%dF_dx(i,j) = dg%dF_dx(i+ido,j+jdo)
64 og%dF_dy(i,j) = dg%dF_dy(i+ido,j+jdo)
65 og%sin_rot(i,j) = dg%sin_rot(i+ido,j+jdo)
66 og%cos_rot(i,j) = dg%cos_rot(i+ido,j+jdo)
67 og%mask2dT(i,j) = dg%mask2dT(i+ido,j+jdo)
68 enddo ; enddo
69
70 do i=isdb,iedb ; do j=jsd,jed
71 og%geoLonCu(i,j) = dg%geoLonCu(i+ido,j+jdo)
72 og%geoLatCu(i,j) = dg%geoLatCu(i+ido,j+jdo)
73 og%dxCu(i,j) = dg%dxCu(i+ido,j+jdo)
74 og%dyCu(i,j) = dg%dyCu(i+ido,j+jdo)
75 og%dy_Cu(i,j) = dg%dy_Cu(i+ido,j+jdo)
76
77 og%porous_DminU(i,j) = dg%porous_DminU(i+ido,j+jdo) - og%Z_ref
78 og%porous_DmaxU(i,j) = dg%porous_DmaxU(i+ido,j+jdo) - og%Z_ref
79 og%porous_DavgU(i,j) = dg%porous_DavgU(i+ido,j+jdo) - og%Z_ref
80
81 og%mask2dCu(i,j) = dg%mask2dCu(i+ido,j+jdo)
82 og%OBCmaskCu(i,j) = dg%OBCmaskCu(i+ido,j+jdo)
83 og%areaCu(i,j) = dg%areaCu(i+ido,j+jdo)
84 og%IareaCu(i,j) = dg%IareaCu(i+ido,j+jdo)
85 enddo ; enddo
86
87 do i=isd,ied ; do j=jsdb,jedb
88 og%geoLonCv(i,j) = dg%geoLonCv(i+ido,j+jdo)
89 og%geoLatCv(i,j) = dg%geoLatCv(i+ido,j+jdo)
90 og%dxCv(i,j) = dg%dxCv(i+ido,j+jdo)
91 og%dyCv(i,j) = dg%dyCv(i+ido,j+jdo)
92 og%dx_Cv(i,j) = dg%dx_Cv(i+ido,j+jdo)
93
94 og%porous_DminV(i,j) = dg%porous_DminV(i+ido,j+jdo) - og%Z_ref
95 og%porous_DmaxV(i,j) = dg%porous_DmaxV(i+ido,j+jdo) - og%Z_ref
96 og%porous_DavgV(i,j) = dg%porous_DavgV(i+ido,j+jdo) - og%Z_ref
97
98 og%mask2dCv(i,j) = dg%mask2dCv(i+ido,j+jdo)
99 og%OBCmaskCv(i,j) = dg%OBCmaskCv(i+ido,j+jdo)
100 og%areaCv(i,j) = dg%areaCv(i+ido,j+jdo)
101 og%IareaCv(i,j) = dg%IareaCv(i+ido,j+jdo)
102 enddo ; enddo
103
104 do i=isdb,iedb ; do j=jsdb,jedb
105 og%geoLonBu(i,j) = dg%geoLonBu(i+ido,j+jdo)
106 og%geoLatBu(i,j) = dg%geoLatBu(i+ido,j+jdo)
107 og%dxBu(i,j) = dg%dxBu(i+ido,j+jdo)
108 og%dyBu(i,j) = dg%dyBu(i+ido,j+jdo)
109 og%areaBu(i,j) = dg%areaBu(i+ido,j+jdo)
110 og%CoriolisBu(i,j) = dg%CoriolisBu(i+ido,j+jdo)
111 og%Coriolis2Bu(i,j) = dg%Coriolis2Bu(i+ido,j+jdo)
112 og%mask2dBu(i,j) = dg%mask2dBu(i+ido,j+jdo)
113 enddo ; enddo
114
115 og%bathymetry_at_vel = dg%bathymetry_at_vel
116 if (og%bathymetry_at_vel) then
117 do i=isdb,iedb ; do j=jsd,jed
118 og%Dblock_u(i,j) = dg%Dblock_u(i+ido,j+jdo) - og%Z_ref
119 og%Dopen_u(i,j) = dg%Dopen_u(i+ido,j+jdo) - og%Z_ref
120 enddo ; enddo
121 do i=isd,ied ; do j=jsdb,jedb
122 og%Dblock_v(i,j) = dg%Dblock_v(i+ido,j+jdo) - og%Z_ref
123 og%Dopen_v(i,j) = dg%Dopen_v(i+ido,j+jdo) - og%Z_ref
124 enddo ; enddo
125 endif
126
127 og%gridLonT(og%isg:og%ieg) = dg%gridLonT(dg%isg:dg%ieg)
128 og%gridLatT(og%jsg:og%jeg) = dg%gridLatT(dg%jsg:dg%jeg)
129 ! The more complicated logic here avoids segmentation faults if one grid uses
130 ! global symmetric memory while the other does not. Because a northeast grid
131 ! convention is being used, the upper bounds for each array correspond.
132 ! Note that the dynamic grid always uses symmetric memory.
133 ido2 = dg%IegB-og%IegB ; igst = max(og%IsgB, (dg%isg-1)-ido2)
134 jdo2 = dg%JegB-og%JegB ; jgst = max(og%JsgB, (dg%jsg-1)-jdo2)
135 do i=igst,og%IegB ; og%gridLonB(i) = dg%gridLonB(i+ido2) ; enddo
136 do j=jgst,og%JegB ; og%gridLatB(j) = dg%gridLatB(j+jdo2) ; enddo
137
138 ! Copy various scalar variables and strings.
139 og%x_axis_units = dg%x_axis_units ; og%y_axis_units = dg%y_axis_units
140 og%x_ax_unit_short = dg%x_ax_unit_short ; og%y_ax_unit_short = dg%y_ax_unit_short
141 og%grid_unit_to_L = dg%grid_unit_to_L
142 og%areaT_global = dg%areaT_global ; og%IareaT_global = dg%IareaT_global
143 og%south_lat = dg%south_lat ; og%west_lon = dg%west_lon
144 og%len_lat = dg%len_lat ; og%len_lon = dg%len_lon
145 og%Rad_Earth_L = dg%Rad_Earth_L
146 og%max_depth = dg%max_depth
147
148! Update the halos in case the dynamic grid has smaller halos than the ocean grid.
149 call pass_var(og%areaT, og%Domain)
150 call pass_var(og%bathyT, og%Domain)
151 call pass_var(og%meanSL, og%Domain)
152 call pass_var(og%geoLonT, og%Domain)
153 call pass_var(og%geoLatT, og%Domain)
154 call pass_vector(og%dxT, og%dyT, og%Domain, to_all+scalar_pair, agrid)
155 call pass_vector(og%dF_dx, og%dF_dy, og%Domain, to_all, agrid)
156 call pass_vector(og%cos_rot, og%sin_rot, og%Domain, to_all, agrid)
157 call pass_var(og%mask2dT, og%Domain)
158
159 call pass_vector(og%areaCu, og%areaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
160 call pass_vector(og%dyCu, og%dxCv, og%Domain, to_all+scalar_pair, cgrid_ne)
161 call pass_vector(og%dxCu, og%dyCv, og%Domain, to_all+scalar_pair, cgrid_ne)
162 call pass_vector(og%dy_Cu, og%dx_Cv, og%Domain, to_all+scalar_pair, cgrid_ne)
163 call pass_vector(og%mask2dCu, og%mask2dCv, og%Domain, to_all+scalar_pair, cgrid_ne)
164 call pass_vector(og%OBCmaskCu, og%OBCmaskCv, og%Domain, to_all+scalar_pair, cgrid_ne)
165 call pass_vector(og%IareaCu, og%IareaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
166 call pass_vector(og%IareaCu, og%IareaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
167 call pass_vector(og%geoLatCu, og%geoLatCv, og%Domain, to_all+scalar_pair, cgrid_ne)
168
169 call pass_var(og%areaBu, og%Domain, position=corner)
170 call pass_var(og%geoLonBu, og%Domain, position=corner, inner_halo=og%isc-isd)
171 call pass_var(og%geoLatBu, og%Domain, position=corner)
172 call pass_vector(og%dxBu, og%dyBu, og%Domain, to_all+scalar_pair, bgrid_ne)
173 call pass_var(og%CoriolisBu, og%Domain, position=corner)
174 call pass_var(og%Coriolis2Bu, og%Domain, position=corner)
175 call pass_var(og%mask2dBu, og%Domain, position=corner)
176
177 if (og%bathymetry_at_vel) then
178 call pass_vector(og%Dblock_u, og%Dblock_v, og%Domain, to_all+scalar_pair, cgrid_ne)
179 call pass_vector(og%Dopen_u, og%Dopen_v, og%Domain, to_all+scalar_pair, cgrid_ne)
180 endif
181
182 call set_derived_metrics(og, us)
183
184end subroutine copy_dyngrid_to_mom_grid
185
186
187!> Copies information from an ocean_grid_type into a dynamic (shared)
188!! horizontal grid type. There may also be a change in the reference
189!! height for topography between the two grids.
190subroutine copy_mom_grid_to_dyngrid(oG, dG, US)
191 type(ocean_grid_type), intent(in) :: og !< Ocean grid type
192 type(dyn_horgrid_type), intent(inout) :: dg !< Common horizontal grid type
193 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
194
195 integer :: isd, ied, jsd, jed ! Common data domains.
196 integer :: isdb, iedb, jsdb, jedb ! Common data domains.
197 integer :: ido, jdo, ido2, jdo2 ! Indexing offsets between the grids.
198 integer :: igst, jgst ! Global starting indices.
199 integer :: i, j
200
201 ! MOM_grid_init and create_dyn_horgrid are called outside of this routine.
202 ! This routine copies over the fields that were set by MOM_initialized_fixed.
203
204 ! Determine the indexing offsets between the grids.
205 ido = og%idG_offset - dg%idG_offset
206 jdo = og%jdG_offset - dg%jdG_offset
207
208 isd = max(dg%isd, og%isd+ido) ; jsd = max(dg%jsd, og%jsd+jdo)
209 ied = min(dg%ied, og%ied+ido) ; jed = min(dg%jed, og%jed+jdo)
210 isdb = max(dg%IsdB, og%IsdB+ido) ; jsdb = max(dg%JsdB, og%JsdB+jdo)
211 iedb = min(dg%IedB, og%IedB+ido) ; jedb = min(dg%JedB, og%JedB+jdo)
212
213 ! Check that the grids conform.
214 if ((isd > dg%isc) .or. (ied < dg%ied) .or. (jsd > dg%jsc) .or. (jed > dg%jed)) &
215 call mom_error(fatal, "copy_dyngrid_to_MOM_grid called with incompatible grids.")
216
217 do i=isd,ied ; do j=jsd,jed
218 dg%geoLonT(i,j) = og%geoLonT(i+ido,j+jdo)
219 dg%geoLatT(i,j) = og%geoLatT(i+ido,j+jdo)
220 dg%dxT(i,j) = og%dxT(i+ido,j+jdo)
221 dg%dyT(i,j) = og%dyT(i+ido,j+jdo)
222 dg%areaT(i,j) = og%areaT(i+ido,j+jdo)
223 dg%bathyT(i,j) = og%bathyT(i+ido,j+jdo) + og%Z_ref
224 dg%meanSL(i,j) = og%meanSL(i+ido,j+jdo) - og%Z_ref
225
226 dg%dF_dx(i,j) = og%dF_dx(i+ido,j+jdo)
227 dg%dF_dy(i,j) = og%dF_dy(i+ido,j+jdo)
228 dg%sin_rot(i,j) = og%sin_rot(i+ido,j+jdo)
229 dg%cos_rot(i,j) = og%cos_rot(i+ido,j+jdo)
230 dg%mask2dT(i,j) = og%mask2dT(i+ido,j+jdo)
231 enddo ; enddo
232
233 do i=isdb,iedb ; do j=jsd,jed
234 dg%geoLonCu(i,j) = og%geoLonCu(i+ido,j+jdo)
235 dg%geoLatCu(i,j) = og%geoLatCu(i+ido,j+jdo)
236 dg%dxCu(i,j) = og%dxCu(i+ido,j+jdo)
237 dg%dyCu(i,j) = og%dyCu(i+ido,j+jdo)
238 dg%dy_Cu(i,j) = og%dy_Cu(i+ido,j+jdo)
239
240 dg%porous_DminU(i,j) = og%porous_DminU(i+ido,j+jdo) + og%Z_ref
241 dg%porous_DmaxU(i,j) = og%porous_DmaxU(i+ido,j+jdo) + og%Z_ref
242 dg%porous_DavgU(i,j) = og%porous_DavgU(i+ido,j+jdo) + og%Z_ref
243
244 dg%mask2dCu(i,j) = og%mask2dCu(i+ido,j+jdo)
245 dg%OBCmaskCu(i,j) = og%OBCmaskCu(i+ido,j+jdo)
246 dg%areaCu(i,j) = og%areaCu(i+ido,j+jdo)
247 dg%IareaCu(i,j) = og%IareaCu(i+ido,j+jdo)
248 enddo ; enddo
249
250 do i=isd,ied ; do j=jsdb,jedb
251 dg%geoLonCv(i,j) = og%geoLonCv(i+ido,j+jdo)
252 dg%geoLatCv(i,j) = og%geoLatCv(i+ido,j+jdo)
253 dg%dxCv(i,j) = og%dxCv(i+ido,j+jdo)
254 dg%dyCv(i,j) = og%dyCv(i+ido,j+jdo)
255 dg%dx_Cv(i,j) = og%dx_Cv(i+ido,j+jdo)
256
257 dg%porous_DminV(i,j) = og%porous_DminU(i+ido,j+jdo) + og%Z_ref
258 dg%porous_DmaxV(i,j) = og%porous_DmaxU(i+ido,j+jdo) + og%Z_ref
259 dg%porous_DavgV(i,j) = og%porous_DavgU(i+ido,j+jdo) + og%Z_ref
260
261 dg%mask2dCv(i,j) = og%mask2dCv(i+ido,j+jdo)
262 dg%OBCmaskCv(i,j) = og%OBCmaskCv(i+ido,j+jdo)
263 dg%areaCv(i,j) = og%areaCv(i+ido,j+jdo)
264 dg%IareaCv(i,j) = og%IareaCv(i+ido,j+jdo)
265 enddo ; enddo
266
267 do i=isdb,iedb ; do j=jsdb,jedb
268 dg%geoLonBu(i,j) = og%geoLonBu(i+ido,j+jdo)
269 dg%geoLatBu(i,j) = og%geoLatBu(i+ido,j+jdo)
270 dg%dxBu(i,j) = og%dxBu(i+ido,j+jdo)
271 dg%dyBu(i,j) = og%dyBu(i+ido,j+jdo)
272 dg%areaBu(i,j) = og%areaBu(i+ido,j+jdo)
273 dg%CoriolisBu(i,j) = og%CoriolisBu(i+ido,j+jdo)
274 dg%Coriolis2Bu(i,j) = og%Coriolis2Bu(i+ido,j+jdo)
275 dg%mask2dBu(i,j) = og%mask2dBu(i+ido,j+jdo)
276 enddo ; enddo
277
278 dg%bathymetry_at_vel = og%bathymetry_at_vel
279 if (dg%bathymetry_at_vel) then
280 do i=isdb,iedb ; do j=jsd,jed
281 dg%Dblock_u(i,j) = og%Dblock_u(i+ido,j+jdo) + og%Z_ref
282 dg%Dopen_u(i,j) = og%Dopen_u(i+ido,j+jdo) + og%Z_ref
283 enddo ; enddo
284 do i=isd,ied ; do j=jsdb,jedb
285 dg%Dblock_v(i,j) = og%Dblock_v(i+ido,j+jdo) + og%Z_ref
286 dg%Dopen_v(i,j) = og%Dopen_v(i+ido,j+jdo) + og%Z_ref
287 enddo ; enddo
288 endif
289
290 dg%gridLonT(dg%isg:dg%ieg) = og%gridLonT(og%isg:og%ieg)
291 dg%gridLatT(dg%jsg:dg%jeg) = og%gridLatT(og%jsg:og%jeg)
292
293 ! The more complicated logic here avoids segmentation faults if one grid uses
294 ! global symmetric memory while the other does not. Because a northeast grid
295 ! convention is being used, the upper bounds for each array correspond.
296 ! Note that the dynamic grid always uses symmetric memory.
297 ido2 = og%IegB-dg%IegB ; igst = max(dg%isg-1, og%IsgB-ido2)
298 jdo2 = og%JegB-dg%JegB ; jgst = max(dg%jsg-1, og%JsgB-jdo2)
299 do i=igst,dg%IegB ; dg%gridLonB(i) = og%gridLonB(i+ido2) ; enddo
300 do j=jgst,dg%JegB ; dg%gridLatB(j) = og%gridLatB(j+jdo2) ; enddo
301
302 ! Copy various scalar variables and strings.
303 dg%x_axis_units = og%x_axis_units ; dg%y_axis_units = og%y_axis_units
304 dg%x_ax_unit_short = og%x_ax_unit_short ; dg%y_ax_unit_short = og%y_ax_unit_short
305 dg%grid_unit_to_L = og%grid_unit_to_L
306 dg%areaT_global = og%areaT_global ; dg%IareaT_global = og%IareaT_global
307 dg%south_lat = og%south_lat ; dg%west_lon = og%west_lon
308 dg%len_lat = og%len_lat ; dg%len_lon = og%len_lon
309 dg%Rad_Earth_L = og%Rad_Earth_L
310 dg%max_depth = og%max_depth
311
312! Update the halos in case the dynamic grid has smaller halos than the ocean grid.
313 call pass_var(dg%areaT, dg%Domain)
314 call pass_var(dg%bathyT, dg%Domain)
315 call pass_var(dg%meanSL, dg%Domain)
316 call pass_var(dg%geoLonT, dg%Domain)
317 call pass_var(dg%geoLatT, dg%Domain)
318 call pass_vector(dg%dxT, dg%dyT, dg%Domain, to_all+scalar_pair, agrid)
319 call pass_vector(dg%dF_dx, dg%dF_dy, dg%Domain, to_all, agrid)
320 call pass_vector(dg%cos_rot, dg%sin_rot, dg%Domain, to_all, agrid)
321 call pass_var(dg%mask2dT, dg%Domain)
322
323 call pass_vector(dg%areaCu, dg%areaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
324 call pass_vector(dg%dyCu, dg%dxCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
325 call pass_vector(dg%dxCu, dg%dyCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
326 call pass_vector(dg%dy_Cu, dg%dx_Cv, dg%Domain, to_all+scalar_pair, cgrid_ne)
327 call pass_vector(dg%mask2dCu, dg%mask2dCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
328 call pass_vector(dg%OBCmaskCu, dg%OBCmaskCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
329 call pass_vector(dg%IareaCu, dg%IareaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
330 call pass_vector(dg%IareaCu, dg%IareaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
331 call pass_vector(dg%geoLatCu, dg%geoLatCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
332
333 call pass_var(dg%areaBu, dg%Domain, position=corner)
334 call pass_var(dg%geoLonBu, dg%Domain, position=corner, inner_halo=dg%isc-isd)
335 call pass_var(dg%geoLatBu, dg%Domain, position=corner)
336 call pass_vector(dg%dxBu, dg%dyBu, dg%Domain, to_all+scalar_pair, bgrid_ne)
337 call pass_var(dg%CoriolisBu, dg%Domain, position=corner)
338 call pass_var(dg%Coriolis2Bu, dg%Domain, position=corner)
339 call pass_var(dg%mask2dBu, dg%Domain, position=corner)
340
341 if (dg%bathymetry_at_vel) then
342 call pass_vector(dg%Dblock_u, dg%Dblock_v, dg%Domain, to_all+scalar_pair, cgrid_ne)
343 call pass_vector(dg%Dopen_u, dg%Dopen_v, dg%Domain, to_all+scalar_pair, cgrid_ne)
344 endif
345
346 call set_derived_dyn_horgrid(dg, us)
347
348end subroutine copy_mom_grid_to_dyngrid
349
350end module mom_transcribe_grid