MOM6
mom_diag_mediator::downsample_field Interface Reference

Detailed Description

Down sample a field.

Definition at line 75 of file MOM_diag_mediator.F90.

Private functions

subroutine downsample_field_2d (field_in, field_out, dl, method, mask, diag_cs, diag, isv_o, jsv_o, isv_d, iev_d, jsv_d, jev_d)
 This subroutine allocates and computes a down sampled 2d array given an input array The down sample method is based on the "cell_methods" for the diagnostics as explained in the above table. More...
 
subroutine downsample_field_3d (field_in, field_out, dl, method, mask, diag_cs, diag, isv_o, jsv_o, isv_d, iev_d, jsv_d, jev_d)
 This subroutine allocates and computes a down sampled 3d array given an input array The down sample method is based on the "cell_methods" for the diagnostics as explained in the above table. More...
 

Functions and subroutines

◆ downsample_field_2d()

subroutine mom_diag_mediator::downsample_field::downsample_field_2d ( real, dimension(:,:), pointer  field_in,
real, dimension(:,:), allocatable  field_out,
integer, intent(in)  dl,
integer, intent(in)  method,
real, dimension(:,:), pointer  mask,
type(diag_ctrl), intent(in)  diag_cs,
type(diag_type), intent(in)  diag,
integer, intent(in)  isv_o,
integer, intent(in)  jsv_o,
integer, intent(in)  isv_d,
integer, intent(in)  iev_d,
integer, intent(in)  jsv_d,
integer, intent(in)  jev_d 
)
private

This subroutine allocates and computes a down sampled 2d array given an input array The down sample method is based on the "cell_methods" for the diagnostics as explained in the above table.

Parameters
field_inOriginal field to be down sampled
field_outDown sampled field
[in]dlLevel of down sampling
[in]methodSampling method
maskMask for field
[in]diag_csStructure used to regulate diagnostic output
[in]diagA structure describing the diagnostic to post
[in]isv_oOriginal i-start index
[in]jsv_oOriginal j-start index
[in]isv_di-start index of down sampled data
[in]iev_di-end index of down sampled data
[in]jsv_dj-start index of down sampled data
[in]jev_dj-end index of down sampled data

Definition at line 4005 of file MOM_diag_mediator.F90.

4005  real, dimension(:,:), pointer :: field_in !< Original field to be down sampled
4006  real, dimension(:,:), allocatable :: field_out !< Down sampled field
4007  integer, intent(in) :: dl !< Level of down sampling
4008  integer, intent(in) :: method !< Sampling method
4009  real, dimension(:,:), pointer :: mask !< Mask for field
4010  type(diag_ctrl), intent(in) :: diag_CS !< Structure used to regulate diagnostic output
4011  type(diag_type), intent(in) :: diag !< A structure describing the diagnostic to post
4012  integer, intent(in) :: isv_o !< Original i-start index
4013  integer, intent(in) :: jsv_o !< Original j-start index
4014  integer, intent(in) :: isv_d !< i-start index of down sampled data
4015  integer, intent(in) :: iev_d !< i-end index of down sampled data
4016  integer, intent(in) :: jsv_d !< j-start index of down sampled data
4017  integer, intent(in) :: jev_d !< j-end index of down sampled data
4018  ! Locals
4019  character(len=240) :: mesg
4020  integer :: i,j,ii,jj,i0,j0,f1,f2,f_in1,f_in2
4021  real :: ave, total_weight, weight
4022  real :: epsilon = 1.0e-20 ! A negligibly small count of weights [nondim]
4023  real :: eps_area ! A negligibly small area [L2 ~> m2]
4024  real :: eps_len ! A negligibly small horizontal length [L ~> m]
4025 
4026  eps_len = 1.0e-20 * diag_cs%G%US%m_to_L
4027  eps_area = 1.0e-20 * diag_cs%G%US%m_to_L**2
4028 
4029  ! Allocate the down sampled field on the down sampled data domain
4030 ! allocate(field_out(diag_cs%dsamp(dl)%isd:diag_cs%dsamp(dl)%ied,diag_cs%dsamp(dl)%jsd:diag_cs%dsamp(dl)%jed))
4031 ! allocate(field_out(1:size(field_in,1)/dl,1:size(field_in,2)/dl))
4032  ! Fill the down sampled field on the down sampled diagnostics (almost always compuate) domain
4033  f_in1 = size(field_in,1)
4034  f_in2 = size(field_in,2)
4035  f1 = f_in1/dl
4036  f2 = f_in2/dl
4037  ! Correction for the symmetric case
4038  if (diag_cs%G%symmetric) then
4039  f1 = f1 + mod(f_in1,dl)
4040  f2 = f2 + mod(f_in2,dl)
4041  endif
4042  allocate(field_out(1:f1,1:f2))
4043 
4044  if (method == mmp) then
4045  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4046  i0 = isv_o+dl*(i-isv_d)
4047  j0 = jsv_o+dl*(j-jsv_d)
4048  ave = 0.0
4049  total_weight = 0.0
4050  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4051 ! do ii=i0,i0+dl-1 ; do jj=j0,j0+dl-1
4052  weight = mask(ii,jj)*diag_cs%G%areaT(ii,jj)
4053  total_weight = total_weight + weight
4054  ave = ave+field_in(ii,jj)*weight
4055  enddo; enddo
4056  field_out(i,j) = ave/(total_weight + eps_area) !Avoid zero mask at all aggregating cells where ave=0.0
4057  enddo; enddo
4058  elseif(method == ssp) then ! e.g., T_dfxy_cont_tendency_2d
4059  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4060  i0 = isv_o+dl*(i-isv_d)
4061  j0 = jsv_o+dl*(j-jsv_d)
4062  ave = 0.0
4063  total_weight = 0.0
4064  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4065 ! do ii=i0,i0+dl-1 ; do jj=j0,j0+dl-1
4066  weight = mask(ii,jj)
4067  total_weight = total_weight + weight
4068  ave=ave+field_in(ii,jj)*weight
4069  enddo; enddo
4070  field_out(i,j) = ave/(total_weight+epsilon) !Avoid zero mask at all aggregating cells where ave=0.0
4071  enddo; enddo
4072  elseif(method == psp) then ! e.g., umo_2d
4073  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4074  i0 = isv_o+dl*(i-isv_d)
4075  j0 = jsv_o+dl*(j-jsv_d)
4076  ave = 0.0
4077  total_weight = 0.0
4078  ii=i0
4079  do jj=j0,j0+dl-1
4080  weight = mask(ii,jj)
4081  total_weight = total_weight +weight
4082  ave=ave+field_in(ii,jj)*weight
4083  enddo
4084  field_out(i,j) = ave/(total_weight+epsilon) !Avoid zero mask at all aggregating cells where ave=0.0
4085  enddo; enddo
4086  elseif(method == spp) then ! e.g., vmo_2d
4087  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4088  i0 = isv_o+dl*(i-isv_d)
4089  j0 = jsv_o+dl*(j-jsv_d)
4090  ave = 0.0
4091  total_weight = 0.0
4092  jj=j0
4093  do ii=i0,i0+dl-1
4094  weight = mask(ii,jj)
4095  total_weight = total_weight +weight
4096  ave=ave+field_in(ii,jj)*weight
4097  enddo
4098  field_out(i,j) = ave/(total_weight+epsilon) !Avoid zero mask at all aggregating cells where ave=0.0
4099  enddo; enddo
4100  elseif(method == pmp) then
4101  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4102  i0 = isv_o+dl*(i-isv_d)
4103  j0 = jsv_o+dl*(j-jsv_d)
4104  ave = 0.0
4105  total_weight = 0.0
4106  ii=i0
4107  do jj=j0,j0+dl-1
4108  weight = mask(ii,jj) * diag_cs%G%dyCu(ii,jj)!*diag_cs%h(ii,jj,1) !Niki?
4109  total_weight = total_weight +weight
4110  ave=ave+field_in(ii,jj)*weight
4111  enddo
4112  field_out(i,j) = ave/(total_weight+eps_len) !Avoid zero mask at all aggregating cells where ave=0.0
4113  enddo; enddo
4114  elseif(method == mpp) then
4115  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4116  i0 = isv_o+dl*(i-isv_d)
4117  j0 = jsv_o+dl*(j-jsv_d)
4118  ave = 0.0
4119  total_weight = 0.0
4120  jj=j0
4121  do ii=i0,i0+dl-1
4122  weight = mask(ii,jj)* diag_cs%G%dxCv(ii,jj)!*diag_cs%h(ii,jj,1) !Niki?
4123  total_weight = total_weight +weight
4124  ave=ave+field_in(ii,jj)*weight
4125  enddo
4126  field_out(i,j) = ave/(total_weight+eps_len) !Avoid zero mask at all aggregating cells where ave=0.0
4127  enddo; enddo
4128  elseif(method == msk) then !The input field is a mask, subsample
4129  field_out(:,:) = 0.0
4130  do j=jsv_d,jev_d ; do i=isv_d,iev_d
4131  i0 = isv_o+dl*(i-isv_d)
4132  j0 = jsv_o+dl*(j-jsv_d)
4133  ave = 0.0
4134  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
4135  ave=ave+field_in(ii,jj)
4136  enddo; enddo
4137  if(ave > 0.0) field_out(i,j)=1.0
4138  enddo; enddo
4139  else
4140  write (mesg,*) " unknown sampling method: ",method
4141  call mom_error(fatal, "downsample_field_2d: "//trim(mesg)//" "//trim(diag%debug_str))
4142  endif
4143 

◆ downsample_field_3d()

subroutine mom_diag_mediator::downsample_field::downsample_field_3d ( real, dimension(:,:,:), pointer  field_in,
real, dimension(:,:,:), allocatable  field_out,
integer, intent(in)  dl,
integer, intent(in)  method,
real, dimension(:,:,:), pointer  mask,
type(diag_ctrl), intent(in)  diag_cs,
type(diag_type), intent(in)  diag,
integer, intent(in)  isv_o,
integer, intent(in)  jsv_o,
integer, intent(in)  isv_d,
integer, intent(in)  iev_d,
integer, intent(in)  jsv_d,
integer, intent(in)  jev_d 
)
private

This subroutine allocates and computes a down sampled 3d array given an input array The down sample method is based on the "cell_methods" for the diagnostics as explained in the above table.

The down sample algorithm

The down sample method could be deduced (before send_data call) from the diagx_cell_method, diagy_cell_method and diagv_cell_method

This is the summary of the down sample algoritm for a diagnostic field f:

\[ f(Id,Jd) = \sum_{i,j} f(Id+i,Jd+j) * weight(Id+i,Jd+j) / [ \sum_{i,j} weight(Id+i,Jd+j)] \]

Here, i and j run from 0 to dl-1 (dl being the down sample level). Id,Jd are the down sampled (coarse grid) indices run over the coarsened compute grid, if and jf are the original (fine grid) indices.

 Example   x_cell y_cell v_cell algorithm_id    implemented weight(if,jf)
 ---------------------------------------------------------------------------------------
 theta     mean   mean   mean   MMM =222        G%areaT(if,jf)*h(if,jf)
 u         point  mean   mean   PMM =022        dyCu(if,jf)*h(if,jf)*delta(if,Id)
 v         mean   point  mean   MPM =202        dxCv(if,jf)*h(if,jf)*delta(jf,Jd)
 ?         point  sum    mean   PSM =012        h(if,jf)*delta(if,Id)
 volcello  sum    sum    sum    SSS =111        1
 T_dfxy_co sum    sum    point  SSP =110        1
 umo       point  sum    sum    PSS =011        1*delta(if,Id)
 vmo       sum    point  sum    SPS =101        1*delta(jf,Jd)
 umo_2d    point  sum    point  PSP =010        1*delta(if,Id)
 vmo_2d    sum    point  point  SPP =100        1*delta(jf,Jd)
 ?         point  mean   point  PMP =020        dyCu(if,jf)*delta(if,Id)
 ?         mean   point  point  MPP =200        dxCv(if,jf)*delta(jf,Jd)
 w         mean   mean   point  MMP =220        G%areaT(if,jf)
 h*theta   mean   mean   sum    MMS =221        G%areaT(if,jf)

 delta is the Kronecker delta
Parameters
field_inOriginal field to be down sampled
field_outdown sampled field
[in]dlLevel of down sampling
[in]methodSampling method
maskMask for field
[in]diag_csStructure used to regulate diagnostic output
[in]diagA structure describing the diagnostic to post
[in]isv_oOriginal i-start index
[in]jsv_oOriginal j-start index
[in]isv_di-start index of down sampled data
[in]iev_di-end index of down sampled data
[in]jsv_dj-start index of down sampled data
[in]jev_dj-end index of down sampled data

Definition at line 3849 of file MOM_diag_mediator.F90.

3849  real, dimension(:,:,:), pointer :: field_in !< Original field to be down sampled
3850  real, dimension(:,:,:), allocatable :: field_out !< down sampled field
3851  integer, intent(in) :: dl !< Level of down sampling
3852  integer, intent(in) :: method !< Sampling method
3853  real, dimension(:,:,:), pointer :: mask !< Mask for field
3854  type(diag_ctrl), intent(in) :: diag_CS !< Structure used to regulate diagnostic output
3855  type(diag_type), intent(in) :: diag !< A structure describing the diagnostic to post
3856  integer, intent(in) :: isv_o !< Original i-start index
3857  integer, intent(in) :: jsv_o !< Original j-start index
3858  integer, intent(in) :: isv_d !< i-start index of down sampled data
3859  integer, intent(in) :: iev_d !< i-end index of down sampled data
3860  integer, intent(in) :: jsv_d !< j-start index of down sampled data
3861  integer, intent(in) :: jev_d !< j-end index of down sampled data
3862  ! Locals
3863  character(len=240) :: mesg
3864  integer :: i,j,ii,jj,i0,j0,f1,f2,f_in1,f_in2
3865  integer :: k,ks,ke
3866  real :: ave,total_weight,weight
3867  real :: eps_vol ! A negligibly small volume or mass [H L2 ~> m3 or kg]
3868  real :: eps_area ! A negligibly small area [L2 ~> m2]
3869  real :: eps_face ! A negligibly small face area [H L ~> m2 or kg m-1]
3870 
3871  ks = 1 ; ke = size(field_in,3)
3872  eps_face = 1.0e-20 * diag_cs%G%US%m_to_L * diag_cs%GV%m_to_H
3873  eps_area = 1.0e-20 * diag_cs%G%US%m_to_L**2
3874  eps_vol = 1.0e-20 * diag_cs%G%US%m_to_L**2 * diag_cs%GV%m_to_H
3875 
3876  ! Allocate the down sampled field on the down sampled data domain
3877 ! allocate(field_out(diag_cs%dsamp(dl)%isd:diag_cs%dsamp(dl)%ied,diag_cs%dsamp(dl)%jsd:diag_cs%dsamp(dl)%jed,ks:ke))
3878 ! allocate(field_out(1:size(field_in,1)/dl,1:size(field_in,2)/dl,ks:ke))
3879  f_in1 = size(field_in,1)
3880  f_in2 = size(field_in,2)
3881  f1 = f_in1/dl
3882  f2 = f_in2/dl
3883  !Correction for the symmetric case
3884  if (diag_cs%G%symmetric) then
3885  f1 = f1 + mod(f_in1,dl)
3886  f2 = f2 + mod(f_in2,dl)
3887  endif
3888  allocate(field_out(1:f1,1:f2,ks:ke))
3889 
3890  ! Fill the down sampled field on the down sampled diagnostics (almost always compuate) domain
3891  if (method == mmm) then
3892  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3893  i0 = isv_o+dl*(i-isv_d)
3894  j0 = jsv_o+dl*(j-jsv_d)
3895  ave = 0.0
3896  total_weight = 0.0
3897  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
3898 ! do ii=i0,i0+dl-1 ; do jj=j0,j0+dl-1 !This seems to be faster!!!!
3899  weight = mask(ii,jj,k) * diag_cs%G%areaT(ii,jj) * diag_cs%h(ii,jj,k)
3900  total_weight = total_weight + weight
3901  ave = ave+field_in(ii,jj,k) * weight
3902  enddo; enddo
3903  field_out(i,j,k) = ave/(total_weight + eps_vol) !Avoid zero mask at all aggregating cells where ave=0.0
3904  enddo; enddo; enddo
3905  elseif (method == sss) then !e.g., volcello
3906  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3907  i0 = isv_o+dl*(i-isv_d)
3908  j0 = jsv_o+dl*(j-jsv_d)
3909  ave = 0.0
3910  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
3911  weight = mask(ii,jj,k)
3912  ave = ave+field_in(ii,jj,k)*weight
3913  enddo; enddo
3914  field_out(i,j,k) = ave !Masked Sum (total_weight=1)
3915  enddo; enddo; enddo
3916  elseif(method == mmp .or. method == mms) then !e.g., T_advection_xy
3917  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3918  i0 = isv_o+dl*(i-isv_d)
3919  j0 = jsv_o+dl*(j-jsv_d)
3920  ave = 0.0
3921  total_weight = 0.0
3922  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
3923 ! do ii=i0,i0+dl-1 ; do jj=j0,j0+dl-1
3924  weight = mask(ii,jj,k) * diag_cs%G%areaT(ii,jj)
3925  total_weight = total_weight + weight
3926  ave = ave+field_in(ii,jj,k)*weight
3927  enddo; enddo
3928  field_out(i,j,k) = ave / (total_weight+eps_area) !Avoid zero mask at all aggregating cells where ave=0.0
3929  enddo; enddo; enddo
3930  elseif(method == pmm) then
3931  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3932  i0 = isv_o+dl*(i-isv_d)
3933  j0 = jsv_o+dl*(j-jsv_d)
3934  ave = 0.0
3935  total_weight = 0.0
3936  ii=i0
3937  do jj=j0,j0+dl-1
3938  weight =mask(ii,jj,k) * diag_cs%G%dyCu(ii,jj) * diag_cs%h(ii,jj,k)
3939  total_weight = total_weight +weight
3940  ave=ave+field_in(ii,jj,k)*weight
3941  enddo
3942  field_out(i,j,k) = ave/(total_weight+eps_face) !Avoid zero mask at all aggregating cells where ave=0.0
3943  enddo; enddo; enddo
3944  elseif(method == pss) then !e.g. umo
3945  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3946  i0 = isv_o+dl*(i-isv_d)
3947  j0 = jsv_o+dl*(j-jsv_d)
3948  ave = 0.0
3949  ii=i0
3950  do jj=j0,j0+dl-1
3951  weight =mask(ii,jj,k)
3952  ave=ave+field_in(ii,jj,k)*weight
3953  enddo
3954  field_out(i,j,k) = ave !Masked Sum (total_weight=1)
3955  enddo; enddo; enddo
3956  elseif(method == sps) then !e.g. vmo
3957  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3958  i0 = isv_o+dl*(i-isv_d)
3959  j0 = jsv_o+dl*(j-jsv_d)
3960  ave = 0.0
3961  jj=j0
3962  do ii=i0,i0+dl-1
3963  weight =mask(ii,jj,k)
3964  ave=ave+field_in(ii,jj,k)*weight
3965  enddo
3966  field_out(i,j,k) = ave !Masked Sum (total_weight=1)
3967  enddo; enddo; enddo
3968  elseif(method == mpm) then
3969  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3970  i0 = isv_o+dl*(i-isv_d)
3971  j0 = jsv_o+dl*(j-jsv_d)
3972  ave = 0.0
3973  total_weight = 0.0
3974  jj=j0
3975  do ii=i0,i0+dl-1
3976  weight = mask(ii,jj,k) * diag_cs%G%dxCv(ii,jj) * diag_cs%h(ii,jj,k)
3977  total_weight = total_weight + weight
3978  ave=ave+field_in(ii,jj,k)*weight
3979  enddo
3980  field_out(i,j,k) = ave/(total_weight+eps_face) !Avoid zero mask at all aggregating cells where ave=0.0
3981  enddo; enddo; enddo
3982  elseif(method == msk) then !The input field is a mask, subsample
3983  field_out(:,:,:) = 0.0
3984  do k= ks,ke ; do j=jsv_d,jev_d ; do i=isv_d,iev_d
3985  i0 = isv_o+dl*(i-isv_d)
3986  j0 = jsv_o+dl*(j-jsv_d)
3987  ave = 0.0
3988  do jj=j0,j0+dl-1 ; do ii=i0,i0+dl-1
3989  ave=ave+field_in(ii,jj,k)
3990  enddo; enddo
3991  if(ave > 0.0) field_out(i,j,k)=1.0
3992  enddo; enddo; enddo
3993  else
3994  write (mesg,*) " unknown sampling method: ",method
3995  call mom_error(fatal, "downsample_field_3d: "//trim(mesg)//" "//trim(diag%debug_str))
3996  endif
3997 

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