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 4139 of file MOM_diag_mediator.F90.

4139  real, dimension(:,:), intent(in) :: field_in !< Original field to be down sampled
4140  real, dimension(:,:), pointer :: field_out !< Down sampled field
4141  integer, intent(in) :: dl !< Level of down sampling
4142  integer, intent(in) :: isc_o !< Original i-start index
4143  integer, intent(in) :: jsc_o !< Original j-start index
4144  integer, intent(in) :: isc_d !< Computational i-start index of down sampled data
4145  integer, intent(in) :: iec_d !< Computational i-end index of down sampled data
4146  integer, intent(in) :: jsc_d !< Computational j-start index of down sampled data
4147  integer, intent(in) :: jec_d !< Computational j-end index of down sampled data
4148  integer, intent(in) :: isd_d !< Computational i-start index of down sampled data
4149  integer, intent(in) :: ied_d !< Computational i-end index of down sampled data
4150  integer, intent(in) :: jsd_d !< Computational j-start index of down sampled data
4151  integer, intent(in) :: jed_d !< Computational j-end index of down sampled data
4152  ! Locals
4153  integer :: i,j,ii,jj,i0,j0
4154  real :: tot_non_zero
4155  ! down sampled mask = 0 unless the mask value of one of the down sampling cells is 1
4156  allocate(field_out(isd_d:ied_d,jsd_d:jed_d))
4157  field_out(:,:) = 0.0
4158  do j=jsc_d,jec_d ; do i=isc_d,iec_d
4159  i0 = isc_o+dl*(i-isc_d)
4160  j0 = jsc_o+dl*(j-jsc_d)
4161  tot_non_zero = 0.0
4162  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4163  tot_non_zero = tot_non_zero + field_in(ii,jj)
4164  enddo;enddo
4165  if(tot_non_zero > 0.0) field_out(i,j)=1.0
4166  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 4174 of file MOM_diag_mediator.F90.

4174  real, dimension(:,:,:), intent(in) :: field_in !< Original field to be down sampled
4175  real, dimension(:,:,:), pointer :: field_out !< down sampled field
4176  integer, intent(in) :: dl !< Level of down sampling
4177  integer, intent(in) :: isc_o !< Original i-start index
4178  integer, intent(in) :: jsc_o !< Original j-start index
4179  integer, intent(in) :: isc_d !< Computational i-start index of down sampled data
4180  integer, intent(in) :: iec_d !< Computational i-end index of down sampled data
4181  integer, intent(in) :: jsc_d !< Computational j-start index of down sampled data
4182  integer, intent(in) :: jec_d !< Computational j-end index of down sampled data
4183  integer, intent(in) :: isd_d !< Computational i-start index of down sampled data
4184  integer, intent(in) :: ied_d !< Computational i-end index of down sampled data
4185  integer, intent(in) :: jsd_d !< Computational j-start index of down sampled data
4186  integer, intent(in) :: jed_d !< Computational j-end index of down sampled data
4187  ! Locals
4188  integer :: i,j,ii,jj,i0,j0,k,ks,ke
4189  real :: tot_non_zero
4190  ! down sampled mask = 0 unless the mask value of one of the down sampling cells is 1
4191  ks = lbound(field_in,3) ; ke = ubound(field_in,3)
4192  allocate(field_out(isd_d:ied_d,jsd_d:jed_d,ks:ke))
4193  field_out(:,:,:) = 0.0
4194  do k= ks,ke ; do j=jsc_d,jec_d ; do i=isc_d,iec_d
4195  i0 = isc_o+dl*(i-isc_d)
4196  j0 = jsc_o+dl*(j-jsc_d)
4197  tot_non_zero = 0.0
4198  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4199  tot_non_zero = tot_non_zero + field_in(ii,jj,k)
4200  enddo;enddo
4201  if(tot_non_zero > 0.0) field_out(i,j,k)=1.0
4202  enddo; enddo; enddo

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