MOM6
mom_diag_mediator::downsample_mask Interface Reference

Detailed Description

Down sample the mask of a field.

Definition at line 80 of file MOM_diag_mediator.F90.

Private functions

subroutine downsample_mask_2d (field_in, field_out, dl, isc_o, jsc_o, isc_d, iec_d, jsc_d, jec_d, isd_d, ied_d, jsd_d, jed_d)
 Allocate and compute the 2d down sampled mask The masks are down sampled based on a minority rule, i.e., a coarse cell is open (1) if at least one of the sub-cells are open, otherwise it's closed (0) More...
 
subroutine downsample_mask_3d (field_in, field_out, dl, isc_o, jsc_o, isc_d, iec_d, jsc_d, jec_d, isd_d, ied_d, jsd_d, jed_d)
 Allocate and compute the 3d down sampled mask The masks are down sampled based on a minority rule, i.e., a coarse cell is open (1) if at least one of the sub-cells are open, otherwise it's closed (0) More...
 

Functions and subroutines

◆ downsample_mask_2d()

subroutine mom_diag_mediator::downsample_mask::downsample_mask_2d ( real, dimension(:,:), intent(in)  field_in,
real, dimension(:,:), pointer  field_out,
integer, intent(in)  dl,
integer, intent(in)  isc_o,
integer, intent(in)  jsc_o,
integer, intent(in)  isc_d,
integer, intent(in)  iec_d,
integer, intent(in)  jsc_d,
integer, intent(in)  jec_d,
integer, intent(in)  isd_d,
integer, intent(in)  ied_d,
integer, intent(in)  jsd_d,
integer, intent(in)  jed_d 
)
private

Allocate and compute the 2d down sampled mask The masks are down sampled based on a minority rule, i.e., a coarse cell is open (1) if at least one of the sub-cells are open, otherwise it's closed (0)

Parameters
[in]field_inOriginal field to be down sampled
field_outDown sampled field
[in]dlLevel of down sampling
[in]isc_oOriginal i-start index
[in]jsc_oOriginal j-start index
[in]isc_dComputational i-start index of down sampled data
[in]iec_dComputational i-end index of down sampled data
[in]jsc_dComputational j-start index of down sampled data
[in]jec_dComputational j-end index of down sampled data
[in]isd_dComputational i-start index of down sampled data
[in]ied_dComputational i-end index of down sampled data
[in]jsd_dComputational j-start index of down sampled data
[in]jed_dComputational j-end index of down sampled data

Definition at line 4151 of file MOM_diag_mediator.F90.

4151  real, dimension(:,:), intent(in) :: field_in !< Original field to be down sampled
4152  real, dimension(:,:), pointer :: field_out !< Down sampled field
4153  integer, intent(in) :: dl !< Level of down sampling
4154  integer, intent(in) :: isc_o !< Original i-start index
4155  integer, intent(in) :: jsc_o !< Original j-start index
4156  integer, intent(in) :: isc_d !< Computational i-start index of down sampled data
4157  integer, intent(in) :: iec_d !< Computational i-end index of down sampled data
4158  integer, intent(in) :: jsc_d !< Computational j-start index of down sampled data
4159  integer, intent(in) :: jec_d !< Computational j-end index of down sampled data
4160  integer, intent(in) :: isd_d !< Computational i-start index of down sampled data
4161  integer, intent(in) :: ied_d !< Computational i-end index of down sampled data
4162  integer, intent(in) :: jsd_d !< Computational j-start index of down sampled data
4163  integer, intent(in) :: jed_d !< Computational j-end index of down sampled data
4164  ! Locals
4165  integer :: i,j,ii,jj,i0,j0
4166  real :: tot_non_zero
4167  ! down sampled mask = 0 unless the mask value of one of the down sampling cells is 1
4168  allocate(field_out(isd_d:ied_d,jsd_d:jed_d))
4169  field_out(:,:) = 0.0
4170  do j=jsc_d,jec_d ; do i=isc_d,iec_d
4171  i0 = isc_o+dl*(i-isc_d)
4172  j0 = jsc_o+dl*(j-jsc_d)
4173  tot_non_zero = 0.0
4174  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4175  tot_non_zero = tot_non_zero + field_in(ii,jj)
4176  enddo;enddo
4177  if(tot_non_zero > 0.0) field_out(i,j)=1.0
4178  enddo; enddo

◆ downsample_mask_3d()

subroutine mom_diag_mediator::downsample_mask::downsample_mask_3d ( real, dimension(:,:,:), intent(in)  field_in,
real, dimension(:,:,:), pointer  field_out,
integer, intent(in)  dl,
integer, intent(in)  isc_o,
integer, intent(in)  jsc_o,
integer, intent(in)  isc_d,
integer, intent(in)  iec_d,
integer, intent(in)  jsc_d,
integer, intent(in)  jec_d,
integer, intent(in)  isd_d,
integer, intent(in)  ied_d,
integer, intent(in)  jsd_d,
integer, intent(in)  jed_d 
)
private

Allocate and compute the 3d down sampled mask The masks are down sampled based on a minority rule, i.e., a coarse cell is open (1) if at least one of the sub-cells are open, otherwise it's closed (0)

Parameters
[in]field_inOriginal field to be down sampled
field_outdown sampled field
[in]dlLevel of down sampling
[in]isc_oOriginal i-start index
[in]jsc_oOriginal j-start index
[in]isc_dComputational i-start index of down sampled data
[in]iec_dComputational i-end index of down sampled data
[in]jsc_dComputational j-start index of down sampled data
[in]jec_dComputational j-end index of down sampled data
[in]isd_dComputational i-start index of down sampled data
[in]ied_dComputational i-end index of down sampled data
[in]jsd_dComputational j-start index of down sampled data
[in]jed_dComputational j-end index of down sampled data

Definition at line 4186 of file MOM_diag_mediator.F90.

4186  real, dimension(:,:,:), intent(in) :: field_in !< Original field to be down sampled
4187  real, dimension(:,:,:), pointer :: field_out !< down sampled field
4188  integer, intent(in) :: dl !< Level of down sampling
4189  integer, intent(in) :: isc_o !< Original i-start index
4190  integer, intent(in) :: jsc_o !< Original j-start index
4191  integer, intent(in) :: isc_d !< Computational i-start index of down sampled data
4192  integer, intent(in) :: iec_d !< Computational i-end index of down sampled data
4193  integer, intent(in) :: jsc_d !< Computational j-start index of down sampled data
4194  integer, intent(in) :: jec_d !< Computational j-end index of down sampled data
4195  integer, intent(in) :: isd_d !< Computational i-start index of down sampled data
4196  integer, intent(in) :: ied_d !< Computational i-end index of down sampled data
4197  integer, intent(in) :: jsd_d !< Computational j-start index of down sampled data
4198  integer, intent(in) :: jed_d !< Computational j-end index of down sampled data
4199  ! Locals
4200  integer :: i,j,ii,jj,i0,j0,k,ks,ke
4201  real :: tot_non_zero
4202  ! down sampled mask = 0 unless the mask value of one of the down sampling cells is 1
4203  ks = lbound(field_in,3) ; ke = ubound(field_in,3)
4204  allocate(field_out(isd_d:ied_d,jsd_d:jed_d,ks:ke))
4205  field_out(:,:,:) = 0.0
4206  do k= ks,ke ; do j=jsc_d,jec_d ; do i=isc_d,iec_d
4207  i0 = isc_o+dl*(i-isc_d)
4208  j0 = jsc_o+dl*(j-jsc_d)
4209  tot_non_zero = 0.0
4210  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4211  tot_non_zero = tot_non_zero + field_in(ii,jj,k)
4212  enddo;enddo
4213  if(tot_non_zero > 0.0) field_out(i,j,k)=1.0
4214  enddo; enddo; enddo

The documentation for this interface was generated from the following file: