MOM6
mom_dyn_horgrid Module Reference

Detailed Description

Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines that work on this type.

Data Types

type  dyn_horgrid_type
 Describes the horizontal ocean grid with only dynamic memory arrays. More...
 

Functions/Subroutines

subroutine, public create_dyn_horgrid (G, HI, bathymetry_at_vel)
 Allocate memory used by the dyn_horgrid_type and related structures. More...
 
subroutine, public rescale_dyn_horgrid_bathymetry (G, m_in_new_units)
 rescale_dyn_horgrid_bathymetry permits a change in the internal units for the bathymetry on the grid, both rescaling the depths and recording the new internal depth units. More...
 
subroutine, public set_derived_dyn_horgrid (G)
 set_derived_dyn_horgrid calculates metric terms that are derived from other metrics. More...
 
real function adcroft_reciprocal (val)
 Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0. More...
 
subroutine, public destroy_dyn_horgrid (G)
 Release memory used by the dyn_horgrid_type and related structures. More...
 

Function/Subroutine Documentation

◆ adcroft_reciprocal()

real function mom_dyn_horgrid::adcroft_reciprocal ( real, intent(in)  val)
private

Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0.

Parameters
[in]valThe value being inverted.
Returns
The Adcroft reciprocal of val.

Definition at line 362 of file MOM_dyn_horgrid.F90.

362  real, intent(in) :: val !< The value being inverted.
363  real :: I_val !< The Adcroft reciprocal of val.
364 
365  i_val = 0.0 ; if (val /= 0.0) i_val = 1.0/val

◆ create_dyn_horgrid()

subroutine, public mom_dyn_horgrid::create_dyn_horgrid ( type(dyn_horgrid_type), pointer  G,
type(hor_index_type), intent(in)  HI,
logical, intent(in), optional  bathymetry_at_vel 
)

Allocate memory used by the dyn_horgrid_type and related structures.

Parameters
gA pointer to the dynamic horizontal grid type
[in]hiA hor_index_type for array extents
[in]bathymetry_at_velIf true, there are separate values for the basin depths at velocity points. Otherwise the effects of topography are entirely determined from thickness points.

Definition at line 175 of file MOM_dyn_horgrid.F90.

175  type(dyn_horgrid_type), pointer :: G !< A pointer to the dynamic horizontal grid type
176  type(hor_index_type), intent(in) :: HI !< A hor_index_type for array extents
177  logical, optional, intent(in) :: bathymetry_at_vel !< If true, there are
178  !! separate values for the basin depths at velocity
179  !! points. Otherwise the effects of topography are
180  !! entirely determined from thickness points.
181  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, isg, ieg, jsg, jeg
182 
183  ! This subroutine allocates the lateral elements of the dyn_horgrid_type that
184  ! are always used and zeros them out.
185 
186  if (associated(g)) then
187  call mom_error(warning, "create_dyn_horgrid called with an associated horgrid_type.")
188  else
189  allocate(g)
190  endif
191 
192  g%HI = hi
193 
194  g%isc = hi%isc ; g%iec = hi%iec ; g%jsc = hi%jsc ; g%jec = hi%jec
195  g%isd = hi%isd ; g%ied = hi%ied ; g%jsd = hi%jsd ; g%jed = hi%jed
196  g%isg = hi%isg ; g%ieg = hi%ieg ; g%jsg = hi%jsg ; g%jeg = hi%jeg
197 
198  g%IscB = hi%IscB ; g%IecB = hi%IecB ; g%JscB = hi%JscB ; g%JecB = hi%JecB
199  g%IsdB = hi%IsdB ; g%IedB = hi%IedB ; g%JsdB = hi%JsdB ; g%JedB = hi%JedB
200  g%IsgB = hi%IsgB ; g%IegB = hi%IegB ; g%JsgB = hi%JsgB ; g%JegB = hi%JegB
201 
202  g%idg_offset = hi%idg_offset ; g%jdg_offset = hi%jdg_offset
203  g%isd_global = g%isd + hi%idg_offset ; g%jsd_global = g%jsd + hi%jdg_offset
204  g%symmetric = hi%symmetric
205 
206  g%bathymetry_at_vel = .false.
207  if (present(bathymetry_at_vel)) g%bathymetry_at_vel = bathymetry_at_vel
208 
209  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
210  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
211  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
212 
213  allocate(g%dxT(isd:ied,jsd:jed)) ; g%dxT(:,:) = 0.0
214  allocate(g%dxCu(isdb:iedb,jsd:jed)) ; g%dxCu(:,:) = 0.0
215  allocate(g%dxCv(isd:ied,jsdb:jedb)) ; g%dxCv(:,:) = 0.0
216  allocate(g%dxBu(isdb:iedb,jsdb:jedb)) ; g%dxBu(:,:) = 0.0
217  allocate(g%IdxT(isd:ied,jsd:jed)) ; g%IdxT(:,:) = 0.0
218  allocate(g%IdxCu(isdb:iedb,jsd:jed)) ; g%IdxCu(:,:) = 0.0
219  allocate(g%IdxCv(isd:ied,jsdb:jedb)) ; g%IdxCv(:,:) = 0.0
220  allocate(g%IdxBu(isdb:iedb,jsdb:jedb)) ; g%IdxBu(:,:) = 0.0
221 
222  allocate(g%dyT(isd:ied,jsd:jed)) ; g%dyT(:,:) = 0.0
223  allocate(g%dyCu(isdb:iedb,jsd:jed)) ; g%dyCu(:,:) = 0.0
224  allocate(g%dyCv(isd:ied,jsdb:jedb)) ; g%dyCv(:,:) = 0.0
225  allocate(g%dyBu(isdb:iedb,jsdb:jedb)) ; g%dyBu(:,:) = 0.0
226  allocate(g%IdyT(isd:ied,jsd:jed)) ; g%IdyT(:,:) = 0.0
227  allocate(g%IdyCu(isdb:iedb,jsd:jed)) ; g%IdyCu(:,:) = 0.0
228  allocate(g%IdyCv(isd:ied,jsdb:jedb)) ; g%IdyCv(:,:) = 0.0
229  allocate(g%IdyBu(isdb:iedb,jsdb:jedb)) ; g%IdyBu(:,:) = 0.0
230 
231  allocate(g%areaT(isd:ied,jsd:jed)) ; g%areaT(:,:) = 0.0
232  allocate(g%IareaT(isd:ied,jsd:jed)) ; g%IareaT(:,:) = 0.0
233  allocate(g%areaBu(isdb:iedb,jsdb:jedb)) ; g%areaBu(:,:) = 0.0
234  allocate(g%IareaBu(isdb:iedb,jsdb:jedb)) ; g%IareaBu(:,:) = 0.0
235 
236  allocate(g%mask2dT(isd:ied,jsd:jed)) ; g%mask2dT(:,:) = 0.0
237  allocate(g%mask2dCu(isdb:iedb,jsd:jed)) ; g%mask2dCu(:,:) = 0.0
238  allocate(g%mask2dCv(isd:ied,jsdb:jedb)) ; g%mask2dCv(:,:) = 0.0
239  allocate(g%mask2dBu(isdb:iedb,jsdb:jedb)) ; g%mask2dBu(:,:) = 0.0
240  allocate(g%geoLatT(isd:ied,jsd:jed)) ; g%geoLatT(:,:) = 0.0
241  allocate(g%geoLatCu(isdb:iedb,jsd:jed)) ; g%geoLatCu(:,:) = 0.0
242  allocate(g%geoLatCv(isd:ied,jsdb:jedb)) ; g%geoLatCv(:,:) = 0.0
243  allocate(g%geoLatBu(isdb:iedb,jsdb:jedb)) ; g%geoLatBu(:,:) = 0.0
244  allocate(g%geoLonT(isd:ied,jsd:jed)) ; g%geoLonT(:,:) = 0.0
245  allocate(g%geoLonCu(isdb:iedb,jsd:jed)) ; g%geoLonCu(:,:) = 0.0
246  allocate(g%geoLonCv(isd:ied,jsdb:jedb)) ; g%geoLonCv(:,:) = 0.0
247  allocate(g%geoLonBu(isdb:iedb,jsdb:jedb)) ; g%geoLonBu(:,:) = 0.0
248 
249  allocate(g%dx_Cv(isd:ied,jsdb:jedb)) ; g%dx_Cv(:,:) = 0.0
250  allocate(g%dy_Cu(isdb:iedb,jsd:jed)) ; g%dy_Cu(:,:) = 0.0
251 
252  allocate(g%areaCu(isdb:iedb,jsd:jed)) ; g%areaCu(:,:) = 0.0
253  allocate(g%areaCv(isd:ied,jsdb:jedb)) ; g%areaCv(:,:) = 0.0
254  allocate(g%IareaCu(isdb:iedb,jsd:jed)) ; g%IareaCu(:,:) = 0.0
255  allocate(g%IareaCv(isd:ied,jsdb:jedb)) ; g%IareaCv(:,:) = 0.0
256 
257  allocate(g%bathyT(isd:ied, jsd:jed)) ; g%bathyT(:,:) = 0.0
258  allocate(g%CoriolisBu(isdb:iedb, jsdb:jedb)) ; g%CoriolisBu(:,:) = 0.0
259  allocate(g%dF_dx(isd:ied, jsd:jed)) ; g%dF_dx(:,:) = 0.0
260  allocate(g%dF_dy(isd:ied, jsd:jed)) ; g%dF_dy(:,:) = 0.0
261 
262  allocate(g%sin_rot(isd:ied,jsd:jed)) ; g%sin_rot(:,:) = 0.0
263  allocate(g%cos_rot(isd:ied,jsd:jed)) ; g%cos_rot(:,:) = 1.0
264 
265  if (g%bathymetry_at_vel) then
266  allocate(g%Dblock_u(isdb:iedb, jsd:jed)) ; g%Dblock_u(:,:) = 0.0
267  allocate(g%Dopen_u(isdb:iedb, jsd:jed)) ; g%Dopen_u(:,:) = 0.0
268  allocate(g%Dblock_v(isd:ied, jsdb:jedb)) ; g%Dblock_v(:,:) = 0.0
269  allocate(g%Dopen_v(isd:ied, jsdb:jedb)) ; g%Dopen_v(:,:) = 0.0
270  endif
271 
272  ! gridLonB and gridLatB are used as edge values in some cases, so they
273  ! always need to use symmetric memory allcoations.
274  allocate(g%gridLonT(isg:ieg)) ; g%gridLonT(:) = 0.0
275  allocate(g%gridLonB(isg-1:ieg)) ; g%gridLonB(:) = 0.0
276  allocate(g%gridLatT(jsg:jeg)) ; g%gridLatT(:) = 0.0
277  allocate(g%gridLatB(jsg-1:jeg)) ; g%gridLatB(:) = 0.0
278 

◆ destroy_dyn_horgrid()

subroutine, public mom_dyn_horgrid::destroy_dyn_horgrid ( type(dyn_horgrid_type), pointer  G)

Release memory used by the dyn_horgrid_type and related structures.

Parameters
gThe dynamic horizontal grid type

Definition at line 371 of file MOM_dyn_horgrid.F90.

371  type(dyn_horgrid_type), pointer :: G !< The dynamic horizontal grid type
372 
373  if (.not.associated(g)) then
374  call mom_error(fatal, "destroy_dyn_horgrid called with an unassociated horgrid_type.")
375  endif
376 
377  deallocate(g%dxT) ; deallocate(g%dxCu) ; deallocate(g%dxCv) ; deallocate(g%dxBu)
378  deallocate(g%IdxT) ; deallocate(g%IdxCu) ; deallocate(g%IdxCv) ; deallocate(g%IdxBu)
379 
380  deallocate(g%dyT) ; deallocate(g%dyCu) ; deallocate(g%dyCv) ; deallocate(g%dyBu)
381  deallocate(g%IdyT) ; deallocate(g%IdyCu) ; deallocate(g%IdyCv) ; deallocate(g%IdyBu)
382 
383  deallocate(g%areaT) ; deallocate(g%IareaT)
384  deallocate(g%areaBu) ; deallocate(g%IareaBu)
385  deallocate(g%areaCu) ; deallocate(g%IareaCu)
386  deallocate(g%areaCv) ; deallocate(g%IareaCv)
387 
388  deallocate(g%mask2dT) ; deallocate(g%mask2dCu)
389  deallocate(g%mask2dCv) ; deallocate(g%mask2dBu)
390 
391  deallocate(g%geoLatT) ; deallocate(g%geoLatCu)
392  deallocate(g%geoLatCv) ; deallocate(g%geoLatBu)
393  deallocate(g%geoLonT) ; deallocate(g%geoLonCu)
394  deallocate(g%geoLonCv) ; deallocate(g%geoLonBu)
395 
396  deallocate(g%dx_Cv) ; deallocate(g%dy_Cu)
397 
398  deallocate(g%bathyT) ; deallocate(g%CoriolisBu)
399  deallocate(g%dF_dx) ; deallocate(g%dF_dy)
400  deallocate(g%sin_rot) ; deallocate(g%cos_rot)
401 
402  if (allocated(g%Dblock_u)) deallocate(g%Dblock_u)
403  if (allocated(g%Dopen_u)) deallocate(g%Dopen_u)
404  if (allocated(g%Dblock_v)) deallocate(g%Dblock_v)
405  if (allocated(g%Dopen_v)) deallocate(g%Dopen_v)
406 
407  deallocate(g%gridLonT) ; deallocate(g%gridLatT)
408  deallocate(g%gridLonB) ; deallocate(g%gridLatB)
409 
410  deallocate(g%Domain%mpp_domain)
411  deallocate(g%Domain)
412 
413  deallocate(g)
414 

◆ rescale_dyn_horgrid_bathymetry()

subroutine, public mom_dyn_horgrid::rescale_dyn_horgrid_bathymetry ( type(dyn_horgrid_type), intent(inout)  G,
real, intent(in)  m_in_new_units 
)

rescale_dyn_horgrid_bathymetry permits a change in the internal units for the bathymetry on the grid, both rescaling the depths and recording the new internal depth units.

Parameters
[in,out]gThe dynamic horizontal grid type
[in]m_in_new_unitsThe new internal representation of 1 m depth.

Definition at line 284 of file MOM_dyn_horgrid.F90.

284  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
285  real, intent(in) :: m_in_new_units !< The new internal representation of 1 m depth.
286 
287  ! Local variables
288  real :: rescale
289  integer :: i, j, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
290 
291  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
292  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
293 
294  if (m_in_new_units == 1.0) return
295  if (m_in_new_units < 0.0) &
296  call mom_error(fatal, "rescale_grid_bathymetry: Negative depth units are not permitted.")
297  if (m_in_new_units == 0.0) &
298  call mom_error(fatal, "rescale_grid_bathymetry: Zero depth units are not permitted.")
299 
300  rescale = 1.0 / m_in_new_units
301  do j=jsd,jed ; do i=isd,ied
302  g%bathyT(i,j) = rescale*g%bathyT(i,j)
303  enddo ; enddo
304  if (g%bathymetry_at_vel) then ; do j=jsd,jed ; do i=isdb,iedb
305  g%Dblock_u(i,j) = rescale*g%Dblock_u(i,j) ; g%Dopen_u(i,j) = rescale*g%Dopen_u(i,j)
306  enddo ; enddo ; endif
307  if (g%bathymetry_at_vel) then ; do j=jsdb,jedb ; do i=isd,ied
308  g%Dblock_v(i,j) = rescale*g%Dblock_v(i,j) ; g%Dopen_v(i,j) = rescale*g%Dopen_v(i,j)
309  enddo ; enddo ; endif
310  g%max_depth = rescale*g%max_depth
311 

◆ set_derived_dyn_horgrid()

subroutine, public mom_dyn_horgrid::set_derived_dyn_horgrid ( type(dyn_horgrid_type), intent(inout)  G)

set_derived_dyn_horgrid calculates metric terms that are derived from other metrics.

Parameters
[in,out]gThe dynamic horizontal grid type

Definition at line 316 of file MOM_dyn_horgrid.F90.

316  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
317 ! Various inverse grid spacings and derived areas are calculated within this
318 ! subroutine.
319  integer :: i, j, isd, ied, jsd, jed
320  integer :: IsdB, IedB, JsdB, JedB
321 
322  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
323  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
324 
325  do j=jsd,jed ; do i=isd,ied
326  if (g%dxT(i,j) < 0.0) g%dxT(i,j) = 0.0
327  if (g%dyT(i,j) < 0.0) g%dyT(i,j) = 0.0
328  g%IdxT(i,j) = adcroft_reciprocal(g%dxT(i,j))
329  g%IdyT(i,j) = adcroft_reciprocal(g%dyT(i,j))
330  g%IareaT(i,j) = adcroft_reciprocal(g%areaT(i,j))
331  enddo ; enddo
332 
333  do j=jsd,jed ; do i=isdb,iedb
334  if (g%dxCu(i,j) < 0.0) g%dxCu(i,j) = 0.0
335  if (g%dyCu(i,j) < 0.0) g%dyCu(i,j) = 0.0
336  g%IdxCu(i,j) = adcroft_reciprocal(g%dxCu(i,j))
337  g%IdyCu(i,j) = adcroft_reciprocal(g%dyCu(i,j))
338  enddo ; enddo
339 
340  do j=jsdb,jedb ; do i=isd,ied
341  if (g%dxCv(i,j) < 0.0) g%dxCv(i,j) = 0.0
342  if (g%dyCv(i,j) < 0.0) g%dyCv(i,j) = 0.0
343  g%IdxCv(i,j) = adcroft_reciprocal(g%dxCv(i,j))
344  g%IdyCv(i,j) = adcroft_reciprocal(g%dyCv(i,j))
345  enddo ; enddo
346 
347  do j=jsdb,jedb ; do i=isdb,iedb
348  if (g%dxBu(i,j) < 0.0) g%dxBu(i,j) = 0.0
349  if (g%dyBu(i,j) < 0.0) g%dyBu(i,j) = 0.0
350 
351  g%IdxBu(i,j) = adcroft_reciprocal(g%dxBu(i,j))
352  g%IdyBu(i,j) = adcroft_reciprocal(g%dyBu(i,j))
353  ! areaBu has usually been set to a positive area elsewhere.
354  if (g%areaBu(i,j) <= 0.0) g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
355  g%IareaBu(i,j) = adcroft_reciprocal(g%areaBu(i,j))
356  enddo ; enddo
357