8 use mom_domains,
only : to_all, scalar_pair, cgrid_ne, agrid, bgrid_ne, corner
13 implicit none ;
private
15 public copy_dyngrid_to_mom_grid, copy_mom_grid_to_dyngrid
21 subroutine copy_dyngrid_to_mom_grid(dG, oG)
25 integer :: isd, ied, jsd, jed
26 integer :: isdb, iedb, jsdb, jedb
27 integer :: ido, jdo, ido2, jdo2
35 ido = dg%idg_offset - og%idg_offset
36 jdo = dg%jdg_offset - og%jdg_offset
38 isd = max(og%isd, dg%isd+ido) ; jsd = max(og%jsd, dg%jsd+jdo)
39 ied = min(og%ied, dg%ied+ido) ; jed = min(og%jed, dg%jed+jdo)
40 isdb = max(og%IsdB, dg%IsdB+ido) ; jsdb = max(og%JsdB, dg%JsdB+jdo)
41 iedb = min(og%IedB, dg%IedB+ido) ; jedb = min(og%JedB, dg%JedB+jdo)
44 if ((isd > og%isc) .or. (ied < og%ied) .or. (jsd > og%jsc) .or. (jed > og%jed)) &
45 call mom_error(fatal,
"copy_dyngrid_to_MOM_grid called with incompatible grids.")
47 do i=isd,ied ;
do j=jsd,jed
48 og%geoLonT(i,j) = dg%geoLonT(i+ido,j+jdo)
49 og%geoLatT(i,j) = dg%geoLatT(i+ido,j+jdo)
50 og%dxT(i,j) = dg%dxT(i+ido,j+jdo)
51 og%dyT(i,j) = dg%dyT(i+ido,j+jdo)
52 og%areaT(i,j) = dg%areaT(i+ido,j+jdo)
53 og%bathyT(i,j) = dg%bathyT(i+ido,j+jdo)
55 og%dF_dx(i,j) = dg%dF_dx(i+ido,j+jdo)
56 og%dF_dy(i,j) = dg%dF_dy(i+ido,j+jdo)
57 og%sin_rot(i,j) = dg%sin_rot(i+ido,j+jdo)
58 og%cos_rot(i,j) = dg%cos_rot(i+ido,j+jdo)
59 og%mask2dT(i,j) = dg%mask2dT(i+ido,j+jdo)
62 do i=isdb,iedb ;
do j=jsd,jed
63 og%geoLonCu(i,j) = dg%geoLonCu(i+ido,j+jdo)
64 og%geoLatCu(i,j) = dg%geoLatCu(i+ido,j+jdo)
65 og%dxCu(i,j) = dg%dxCu(i+ido,j+jdo)
66 og%dyCu(i,j) = dg%dyCu(i+ido,j+jdo)
67 og%dy_Cu(i,j) = dg%dy_Cu(i+ido,j+jdo)
69 og%mask2dCu(i,j) = dg%mask2dCu(i+ido,j+jdo)
70 og%areaCu(i,j) = dg%areaCu(i+ido,j+jdo)
71 og%IareaCu(i,j) = dg%IareaCu(i+ido,j+jdo)
74 do i=isd,ied ;
do j=jsdb,jedb
75 og%geoLonCv(i,j) = dg%geoLonCv(i+ido,j+jdo)
76 og%geoLatCv(i,j) = dg%geoLatCv(i+ido,j+jdo)
77 og%dxCv(i,j) = dg%dxCv(i+ido,j+jdo)
78 og%dyCv(i,j) = dg%dyCv(i+ido,j+jdo)
79 og%dx_Cv(i,j) = dg%dx_Cv(i+ido,j+jdo)
81 og%mask2dCv(i,j) = dg%mask2dCv(i+ido,j+jdo)
82 og%areaCv(i,j) = dg%areaCv(i+ido,j+jdo)
83 og%IareaCv(i,j) = dg%IareaCv(i+ido,j+jdo)
86 do i=isdb,iedb ;
do j=jsdb,jedb
87 og%geoLonBu(i,j) = dg%geoLonBu(i+ido,j+jdo)
88 og%geoLatBu(i,j) = dg%geoLatBu(i+ido,j+jdo)
89 og%dxBu(i,j) = dg%dxBu(i+ido,j+jdo)
90 og%dyBu(i,j) = dg%dyBu(i+ido,j+jdo)
91 og%areaBu(i,j) = dg%areaBu(i+ido,j+jdo)
92 og%CoriolisBu(i,j) = dg%CoriolisBu(i+ido,j+jdo)
93 og%mask2dBu(i,j) = dg%mask2dBu(i+ido,j+jdo)
96 og%bathymetry_at_vel = dg%bathymetry_at_vel
97 if (og%bathymetry_at_vel)
then
98 do i=isdb,iedb ;
do j=jsd,jed
99 og%Dblock_u(i,j) = dg%Dblock_u(i+ido,j+jdo)
100 og%Dopen_u(i,j) = dg%Dopen_u(i+ido,j+jdo)
102 do i=isd,ied ;
do j=jsdb,jedb
103 og%Dblock_v(i,j) = dg%Dblock_v(i+ido,j+jdo)
104 og%Dopen_v(i,j) = dg%Dopen_v(i+ido,j+jdo)
108 og%gridLonT(og%isg:og%ieg) = dg%gridLonT(dg%isg:dg%ieg)
109 og%gridLatT(og%jsg:og%jeg) = dg%gridLatT(dg%jsg:dg%jeg)
114 ido2 = dg%IegB-og%IegB ; igst = max(og%IsgB, (dg%isg-1)-ido2)
115 jdo2 = dg%JegB-og%JegB ; jgst = max(og%JsgB, (dg%jsg-1)-jdo2)
116 do i=igst,og%IegB ; og%gridLonB(i) = dg%gridLonB(i+ido2) ;
enddo
117 do j=jgst,og%JegB ; og%gridLatB(j) = dg%gridLatB(j+jdo2) ;
enddo
120 og%x_axis_units = dg%x_axis_units ; og%y_axis_units = dg%y_axis_units
121 og%areaT_global = dg%areaT_global ; og%IareaT_global = dg%IareaT_global
122 og%south_lat = dg%south_lat ; og%west_lon = dg%west_lon
123 og%len_lat = dg%len_lat ; og%len_lon = dg%len_lon
124 og%Rad_Earth = dg%Rad_Earth ; og%max_depth = dg%max_depth
129 call pass_var(og%geoLonT, og%Domain)
130 call pass_var(og%geoLatT, og%Domain)
131 call pass_vector(og%dxT, og%dyT, og%Domain, to_all+scalar_pair, agrid)
132 call pass_vector(og%dF_dx, og%dF_dy, og%Domain, to_all, agrid)
133 call pass_vector(og%cos_rot, og%sin_rot, og%Domain, to_all, agrid)
134 call pass_var(og%mask2dT, og%Domain)
136 call pass_vector(og%areaCu, og%areaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
137 call pass_vector(og%dyCu, og%dxCv, og%Domain, to_all+scalar_pair, cgrid_ne)
138 call pass_vector(og%dxCu, og%dyCv, og%Domain, to_all+scalar_pair, cgrid_ne)
139 call pass_vector(og%dy_Cu, og%dx_Cv, og%Domain, to_all+scalar_pair, cgrid_ne)
140 call pass_vector(og%mask2dCu, og%mask2dCv, og%Domain, to_all+scalar_pair, cgrid_ne)
141 call pass_vector(og%IareaCu, og%IareaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
142 call pass_vector(og%IareaCu, og%IareaCv, og%Domain, to_all+scalar_pair, cgrid_ne)
143 call pass_vector(og%geoLatCu, og%geoLatCv, og%Domain, to_all+scalar_pair, cgrid_ne)
145 call pass_var(og%areaBu, og%Domain, position=corner)
146 call pass_var(og%geoLonBu, og%Domain, position=corner, inner_halo=og%isc-isd)
147 call pass_var(og%geoLatBu, og%Domain, position=corner)
148 call pass_vector(og%dxBu, og%dyBu, og%Domain, to_all+scalar_pair, bgrid_ne)
149 call pass_var(og%CoriolisBu, og%Domain, position=corner)
150 call pass_var(og%mask2dBu, og%Domain, position=corner)
152 if (og%bathymetry_at_vel)
then
153 call pass_vector(og%Dblock_u, og%Dblock_v, og%Domain, to_all+scalar_pair, cgrid_ne)
154 call pass_vector(og%Dopen_u, og%Dopen_v, og%Domain, to_all+scalar_pair, cgrid_ne)
157 call set_derived_metrics(og)
159 end subroutine copy_dyngrid_to_mom_grid
164 subroutine copy_mom_grid_to_dyngrid(oG, dG)
168 integer :: isd, ied, jsd, jed
169 integer :: isdb, iedb, jsdb, jedb
170 integer :: ido, jdo, ido2, jdo2
171 integer :: igst, jgst
178 ido = og%idG_offset - dg%idG_offset
179 jdo = og%jdG_offset - dg%jdG_offset
181 isd = max(dg%isd, og%isd+ido) ; jsd = max(dg%jsd, og%jsd+jdo)
182 ied = min(dg%ied, og%ied+ido) ; jed = min(dg%jed, og%jed+jdo)
183 isdb = max(dg%IsdB, og%IsdB+ido) ; jsdb = max(dg%JsdB, og%JsdB+jdo)
184 iedb = min(dg%IedB, og%IedB+ido) ; jedb = min(dg%JedB, og%JedB+jdo)
187 if ((isd > dg%isc) .or. (ied < dg%ied) .or. (jsd > dg%jsc) .or. (jed > dg%jed)) &
188 call mom_error(fatal,
"copy_dyngrid_to_MOM_grid called with incompatible grids.")
190 do i=isd,ied ;
do j=jsd,jed
191 dg%geoLonT(i,j) = og%geoLonT(i+ido,j+jdo)
192 dg%geoLatT(i,j) = og%geoLatT(i+ido,j+jdo)
193 dg%dxT(i,j) = og%dxT(i+ido,j+jdo)
194 dg%dyT(i,j) = og%dyT(i+ido,j+jdo)
195 dg%areaT(i,j) = og%areaT(i+ido,j+jdo)
196 dg%bathyT(i,j) = og%bathyT(i+ido,j+jdo)
198 dg%dF_dx(i,j) = og%dF_dx(i+ido,j+jdo)
199 dg%dF_dy(i,j) = og%dF_dy(i+ido,j+jdo)
200 dg%sin_rot(i,j) = og%sin_rot(i+ido,j+jdo)
201 dg%cos_rot(i,j) = og%cos_rot(i+ido,j+jdo)
202 dg%mask2dT(i,j) = og%mask2dT(i+ido,j+jdo)
205 do i=isdb,iedb ;
do j=jsd,jed
206 dg%geoLonCu(i,j) = og%geoLonCu(i+ido,j+jdo)
207 dg%geoLatCu(i,j) = og%geoLatCu(i+ido,j+jdo)
208 dg%dxCu(i,j) = og%dxCu(i+ido,j+jdo)
209 dg%dyCu(i,j) = og%dyCu(i+ido,j+jdo)
210 dg%dy_Cu(i,j) = og%dy_Cu(i+ido,j+jdo)
212 dg%mask2dCu(i,j) = og%mask2dCu(i+ido,j+jdo)
213 dg%areaCu(i,j) = og%areaCu(i+ido,j+jdo)
214 dg%IareaCu(i,j) = og%IareaCu(i+ido,j+jdo)
217 do i=isd,ied ;
do j=jsdb,jedb
218 dg%geoLonCv(i,j) = og%geoLonCv(i+ido,j+jdo)
219 dg%geoLatCv(i,j) = og%geoLatCv(i+ido,j+jdo)
220 dg%dxCv(i,j) = og%dxCv(i+ido,j+jdo)
221 dg%dyCv(i,j) = og%dyCv(i+ido,j+jdo)
222 dg%dx_Cv(i,j) = og%dx_Cv(i+ido,j+jdo)
224 dg%mask2dCv(i,j) = og%mask2dCv(i+ido,j+jdo)
225 dg%areaCv(i,j) = og%areaCv(i+ido,j+jdo)
226 dg%IareaCv(i,j) = og%IareaCv(i+ido,j+jdo)
229 do i=isdb,iedb ;
do j=jsdb,jedb
230 dg%geoLonBu(i,j) = og%geoLonBu(i+ido,j+jdo)
231 dg%geoLatBu(i,j) = og%geoLatBu(i+ido,j+jdo)
232 dg%dxBu(i,j) = og%dxBu(i+ido,j+jdo)
233 dg%dyBu(i,j) = og%dyBu(i+ido,j+jdo)
234 dg%areaBu(i,j) = og%areaBu(i+ido,j+jdo)
235 dg%CoriolisBu(i,j) = og%CoriolisBu(i+ido,j+jdo)
236 dg%mask2dBu(i,j) = og%mask2dBu(i+ido,j+jdo)
239 dg%bathymetry_at_vel = og%bathymetry_at_vel
240 if (dg%bathymetry_at_vel)
then
241 do i=isdb,iedb ;
do j=jsd,jed
242 dg%Dblock_u(i,j) = og%Dblock_u(i+ido,j+jdo)
243 dg%Dopen_u(i,j) = og%Dopen_u(i+ido,j+jdo)
245 do i=isd,ied ;
do j=jsdb,jedb
246 dg%Dblock_v(i,j) = og%Dblock_v(i+ido,j+jdo)
247 dg%Dopen_v(i,j) = og%Dopen_v(i+ido,j+jdo)
251 dg%gridLonT(dg%isg:dg%ieg) = og%gridLonT(og%isg:og%ieg)
252 dg%gridLatT(dg%jsg:dg%jeg) = og%gridLatT(og%jsg:og%jeg)
258 ido2 = og%IegB-dg%IegB ; igst = max(dg%isg-1, og%IsgB-ido2)
259 jdo2 = og%JegB-dg%JegB ; jgst = max(dg%jsg-1, og%JsgB-jdo2)
260 do i=igst,dg%IegB ; dg%gridLonB(i) = og%gridLonB(i+ido2) ;
enddo
261 do j=jgst,dg%JegB ; dg%gridLatB(j) = og%gridLatB(j+jdo2) ;
enddo
264 dg%x_axis_units = og%x_axis_units ; dg%y_axis_units = og%y_axis_units
265 dg%areaT_global = og%areaT_global ; dg%IareaT_global = og%IareaT_global
266 dg%south_lat = og%south_lat ; dg%west_lon = og%west_lon
267 dg%len_lat = og%len_lat ; dg%len_lon = og%len_lon
268 dg%Rad_Earth = og%Rad_Earth ; dg%max_depth = og%max_depth
273 call pass_var(dg%geoLonT, dg%Domain)
274 call pass_var(dg%geoLatT, dg%Domain)
275 call pass_vector(dg%dxT, dg%dyT, dg%Domain, to_all+scalar_pair, agrid)
276 call pass_vector(dg%dF_dx, dg%dF_dy, dg%Domain, to_all, agrid)
277 call pass_vector(dg%cos_rot, dg%sin_rot, dg%Domain, to_all, agrid)
278 call pass_var(dg%mask2dT, dg%Domain)
280 call pass_vector(dg%areaCu, dg%areaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
281 call pass_vector(dg%dyCu, dg%dxCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
282 call pass_vector(dg%dxCu, dg%dyCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
283 call pass_vector(dg%dy_Cu, dg%dx_Cv, dg%Domain, to_all+scalar_pair, cgrid_ne)
284 call pass_vector(dg%mask2dCu, dg%mask2dCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
285 call pass_vector(dg%IareaCu, dg%IareaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
286 call pass_vector(dg%IareaCu, dg%IareaCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
287 call pass_vector(dg%geoLatCu, dg%geoLatCv, dg%Domain, to_all+scalar_pair, cgrid_ne)
289 call pass_var(dg%areaBu, dg%Domain, position=corner)
290 call pass_var(dg%geoLonBu, dg%Domain, position=corner, inner_halo=dg%isc-isd)
291 call pass_var(dg%geoLatBu, dg%Domain, position=corner)
292 call pass_vector(dg%dxBu, dg%dyBu, dg%Domain, to_all+scalar_pair, bgrid_ne)
293 call pass_var(dg%CoriolisBu, dg%Domain, position=corner)
294 call pass_var(dg%mask2dBu, dg%Domain, position=corner)
296 if (dg%bathymetry_at_vel)
then
297 call pass_vector(dg%Dblock_u, dg%Dblock_v, dg%Domain, to_all+scalar_pair, cgrid_ne)
298 call pass_vector(dg%Dopen_u, dg%Dopen_v, dg%Domain, to_all+scalar_pair, cgrid_ne)
301 call set_derived_dyn_horgrid(dg)
303 end subroutine copy_mom_grid_to_dyngrid