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, US)
 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 368 of file MOM_dyn_horgrid.F90.

368  real, intent(in) :: val !< The value being inverted.
369  real :: I_val !< The Adcroft reciprocal of val.
370 
371  i_val = 0.0 ; if (val /= 0.0) i_val = 1.0/val

Referenced by set_derived_dyn_horgrid().

Here is the caller graph for this function:

◆ 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 176 of file MOM_dyn_horgrid.F90.

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

References mom_error_handler::mom_error().

Referenced by mom_oda_driver_mod::init_oda(), and mom::initialize_mom().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 377 of file MOM_dyn_horgrid.F90.

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

References mom_error_handler::mom_error().

Referenced by mom::initialize_mom().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 285 of file MOM_dyn_horgrid.F90.

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

References mom_error_handler::mom_error().

Referenced by mom_ice_shelf::initialize_ice_shelf().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_derived_dyn_horgrid()

subroutine, public mom_dyn_horgrid::set_derived_dyn_horgrid ( type(dyn_horgrid_type), intent(inout)  G,
type(unit_scale_type), intent(in), optional  US 
)

set_derived_dyn_horgrid calculates metric terms that are derived from other metrics.

Parameters
[in,out]gThe dynamic horizontal grid type
[in]usA dimensional unit scaling type

Definition at line 317 of file MOM_dyn_horgrid.F90.

317  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
318  type(unit_scale_type), optional, intent(in) :: US !< A dimensional unit scaling type
319 ! Various inverse grid spacings and derived areas are calculated within this
320 ! subroutine.
321  real :: m_to_L ! A unit conversion factor [L m-1 ~> nondim]
322  real :: L_to_m ! A unit conversion factor [L m-1 ~> nondim]
323  integer :: i, j, isd, ied, jsd, jed
324  integer :: IsdB, IedB, JsdB, JedB
325  m_to_l = 1.0 ; if (present(us)) m_to_l = us%m_to_L
326  l_to_m = 1.0 ; if (present(us)) l_to_m = us%L_to_m
327 
328  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
329  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
330 
331  do j=jsd,jed ; do i=isd,ied
332  if (g%dxT(i,j) < 0.0) g%dxT(i,j) = 0.0
333  if (g%dyT(i,j) < 0.0) g%dyT(i,j) = 0.0
334  g%IdxT(i,j) = adcroft_reciprocal(g%dxT(i,j))
335  g%IdyT(i,j) = adcroft_reciprocal(g%dyT(i,j))
336  g%IareaT(i,j) = adcroft_reciprocal(g%areaT(i,j))
337  enddo ; enddo
338 
339  do j=jsd,jed ; do i=isdb,iedb
340  if (g%dxCu(i,j) < 0.0) g%dxCu(i,j) = 0.0
341  if (g%dyCu(i,j) < 0.0) g%dyCu(i,j) = 0.0
342  g%IdxCu(i,j) = adcroft_reciprocal(g%dxCu(i,j))
343  g%IdyCu(i,j) = adcroft_reciprocal(g%dyCu(i,j))
344  enddo ; enddo
345 
346  do j=jsdb,jedb ; do i=isd,ied
347  if (g%dxCv(i,j) < 0.0) g%dxCv(i,j) = 0.0
348  if (g%dyCv(i,j) < 0.0) g%dyCv(i,j) = 0.0
349  g%IdxCv(i,j) = adcroft_reciprocal(g%dxCv(i,j))
350  g%IdyCv(i,j) = adcroft_reciprocal(g%dyCv(i,j))
351  enddo ; enddo
352 
353  do j=jsdb,jedb ; do i=isdb,iedb
354  if (g%dxBu(i,j) < 0.0) g%dxBu(i,j) = 0.0
355  if (g%dyBu(i,j) < 0.0) g%dyBu(i,j) = 0.0
356 
357  g%IdxBu(i,j) = adcroft_reciprocal(g%dxBu(i,j))
358  g%IdyBu(i,j) = adcroft_reciprocal(g%dyBu(i,j))
359  ! areaBu has usually been set to a positive area elsewhere.
360  if (g%areaBu(i,j) <= 0.0) g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
361  g%IareaBu(i,j) = adcroft_reciprocal(g%areaBu(i,j))
362  enddo ; enddo
363 

References adcroft_reciprocal().

Referenced by mom_transcribe_grid::copy_mom_grid_to_dyngrid(), and mom_grid_initialize::set_grid_metrics().

Here is the call graph for this function:
Here is the caller graph for this function: