MOM6
mom_horizontal_regridding Module Reference

Detailed Description

Horizontal interpolation.

Data Types

interface  fill_boundaries
 Fill grid edges. More...
 
interface  horiz_interp_and_extrap_tracer
 Extrapolate and interpolate data. More...
 

Functions/Subroutines

subroutine fill_miss_2d (aout, good, fill, prev, G, smooth, num_pass, relc, crit, keep_bug, debug)
 Use ICE-9 algorithm to populate points (fill=1) with valid data (good=1). If no information is available, Then use a previous guess (prev). Optionally (smooth) blend the filled points to achieve a more desirable result. More...
 
subroutine horiz_interp_and_extrap_tracer_record (filename, varnam, conversion, recnum, G, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, homogenize, m_to_Z)
 Extrapolate and interpolate from a file record. More...
 
subroutine, public mystats (array, missing, is, ie, js, je, k, mesg)
 Write to the terminal some basic statistics about the k-th level of an array. More...
 
subroutine horiz_interp_and_extrap_tracer_fms_id (fms_id, Time, conversion, G, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, homogenize, m_to_Z)
 Extrapolate and interpolate using a FMS time interpolation handle. More...
 
subroutine meshgrid (x, y, x_T, y_T)
 Create a 2d-mesh of grid coordinates from 1-d arrays. More...
 
integer function, dimension(0:size(m, 1)+1, 0:size(m, 2)+1) fill_boundaries_int (m, cyclic_x, tripolar_n)
 Fill grid edges for integer data. More...
 
real function, dimension(0:size(m, 1)+1, 0:size(m, 2)+1) fill_boundaries_real (m, cyclic_x, tripolar_n)
 Fill grid edges for real data. More...
 
subroutine smooth_heights (zi, fill, bad, sor, niter, cyclic_x, tripolar_n)
 Solve del2 (zi) = 0 using successive iterations with a 5 point stencil. Only points fill==1 are modified. Except where bad==1, information propagates isotropically in index space. The resulting solution in each region is an approximation to del2(zi)=0 subject to boundary conditions along the valid points curve bounding this region. More...
 

Function/Subroutine Documentation

◆ fill_boundaries_int()

integer function, dimension(0:size(m,1)+1,0:size(m,2)+1) mom_horizontal_regridding::fill_boundaries_int ( integer, dimension(:,:), intent(in)  m,
logical, intent(in)  cyclic_x,
logical, intent(in)  tripolar_n 
)
private

Fill grid edges for integer data.

Parameters
[in]minput array (ND)
[in]cyclic_xTrue if domain is zonally re-entrant
[in]tripolar_nTrue if domain has an Arctic fold

Definition at line 903 of file MOM_horizontal_regridding.F90.

903  integer, dimension(:,:), intent(in) :: m !< input array (ND)
904  logical, intent(in) :: cyclic_x !< True if domain is zonally re-entrant
905  logical, intent(in) :: tripolar_n !< True if domain has an Arctic fold
906  integer, dimension(0:size(m,1)+1,0:size(m,2)+1) :: mp
907 
908  real, dimension(size(m,1),size(m,2)) :: m_real
909  real, dimension(0:size(m,1)+1,0:size(m,2)+1) :: mp_real
910 
911  m_real = real(m)
912 
913  mp_real = fill_boundaries_real(m_real,cyclic_x,tripolar_n)
914 
915  mp = int(mp_real)
916 

◆ fill_boundaries_real()

real function, dimension(0:size(m,1)+1,0:size(m,2)+1) mom_horizontal_regridding::fill_boundaries_real ( real, dimension(:,:), intent(in)  m,
logical, intent(in)  cyclic_x,
logical, intent(in)  tripolar_n 
)
private

Fill grid edges for real data.

Parameters
[in]minput array (ND)
[in]cyclic_xTrue if domain is zonally re-entrant
[in]tripolar_nTrue if domain has an Arctic fold

Definition at line 921 of file MOM_horizontal_regridding.F90.

921  real, dimension(:,:), intent(in) :: m !< input array (ND)
922  logical, intent(in) :: cyclic_x !< True if domain is zonally re-entrant
923  logical, intent(in) :: tripolar_n !< True if domain has an Arctic fold
924  real, dimension(0:size(m,1)+1,0:size(m,2)+1) :: mp
925 
926  integer :: ni,nj,i,j
927 
928  ni=size(m,1); nj=size(m,2)
929 
930  mp(1:ni,1:nj)=m(:,:)
931 
932  if (cyclic_x) then
933  mp(0,1:nj)=m(ni,1:nj)
934  mp(ni+1,1:nj)=m(1,1:nj)
935  else
936  mp(0,1:nj)=m(1,1:nj)
937  mp(ni+1,1:nj)=m(ni,1:nj)
938  endif
939 
940  mp(1:ni,0)=m(1:ni,1)
941  if (tripolar_n) then
942  do i=1,ni
943  mp(i,nj+1)=m(ni-i+1,nj)
944  enddo
945  else
946  mp(1:ni,nj+1)=m(1:ni,nj)
947  endif
948 

◆ fill_miss_2d()

subroutine mom_horizontal_regridding::fill_miss_2d ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout)  aout,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  good,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  fill,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in), optional  prev,
type(ocean_grid_type), intent(inout)  G,
logical, intent(in), optional  smooth,
integer, intent(in), optional  num_pass,
real, intent(in), optional  relc,
real, intent(in), optional  crit,
logical, intent(in), optional  keep_bug,
logical, intent(in), optional  debug 
)
private

Use ICE-9 algorithm to populate points (fill=1) with valid data (good=1). If no information is available, Then use a previous guess (prev). Optionally (smooth) blend the filled points to achieve a more desirable result.

Parameters
[in,out]gThe ocean's grid structure.
[in,out]aoutThe array with missing values to fill
[in]goodValid data mask for incoming array
[in]fillSame shape array of points which need
[in]prevFirst guess where isolated holes exist.
[in]smoothIf present and true, apply a number of Laplacian iterations to the interpolated data
[in]num_passThe maximum number of iterations
[in]relcA relaxation coefficient for Laplacian (ND)
[in]critA minimal value for deltas between iterations.
[in]keep_bugUse an algorithm with a bug that dates to the "sienna" code release.
[in]debugIf true, write verbose debugging messages.

Definition at line 105 of file MOM_horizontal_regridding.F90.

105  use mom_coms, only : sum_across_pes
106 
107  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
108  real, dimension(SZI_(G),SZJ_(G)), &
109  intent(inout) :: aout !< The array with missing values to fill
110  real, dimension(SZI_(G),SZJ_(G)), &
111  intent(in) :: good !< Valid data mask for incoming array
112  !! (1==good data; 0==missing data).
113  real, dimension(SZI_(G),SZJ_(G)), &
114  intent(in) :: fill !< Same shape array of points which need
115  !! filling (1==fill;0==dont fill)
116  real, dimension(SZI_(G),SZJ_(G)), &
117  optional, intent(in) :: prev !< First guess where isolated holes exist.
118  logical, optional, intent(in) :: smooth !< If present and true, apply a number of
119  !! Laplacian iterations to the interpolated data
120  integer, optional, intent(in) :: num_pass !< The maximum number of iterations
121  real, optional, intent(in) :: relc !< A relaxation coefficient for Laplacian (ND)
122  real, optional, intent(in) :: crit !< A minimal value for deltas between iterations.
123  logical, optional, intent(in) :: keep_bug !< Use an algorithm with a bug that dates
124  !! to the "sienna" code release.
125  logical, optional, intent(in) :: debug !< If true, write verbose debugging messages.
126 
127 
128  real, dimension(SZI_(G),SZJ_(G)) :: b,r
129  real, dimension(SZI_(G),SZJ_(G)) :: fill_pts, good_, good_new
130 
131  character(len=256) :: mesg ! The text of an error message
132  integer :: i,j,k
133  real :: east,west,north,south,sor
134  real :: ge,gw,gn,gs,ngood
135  logical :: do_smooth,siena_bug
136  real :: nfill, nfill_prev
137  integer, parameter :: num_pass_default = 10000
138  real, parameter :: relc_default = 0.25, crit_default = 1.e-3
139 
140  integer :: npass
141  integer :: is, ie, js, je
142  real :: relax_coeff, acrit, ares
143  logical :: debug_it
144 
145  debug_it=.false.
146  if (PRESENT(debug)) debug_it=debug
147 
148  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
149 
150  npass = num_pass_default
151  if (PRESENT(num_pass)) npass = num_pass
152 
153  relax_coeff = relc_default
154  if (PRESENT(relc)) relax_coeff = relc
155 
156  acrit = crit_default
157  if (PRESENT(crit)) acrit = crit
158 
159  siena_bug=.false.
160  if (PRESENT(keep_bug)) siena_bug = keep_bug
161 
162  do_smooth=.false.
163  if (PRESENT(smooth)) do_smooth=smooth
164 
165  fill_pts(:,:) = fill(:,:)
166 
167  nfill = sum(fill(is:ie,js:je))
168  call sum_across_pes(nfill)
169 
170  nfill_prev = nfill
171  good_(:,:) = good(:,:)
172  r(:,:) = 0.0
173 
174  do while (nfill > 0.0)
175 
176  call pass_var(good_,g%Domain)
177  call pass_var(aout,g%Domain)
178 
179  b(:,:)=aout(:,:)
180  good_new(:,:)=good_(:,:)
181 
182  do j=js,je ; do i=is,ie
183 
184  if (good_(i,j) == 1.0 .or. fill(i,j) == 0.) cycle
185 
186  ge=good_(i+1,j) ; gw=good_(i-1,j)
187  gn=good_(i,j+1) ; gs=good_(i,j-1)
188  east=0.0 ; west=0.0 ; north=0.0 ; south=0.0
189  if (ge == 1.0) east = aout(i+1,j)*ge
190  if (gw == 1.0) west = aout(i-1,j)*gw
191  if (gn == 1.0) north = aout(i,j+1)*gn
192  if (gs == 1.0) south = aout(i,j-1)*gs
193 
194  ngood = ge+gw+gn+gs
195  if (ngood > 0.) then
196  b(i,j)=(east+west+north+south)/ngood
197  !### Replace this with
198  ! b(i,j) = ((east+west) + (north+south))/ngood
199  fill_pts(i,j) = 0.0
200  good_new(i,j) = 1.0
201  endif
202  enddo ; enddo
203 
204  aout(is:ie,js:je) = b(is:ie,js:je)
205  good_(is:ie,js:je) = good_new(is:ie,js:je)
206  nfill_prev = nfill
207  nfill = sum(fill_pts(is:ie,js:je))
208  call sum_across_pes(nfill)
209 
210  if (nfill == nfill_prev .and. PRESENT(prev)) then
211  do j=js,je ; do i=is,ie ; if (fill_pts(i,j) == 1.0) then
212  aout(i,j) = prev(i,j)
213  fill_pts(i,j) = 0.0
214  endif ; enddo ; enddo
215  elseif (nfill == nfill_prev) then
216  call mom_error(warning, &
217  'Unable to fill missing points using either data at the same vertical level from a connected basin'//&
218  'or using a point from a previous vertical level. Make sure that the original data has some valid'//&
219  'data in all basins.', .true.)
220  write(mesg,*) 'nfill=',nfill
221  call mom_error(warning, mesg, .true.)
222  endif
223 
224  nfill = sum(fill_pts(is:ie,js:je))
225  call sum_across_pes(nfill)
226 
227  enddo
228 
229  if (do_smooth) then ; do k=1,npass
230  call pass_var(aout,g%Domain)
231  do j=js,je ; do i=is,ie
232  if (fill(i,j) == 1) then
233  east = max(good(i+1,j),fill(i+1,j)) ; west = max(good(i-1,j),fill(i-1,j))
234  north = max(good(i,j+1),fill(i,j+1)) ; south = max(good(i,j-1),fill(i,j-1))
235  r(i,j) = relax_coeff*(south*aout(i,j-1)+north*aout(i,j+1) + &
236  west*aout(i-1,j)+east*aout(i+1,j) - &
237  (south+north+west+east)*aout(i,j))
238  !### Appropriate parentheses should be added here, but they will change answers.
239  ! r(i,j) = relax_coeff*( ((south*aout(i,j-1) + north*aout(i,j+1)) + &
240  ! (west*aout(i-1,j)+east*aout(i+1,j))) - &
241  ! ((south+north)+(west+east))*aout(i,j) )
242  else
243  r(i,j) = 0.
244  endif
245  enddo ; enddo
246  ares = 0.0
247  do j=js,je ; do i=is,ie
248  aout(i,j) = r(i,j) + aout(i,j)
249  ares = max(ares, abs(r(i,j)))
250  enddo ; enddo
251  call max_across_pes(ares)
252  if (ares <= acrit) exit
253  enddo ; endif
254 
255  do j=js,je ; do i=is,ie
256  if (good_(i,j) == 0.0 .and. fill_pts(i,j) == 1.0) then
257  write(mesg,*) 'In fill_miss, fill, good,i,j= ',fill_pts(i,j),good_(i,j),i,j
258  call mom_error(warning, mesg, .true.)
259  call mom_error(fatal,"MOM_initialize: "// &
260  "fill is true and good is false after fill_miss, how did this happen? ")
261  endif
262  enddo ; enddo
263 

◆ horiz_interp_and_extrap_tracer_fms_id()

subroutine mom_horizontal_regridding::horiz_interp_and_extrap_tracer_fms_id ( integer, intent(in)  fms_id,
type(time_type), intent(in)  Time,
real, intent(in)  conversion,
type(ocean_grid_type), intent(inout)  G,
real, dimension(:,:,:), allocatable  tr_z,
real, dimension(:,:,:), allocatable  mask_z,
real, dimension(:), allocatable  z_in,
real, dimension(:), allocatable  z_edges_in,
real, intent(out)  missing_value,
logical, intent(in)  reentrant_x,
logical, intent(in)  tripolar_n,
logical, intent(in), optional  homogenize,
real, intent(in), optional  m_to_Z 
)
private

Extrapolate and interpolate using a FMS time interpolation handle.

Parameters
[in]fms_idA unique id used by the FMS time interpolator
[in]timeA FMS time type
[in]conversionConversion factor for tracer.
[in,out]gGrid object
tr_zpointer to allocatable tracer array on local model grid and native vertical levels.
mask_zpointer to allocatable tracer mask array on local model grid and native vertical levels.
z_inCell grid values for input data.
z_edges_inCell grid edge values for input data. (Intent out)
[out]missing_valueThe missing value in the returned array.
[in]reentrant_xIf true, this grid is reentrant in the x-direction
[in]tripolar_nIf true, this is a northern tripolar grid
[in]homogenizeIf present and true, horizontally homogenize data to produce perfectly "flat" initial conditions
[in]m_to_zA conversion factor from meters to the units of depth. If missing, GbathyT must be in m.

Definition at line 604 of file MOM_horizontal_regridding.F90.

604 
605  integer, intent(in) :: fms_id !< A unique id used by the FMS time interpolator
606  type(time_type), intent(in) :: Time !< A FMS time type
607  real, intent(in) :: conversion !< Conversion factor for tracer.
608  type(ocean_grid_type), intent(inout) :: G !< Grid object
609  real, allocatable, dimension(:,:,:) :: tr_z !< pointer to allocatable tracer array on local
610  !! model grid and native vertical levels.
611  real, allocatable, dimension(:,:,:) :: mask_z !< pointer to allocatable tracer mask array on
612  !! local model grid and native vertical levels.
613  real, allocatable, dimension(:) :: z_in !< Cell grid values for input data.
614  real, allocatable, dimension(:) :: z_edges_in !< Cell grid edge values for input data. (Intent out)
615  real, intent(out) :: missing_value !< The missing value in the returned array.
616  logical, intent(in) :: reentrant_x !< If true, this grid is reentrant in the x-direction
617  logical, intent(in) :: tripolar_n !< If true, this is a northern tripolar grid
618  logical, optional, intent(in) :: homogenize !< If present and true, horizontally homogenize data
619  !! to produce perfectly "flat" initial conditions
620  real, optional, intent(in) :: m_to_Z !< A conversion factor from meters to the units
621  !! of depth. If missing, G%bathyT must be in m.
622 
623  ! Local variables
624  real, dimension(:,:), allocatable :: tr_in,tr_inp !< A 2-d array for holding input data on
625  !! native horizontal grid and extended grid
626  !! with poles.
627  real, dimension(:,:,:), allocatable :: data_in !< A buffer for storing the full 3-d time-interpolated array
628  !! on the original grid
629  real, dimension(:,:), allocatable :: mask_in !< A 2-d mask for extended input grid.
630 
631  real :: PI_180
632  integer :: rcode, ncid, varid, ndims, id, jd, kd, jdp
633  integer :: i,j,k
634  integer, dimension(4) :: start, count, dims, dim_id
635  real, dimension(:,:), allocatable :: x_in, y_in
636  real, dimension(:), allocatable :: lon_in, lat_in
637  real, dimension(:), allocatable :: lat_inp, last_row
638  real :: max_lat, min_lat, pole, max_depth, npole
639  real :: roundoff ! The magnitude of roundoff, usually ~2e-16.
640  logical :: add_np
641  character(len=8) :: laynum
642  type(horiz_interp_type) :: Interp
643  type(axistype), dimension(4) :: axes_data
644  integer :: is, ie, js, je ! compute domain indices
645  integer :: isc,iec,jsc,jec ! global compute domain indices
646  integer :: isg, ieg, jsg, jeg ! global extent
647  integer :: isd, ied, jsd, jed ! data domain indices
648  integer :: id_clock_read
649  integer, dimension(4) :: fld_sz
650  character(len=12) :: dim_name(4)
651  logical :: debug=.false.
652  real :: npoints,varAvg
653  real, dimension(SZI_(G),SZJ_(G)) :: lon_out, lat_out, tr_out, mask_out
654  real, dimension(SZI_(G),SZJ_(G)) :: good, fill
655  real, dimension(SZI_(G),SZJ_(G)) :: tr_outf,tr_prev
656  real, dimension(SZI_(G),SZJ_(G)) :: good2,fill2
657  real, dimension(SZI_(G),SZJ_(G)) :: nlevs
658 
659  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
660  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
661  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
662 
663  id_clock_read = cpu_clock_id('(Initialize tracer from Z) read', grain=clock_loop)
664 
665 
666  pi_180=atan(1.0)/45.
667 
668  ! Open NetCDF file and if present, extract data and spatial coordinate information
669  ! The convention adopted here requires that the data be written in (i,j,k) ordering.
670 
671  call cpu_clock_begin(id_clock_read)
672 
673  fld_sz = get_external_field_size(fms_id)
674 
675  if (allocated(lon_in)) deallocate(lon_in)
676  if (allocated(lat_in)) deallocate(lat_in)
677  if (allocated(z_in)) deallocate(z_in)
678  if (allocated(z_edges_in)) deallocate(z_edges_in)
679  if (allocated(tr_z)) deallocate(tr_z)
680  if (allocated(mask_z)) deallocate(mask_z)
681 
682  axes_data = get_external_field_axes(fms_id)
683 
684  id = fld_sz(1) ; jd = fld_sz(2) ; kd = fld_sz(3)
685  allocate(lon_in(id),lat_in(jd),z_in(kd),z_edges_in(kd+1))
686  allocate(tr_z(isd:ied,jsd:jed,kd), mask_z(isd:ied,jsd:jed,kd))
687 
688  call mpp_get_axis_data(axes_data(1), lon_in)
689  call mpp_get_axis_data(axes_data(2), lat_in)
690  call mpp_get_axis_data(axes_data(3), z_in)
691 
692  if (present(m_to_z)) then ; do k=1,kd ; z_in(k) = m_to_z * z_in(k) ; enddo ; endif
693 
694  call cpu_clock_end(id_clock_read)
695 
696  missing_value = get_external_field_missing(fms_id)
697 
698 
699 ! extrapolate the input data to the north pole using the northerm-most latitude
700 
701  max_lat = maxval(lat_in)
702  add_np=.false.
703  if (max_lat < 90.0) then
704  add_np = .true.
705  jdp = jd+1
706  allocate(lat_inp(jdp))
707  lat_inp(1:jd) = lat_in(:)
708  lat_inp(jd+1) = 90.0
709  deallocate(lat_in)
710  allocate(lat_in(1:jdp))
711  lat_in(:) = lat_inp(:)
712  else
713  jdp=jd
714  endif
715 
716 ! construct level cell boundaries as the mid-point between adjacent centers
717 
718  z_edges_in(1) = 0.0
719  do k=2,kd
720  z_edges_in(k) = 0.5*(z_in(k-1)+z_in(k))
721  enddo
722  z_edges_in(kd+1) = 2.0*z_in(kd) - z_in(kd-1)
723 
724  call horiz_interp_init()
725 
726  lon_in = lon_in*pi_180
727  lat_in = lat_in*pi_180
728  allocate(x_in(id,jdp), y_in(id,jdp))
729  call meshgrid(lon_in, lat_in, x_in, y_in)
730 
731  lon_out(:,:) = g%geoLonT(:,:)*pi_180
732  lat_out(:,:) = g%geoLatT(:,:)*pi_180
733 
734  allocate(data_in(id,jd,kd)) ; data_in(:,:,:)=0.0
735  allocate(tr_in(id,jd)) ; tr_in(:,:)=0.0
736  allocate(tr_inp(id,jdp)) ; tr_inp(:,:)=0.0
737  allocate(mask_in(id,jdp)) ; mask_in(:,:)=0.0
738  allocate(last_row(id)) ; last_row(:)=0.0
739 
740  max_depth = maxval(g%bathyT)
741  call mpp_max(max_depth)
742 
743  if (z_edges_in(kd+1)<max_depth) z_edges_in(kd+1)=max_depth
744 
745  if (is_root_pe()) &
746  call time_interp_external(fms_id, time, data_in, verbose=.true.)
747 
748 ! loop through each data level and interpolate to model grid.
749 ! after interpolating, fill in points which will be needed
750 ! to define the layers
751 
752 ! roundoff = 3.0*EPSILON(missing_value)
753  roundoff = 1.e-4
754 
755  do k=1,kd
756  write(laynum,'(I8)') k ; laynum = adjustl(laynum)
757 
758  if (is_root_pe()) then
759  tr_in(1:id,1:jd) = data_in(1:id,1:jd,k)
760  if (add_np) then
761  last_row(:)=tr_in(:,jd); pole=0.0;npole=0.0
762  do i=1,id
763  if (abs(tr_in(i,jd)-missing_value) > abs(roundoff*missing_value)) then
764  pole = pole+last_row(i)
765  npole = npole+1.0
766  endif
767  enddo
768  if (npole > 0) then
769  pole=pole/npole
770  else
771  pole=missing_value
772  endif
773  tr_inp(:,1:jd) = tr_in(:,:)
774  tr_inp(:,jdp) = pole
775  else
776  tr_inp(:,:) = tr_in(:,:)
777  endif
778 
779  endif
780 
781  call mpp_sync()
782  call mpp_broadcast(tr_inp, id*jdp, root_pe())
783  call mpp_sync_self()
784 
785  mask_in=0.0
786 
787  do j=1,jdp ; do i=1,id
788  if (abs(tr_inp(i,j)-missing_value) > abs(roundoff*missing_value)) then
789  mask_in(i,j)=1.0
790  tr_inp(i,j) = tr_inp(i,j) * conversion
791  else
792  tr_inp(i,j) = missing_value
793  endif
794  enddo ; enddo
795 
796  ! call fms routine horiz_interp to interpolate input level data to model horizontal grid
797  if (k == 1) then
798  call horiz_interp_new(interp, x_in, y_in, lon_out(is:ie,js:je), lat_out(is:ie,js:je), &
799  interp_method='bilinear', src_modulo=reentrant_x)
800  endif
801 
802  if (debug) then
803  call mystats(tr_in, missing_value, 1, id, 1, jd, k, 'Tracer from file')
804  endif
805 
806  tr_out(:,:) = 0.0
807 
808  call horiz_interp(interp, tr_inp, tr_out(is:ie,js:je), missing_value=missing_value, &
809  new_missing_handle=.true.)
810 
811  mask_out(:,:) = 1.0
812  do j=js,je ; do i=is,ie
813  if (abs(tr_out(i,j)-missing_value) < abs(roundoff*missing_value)) mask_out(i,j) = 0.
814  enddo ; enddo
815 
816  fill(:,:) = 0.0 ; good(:,:) = 0.0
817 
818  npoints = 0 ; varavg = 0.
819  do j=js,je ; do i=is,ie
820  if (mask_out(i,j) < 1.0) then
821  tr_out(i,j) = missing_value
822  else
823  good(i,j) = 1.0
824  npoints = npoints + 1
825  varavg = varavg + tr_out(i,j)
826  endif
827  if ((g%mask2dT(i,j) == 1.0) .and. (z_edges_in(k) <= g%bathyT(i,j)) .and. &
828  (mask_out(i,j) < 1.0)) &
829  fill(i,j)=1.0
830  enddo ; enddo
831  call pass_var(fill, g%Domain)
832  call pass_var(good, g%Domain)
833 
834  if (debug) then
835  call mystats(tr_out, missing_value, is, ie, js, je, k, 'variable from horiz_interp()')
836  endif
837 
838  ! Horizontally homogenize data to produce perfectly "flat" initial conditions
839  if (PRESENT(homogenize)) then ; if (homogenize) then
840  call sum_across_pes(npoints)
841  call sum_across_pes(varavg)
842  if (npoints>0) then
843  varavg = varavg/real(npoints)
844  endif
845  tr_out(:,:) = varavg
846  endif ; endif
847 
848  ! tr_out contains input z-space data on the model grid with missing values
849  ! now fill in missing values using "ICE-nine" algorithm.
850 
851  tr_outf(:,:) = tr_out(:,:)
852  if (k==1) tr_prev(:,:) = tr_outf(:,:)
853  good2(:,:) = good(:,:)
854  fill2(:,:) = fill(:,:)
855 
856  call fill_miss_2d(tr_outf, good2, fill2, tr_prev, g, smooth=.true.)
857 
858 ! if (debug) then
859 ! call hchksum(tr_outf, 'field from fill_miss_2d ', G%HI)
860 ! endif
861 
862 ! call myStats(tr_outf, missing_value, is, ie, js, je, k, 'field from fill_miss_2d()')
863 
864  tr_z(:,:,k) = tr_outf(:,:)*g%mask2dT(:,:)
865  mask_z(:,:,k) = good2(:,:) + fill2(:,:)
866  tr_prev(:,:) = tr_z(:,:,k)
867 
868  if (debug) then
869  call hchksum(tr_prev,'field after fill ',g%HI)
870  endif
871 
872  enddo ! kd
873 

◆ horiz_interp_and_extrap_tracer_record()

subroutine mom_horizontal_regridding::horiz_interp_and_extrap_tracer_record ( character(len=*), intent(in)  filename,
character(len=*), intent(in)  varnam,
real, intent(in)  conversion,
integer, intent(in)  recnum,
type(ocean_grid_type), intent(inout)  G,
real, dimension(:,:,:), allocatable  tr_z,
real, dimension(:,:,:), allocatable  mask_z,
real, dimension(:), allocatable  z_in,
real, dimension(:), allocatable  z_edges_in,
real, intent(out)  missing_value,
logical, intent(in)  reentrant_x,
logical, intent(in)  tripolar_n,
logical, intent(in), optional  homogenize,
real, intent(in), optional  m_to_Z 
)

Extrapolate and interpolate from a file record.

Parameters
[in]filenamePath to file containing tracer to be interpolated.
[in]varnamName of tracer in filee.
[in]conversionConversion factor for tracer.
[in]recnumRecord number of tracer to be read.
[in,out]gGrid object
tr_zpointer to allocatable tracer array on local model grid and input-file vertical levels.
mask_zpointer to allocatable tracer mask array on local model grid and input-file vertical levels.
z_inCell grid values for input data.
z_edges_inCell grid edge values for input data.
[out]missing_valueThe missing value in the returned array.
[in]reentrant_xIf true, this grid is reentrant in the x-direction
[in]tripolar_nIf true, this is a northern tripolar grid
[in]homogenizeIf present and true, horizontally homogenize data to produce perfectly "flat" initial conditions
[in]m_to_zA conversion factor from meters to the units of depth. If missing, GbathyT must be in m.

Definition at line 270 of file MOM_horizontal_regridding.F90.

270 
271  character(len=*), intent(in) :: filename !< Path to file containing tracer to be
272  !! interpolated.
273  character(len=*), intent(in) :: varnam !< Name of tracer in filee.
274  real, intent(in) :: conversion !< Conversion factor for tracer.
275  integer, intent(in) :: recnum !< Record number of tracer to be read.
276  type(ocean_grid_type), intent(inout) :: G !< Grid object
277  real, allocatable, dimension(:,:,:) :: tr_z !< pointer to allocatable tracer array on local
278  !! model grid and input-file vertical levels.
279  real, allocatable, dimension(:,:,:) :: mask_z !< pointer to allocatable tracer mask array on
280  !! local model grid and input-file vertical levels.
281  real, allocatable, dimension(:) :: z_in !< Cell grid values for input data.
282  real, allocatable, dimension(:) :: z_edges_in !< Cell grid edge values for input data.
283  real, intent(out) :: missing_value !< The missing value in the returned array.
284  logical, intent(in) :: reentrant_x !< If true, this grid is reentrant in the x-direction
285  logical, intent(in) :: tripolar_n !< If true, this is a northern tripolar grid
286  logical, optional, intent(in) :: homogenize !< If present and true, horizontally homogenize data
287  !! to produce perfectly "flat" initial conditions
288  real, optional, intent(in) :: m_to_Z !< A conversion factor from meters to the units
289  !! of depth. If missing, G%bathyT must be in m.
290 
291  ! Local variables
292  real, dimension(:,:), allocatable :: tr_in, tr_inp ! A 2-d array for holding input data on
293  ! native horizontal grid and extended grid
294  ! with poles.
295  real, dimension(:,:), allocatable :: mask_in ! A 2-d mask for extended input grid.
296 
297  real :: PI_180
298  integer :: rcode, ncid, varid, ndims, id, jd, kd, jdp
299  integer :: i,j,k
300  integer, dimension(4) :: start, count, dims, dim_id
301  real, dimension(:,:), allocatable :: x_in, y_in
302  real, dimension(:), allocatable :: lon_in, lat_in
303  real, dimension(:), allocatable :: lat_inp, last_row
304  real :: max_lat, min_lat, pole, max_depth, npole
305  real :: roundoff ! The magnitude of roundoff, usually ~2e-16.
306  real :: add_offset, scale_factor
307  logical :: add_np
308  character(len=8) :: laynum
309  type(horiz_interp_type) :: Interp
310  integer :: is, ie, js, je ! compute domain indices
311  integer :: isc,iec,jsc,jec ! global compute domain indices
312  integer :: isg, ieg, jsg, jeg ! global extent
313  integer :: isd, ied, jsd, jed ! data domain indices
314  integer :: id_clock_read
315  character(len=12) :: dim_name(4)
316  logical :: debug=.false.
317  real :: npoints,varAvg
318  real, dimension(SZI_(G),SZJ_(G)) :: lon_out, lat_out, tr_out, mask_out
319  real, dimension(SZI_(G),SZJ_(G)) :: good, fill
320  real, dimension(SZI_(G),SZJ_(G)) :: tr_outf,tr_prev
321  real, dimension(SZI_(G),SZJ_(G)) :: good2,fill2
322  real, dimension(SZI_(G),SZJ_(G)) :: nlevs
323 
324  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
325  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
326  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
327 
328  id_clock_read = cpu_clock_id('(Initialize tracer from Z) read', grain=clock_loop)
329 
330 
331  if (allocated(tr_z)) deallocate(tr_z)
332  if (allocated(mask_z)) deallocate(mask_z)
333  if (allocated(z_edges_in)) deallocate(z_edges_in)
334 
335  pi_180=atan(1.0)/45.
336 
337  ! Open NetCDF file and if present, extract data and spatial coordinate information
338  ! The convention adopted here requires that the data be written in (i,j,k) ordering.
339 
340  call cpu_clock_begin(id_clock_read)
341 
342 
343  rcode = nf90_open(filename, nf90_nowrite, ncid)
344  if (rcode /= 0) call mom_error(fatal,"error opening file "//trim(filename)//&
345  " in hinterp_extrap")
346  rcode = nf90_inq_varid(ncid, varnam, varid)
347  if (rcode /= 0) call mom_error(fatal,"error finding variable "//trim(varnam)//&
348  " in file "//trim(filename)//" in hinterp_extrap")
349 
350  rcode = nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dims)
351  if (rcode /= 0) call mom_error(fatal,'error inquiring dimensions hinterp_extrap')
352  if (ndims < 3) call mom_error(fatal,"Variable "//trim(varnam)//" in file "// &
353  trim(filename)//" has too few dimensions.")
354 
355  rcode = nf90_inquire_dimension(ncid, dims(1), dim_name(1), len=id)
356  if (rcode /= 0) call mom_error(fatal,"error reading dimension 1 data for "// &
357  trim(varnam)//" in file "// trim(filename)//" in hinterp_extrap")
358  rcode = nf90_inq_varid(ncid, dim_name(1), dim_id(1))
359  if (rcode /= 0) call mom_error(fatal,"error finding variable "//trim(dim_name(1))//&
360  " in file "//trim(filename)//" in hinterp_extrap")
361  rcode = nf90_inquire_dimension(ncid, dims(2), dim_name(2), len=jd)
362  if (rcode /= 0) call mom_error(fatal,"error reading dimension 2 data for "// &
363  trim(varnam)//" in file "// trim(filename)//" in hinterp_extrap")
364  rcode = nf90_inq_varid(ncid, dim_name(2), dim_id(2))
365  if (rcode /= 0) call mom_error(fatal,"error finding variable "//trim(dim_name(2))//&
366  " in file "//trim(filename)//" in hinterp_extrap")
367  rcode = nf90_inquire_dimension(ncid, dims(3), dim_name(3), len=kd)
368  if (rcode /= 0) call mom_error(fatal,"error reading dimension 3 data for "// &
369  trim(varnam)//" in file "// trim(filename)//" in hinterp_extrap")
370  rcode = nf90_inq_varid(ncid, dim_name(3), dim_id(3))
371  if (rcode /= 0) call mom_error(fatal,"error finding variable "//trim(dim_name(3))//&
372  " in file "//trim(filename)//" in hinterp_extrap")
373 
374 
375  missing_value=0.0
376  rcode = nf90_get_att(ncid, varid, "_FillValue", missing_value)
377  if (rcode /= 0) call mom_error(fatal,"error finding missing value for "//&
378  trim(varnam)//" in file "// trim(filename)//" in hinterp_extrap")
379 
380  rcode = nf90_get_att(ncid, varid, "add_offset", add_offset)
381  if (rcode /= 0) add_offset = 0.0
382 
383  rcode = nf90_get_att(ncid, varid, "scale_factor", scale_factor)
384  if (rcode /= 0) scale_factor = 1.0
385 
386 
387  if (allocated(lon_in)) deallocate(lon_in)
388  if (allocated(lat_in)) deallocate(lat_in)
389  if (allocated(z_in)) deallocate(z_in)
390  if (allocated(z_edges_in)) deallocate(z_edges_in)
391  if (allocated(tr_z)) deallocate(tr_z)
392  if (allocated(mask_z)) deallocate(mask_z)
393 
394 
395  allocate(lon_in(id),lat_in(jd),z_in(kd),z_edges_in(kd+1))
396  allocate(tr_z(isd:ied,jsd:jed,kd), mask_z(isd:ied,jsd:jed,kd))
397 
398  start = 1; count = 1; count(1) = id
399  rcode = nf90_get_var(ncid, dim_id(1), lon_in, start, count)
400  if (rcode /= 0) call mom_error(fatal,"error reading dimension 1 values for var_name "// &
401  trim(varnam)//",dim_name "//trim(dim_name(1))//" in file "// trim(filename)//" in hinterp_extrap")
402  start = 1; count = 1; count(1) = jd
403  rcode = nf90_get_var(ncid, dim_id(2), lat_in, start, count)
404  if (rcode /= 0) call mom_error(fatal,"error reading dimension 2 values for var_name "// &
405  trim(varnam)//",dim_name "//trim(dim_name(2))//" in file "// trim(filename)//" in hinterp_extrap")
406  start = 1; count = 1; count(1) = kd
407  rcode = nf90_get_var(ncid, dim_id(3), z_in, start, count)
408  if (rcode /= 0) call mom_error(fatal,"error reading dimension 3 values for var_name "// &
409  trim(varnam//",dim_name "//trim(dim_name(3)))//" in file "// trim(filename)//" in hinterp_extrap")
410 
411  call cpu_clock_end(id_clock_read)
412 
413  if (present(m_to_z)) then ; do k=1,kd ; z_in(k) = m_to_z * z_in(k) ; enddo ; endif
414 
415 ! extrapolate the input data to the north pole using the northerm-most latitude
416 
417  max_lat = maxval(lat_in)
418  add_np=.false.
419  if (max_lat < 90.0) then
420  add_np=.true.
421  jdp=jd+1
422  allocate(lat_inp(jdp))
423  lat_inp(1:jd)=lat_in(:)
424  lat_inp(jd+1)=90.0
425  deallocate(lat_in)
426  allocate(lat_in(1:jdp))
427  lat_in(:)=lat_inp(:)
428  else
429  jdp=jd
430  endif
431 
432 ! construct level cell boundaries as the mid-point between adjacent centers
433 
434  z_edges_in(1) = 0.0
435  do k=2,kd
436  z_edges_in(k)=0.5*(z_in(k-1)+z_in(k))
437  enddo
438  z_edges_in(kd+1)=2.0*z_in(kd) - z_in(kd-1)
439 
440  call horiz_interp_init()
441 
442  lon_in = lon_in*pi_180
443  lat_in = lat_in*pi_180
444  allocate(x_in(id,jdp),y_in(id,jdp))
445  call meshgrid(lon_in,lat_in, x_in, y_in)
446 
447  lon_out(:,:) = g%geoLonT(:,:)*pi_180
448  lat_out(:,:) = g%geoLatT(:,:)*pi_180
449 
450 
451  allocate(tr_in(id,jd)) ; tr_in(:,:)=0.0
452  allocate(tr_inp(id,jdp)) ; tr_inp(:,:)=0.0
453  allocate(mask_in(id,jdp)) ; mask_in(:,:)=0.0
454  allocate(last_row(id)) ; last_row(:)=0.0
455 
456  max_depth = maxval(g%bathyT)
457  call mpp_max(max_depth)
458 
459  if (z_edges_in(kd+1)<max_depth) z_edges_in(kd+1)=max_depth
460 
461 
462 ! loop through each data level and interpolate to model grid.
463 ! after interpolating, fill in points which will be needed
464 ! to define the layers
465 
466  roundoff = 3.0*epsilon(missing_value)
467 
468 
469  do k=1,kd
470  write(laynum,'(I8)') k ; laynum = adjustl(laynum)
471 
472  if (is_root_pe()) then
473  start = 1 ; start(3) = k ; count(:) = 1 ; count(1) = id ; count(2) = jd
474  rcode = nf90_get_var(ncid,varid, tr_in, start, count)
475  if (rcode /= 0) call mom_error(fatal,"hinterp_and_extract_from_Fie: "//&
476  "error reading level "//trim(laynum)//" of variable "//&
477  trim(varnam)//" in file "// trim(filename))
478 
479  if (add_np) then
480  last_row(:)=tr_in(:,jd); pole=0.0;npole=0.0
481  do i=1,id
482  if (abs(tr_in(i,jd)-missing_value) > abs(roundoff*missing_value)) then
483  pole = pole+last_row(i)
484  npole = npole+1.0
485  endif
486  enddo
487  if (npole > 0) then
488  pole=pole/npole
489  else
490  pole=missing_value
491  endif
492  tr_inp(:,1:jd) = tr_in(:,:)
493  tr_inp(:,jdp) = pole
494  else
495  tr_inp(:,:) = tr_in(:,:)
496  endif
497 
498  endif
499 
500  call mpp_sync()
501  call mpp_broadcast(tr_inp, id*jdp, root_pe())
502  call mpp_sync_self()
503 
504  mask_in=0.0
505 
506  do j=1,jdp
507  do i=1,id
508  if (abs(tr_inp(i,j)-missing_value) > abs(roundoff*missing_value)) then
509  mask_in(i,j) = 1.0
510  tr_inp(i,j) = (tr_inp(i,j)*scale_factor+add_offset) * conversion
511  else
512  tr_inp(i,j) = missing_value
513  endif
514  enddo
515  enddo
516 
517 
518 ! call fms routine horiz_interp to interpolate input level data to model horizontal grid
519 
520 
521  if (k == 1) then
522  call horiz_interp_new(interp,x_in,y_in,lon_out(is:ie,js:je),lat_out(is:ie,js:je), &
523  interp_method='bilinear',src_modulo=reentrant_x)
524  endif
525 
526  if (debug) then
527  call mystats(tr_inp,missing_value, is,ie,js,je,k,'Tracer from file')
528  endif
529 
530  tr_out(:,:) = 0.0
531 
532  call horiz_interp(interp,tr_inp,tr_out(is:ie,js:je), missing_value=missing_value, new_missing_handle=.true.)
533 
534  mask_out=1.0
535  do j=js,je
536  do i=is,ie
537  if (abs(tr_out(i,j)-missing_value) < abs(roundoff*missing_value)) mask_out(i,j)=0.
538  enddo
539  enddo
540 
541  fill = 0.0; good = 0.0
542 
543  npoints = 0 ; varavg = 0.
544  do j=js,je
545  do i=is,ie
546  if (mask_out(i,j) < 1.0) then
547  tr_out(i,j)=missing_value
548  else
549  good(i,j)=1.0
550  npoints = npoints + 1
551  varavg = varavg + tr_out(i,j)
552  endif
553  if (g%mask2dT(i,j) == 1.0 .and. z_edges_in(k) <= g%bathyT(i,j) .and. mask_out(i,j) < 1.0) &
554  fill(i,j)=1.0
555  enddo
556  enddo
557  call pass_var(fill,g%Domain)
558  call pass_var(good,g%Domain)
559 
560  if (debug) then
561  call mystats(tr_out,missing_value, is,ie,js,je,k,'variable from horiz_interp()')
562  endif
563 
564  ! Horizontally homogenize data to produce perfectly "flat" initial conditions
565  if (PRESENT(homogenize)) then
566  if (homogenize) then
567  call sum_across_pes(npoints)
568  call sum_across_pes(varavg)
569  if (npoints>0) then
570  varavg = varavg/real(npoints)
571  endif
572  tr_out(:,:) = varavg
573  endif
574  endif
575 
576 ! tr_out contains input z-space data on the model grid with missing values
577 ! now fill in missing values using "ICE-nine" algorithm.
578 
579  tr_outf(:,:) = tr_out(:,:)
580  if (k==1) tr_prev(:,:) = tr_outf(:,:)
581  good2(:,:) = good(:,:)
582  fill2(:,:) = fill(:,:)
583 
584  call fill_miss_2d(tr_outf, good2, fill2, tr_prev, g, smooth=.true.)
585  call mystats(tr_outf, missing_value, is, ie, js, je, k, 'field from fill_miss_2d()')
586 
587  tr_z(:,:,k) = tr_outf(:,:) * g%mask2dT(:,:)
588  mask_z(:,:,k) = good2(:,:) + fill2(:,:)
589 
590  tr_prev(:,:) = tr_z(:,:,k)
591 
592  if (debug) then
593  call hchksum(tr_prev,'field after fill ',g%HI)
594  endif
595 
596  enddo ! kd
597 

◆ meshgrid()

subroutine mom_horizontal_regridding::meshgrid ( real, dimension(:), intent(in)  x,
real, dimension(:), intent(in)  y,
real, dimension(size(x,1),size(y,1)), intent(inout)  x_T,
real, dimension(size(x,1),size(y,1)), intent(inout)  y_T 
)
private

Create a 2d-mesh of grid coordinates from 1-d arrays.

Parameters
[in]xinput 1-dimensional vector
[in]yinput 1-dimensional vector
[in,out]x_toutput 2-dimensional array
[in,out]y_toutput 2-dimensional array

Definition at line 880 of file MOM_horizontal_regridding.F90.

880  real, dimension(:), intent(in) :: x !< input 1-dimensional vector
881  real, dimension(:), intent(in) :: y !< input 1-dimensional vector
882  real, dimension(size(x,1),size(y,1)), intent(inout) :: x_T !< output 2-dimensional array
883  real, dimension(size(x,1),size(y,1)), intent(inout) :: y_T !< output 2-dimensional array
884 
885  integer :: ni,nj,i,j
886 
887  ni=size(x,1) ; nj=size(y,1)
888 
889  do j=1,nj ; do i=1,ni
890  x_t(i,j) = x(i)
891  enddo ; enddo
892 
893  do j=1,nj ; do i=1,ni
894  y_t(i,j) = y(j)
895  enddo ; enddo
896 

◆ mystats()

subroutine, public mom_horizontal_regridding::mystats ( real, dimension(:,:), intent(in)  array,
real, intent(in)  missing,
integer  is,
integer  ie,
integer  js,
integer  je,
integer  k,
character(len=*)  mesg 
)

Write to the terminal some basic statistics about the k-th level of an array.

Parameters
[in]arrayinput array (ND)
[in]missingmissing value (ND)
isHorizontal loop bounds to calculate statistics for
ieHorizontal loop bounds to calculate statistics for
jsHorizontal loop bounds to calculate statistics for
jeHorizontal loop bounds to calculate statistics for
kLevel to calculate statistics for
mesgLabel to use in message

Definition at line 61 of file MOM_horizontal_regridding.F90.

61  real, dimension(:,:), intent(in) :: array !< input array (ND)
62  real, intent(in) :: missing !< missing value (ND)
63  !!@{
64  !> Horizontal loop bounds to calculate statistics for
65  integer :: is,ie,js,je
66  !!@}
67  integer :: k !< Level to calculate statistics for
68  character(len=*) :: mesg !< Label to use in message
69  ! Local variables
70  real :: minA, maxA
71  integer :: i,j
72  logical :: found
73  character(len=120) :: lMesg
74  mina = 9.e24 ; maxa = -9.e24 ; found = .false.
75 
76  do j=js,je ; do i=is,ie
77  if (array(i,j) /= array(i,j)) stop 'Nan!'
78  if (abs(array(i,j)-missing) > 1.e-6*abs(missing)) then
79  if (found) then
80  mina = min(mina, array(i,j))
81  maxa = max(maxa, array(i,j))
82  else
83  found = .true.
84  mina = array(i,j)
85  maxa = array(i,j)
86  endif
87  endif
88  enddo ; enddo
89  call min_across_pes(mina)
90  call max_across_pes(maxa)
91  if (is_root_pe()) then
92  write(lmesg(1:120),'(2(a,es12.4),a,i3,x,a)') &
93  'init_from_Z: min=',mina,' max=',maxa,' Level=',k,trim(mesg)
94  call mom_mesg(lmesg,2)
95  endif
96 

◆ smooth_heights()

subroutine mom_horizontal_regridding::smooth_heights ( real, dimension(:,:), intent(inout)  zi,
integer, dimension(size(zi,1),size(zi,2)), intent(in)  fill,
integer, dimension(size(zi,1),size(zi,2)), intent(in)  bad,
real, intent(in)  sor,
integer, intent(in)  niter,
logical, intent(in)  cyclic_x,
logical, intent(in)  tripolar_n 
)
private

Solve del2 (zi) = 0 using successive iterations with a 5 point stencil. Only points fill==1 are modified. Except where bad==1, information propagates isotropically in index space. The resulting solution in each region is an approximation to del2(zi)=0 subject to boundary conditions along the valid points curve bounding this region.

Parameters
[in,out]ziinput and output array (ND)
[in]fillsame shape as zi, 1=fill
[in]badsame shape as zi, 1=bad data
[in]sorrelaxation coefficient (ND)
[in]nitermaximum number of iterations
[in]cyclic_xtrue if domain is zonally reentrant
[in]tripolar_ntrue if domain has an Arctic fold

Definition at line 958 of file MOM_horizontal_regridding.F90.

958  real, dimension(:,:), intent(inout) :: zi !< input and output array (ND)
959  integer, dimension(size(zi,1),size(zi,2)), intent(in) :: fill !< same shape as zi, 1=fill
960  integer, dimension(size(zi,1),size(zi,2)), intent(in) :: bad !< same shape as zi, 1=bad data
961  real, intent(in) :: sor !< relaxation coefficient (ND)
962  integer, intent(in) :: niter !< maximum number of iterations
963  logical, intent(in) :: cyclic_x !< true if domain is zonally reentrant
964  logical, intent(in) :: tripolar_n !< true if domain has an Arctic fold
965 
966  ! Local variables
967  real, dimension(size(zi,1),size(zi,2)) :: res, m
968  integer, dimension(size(zi,1),size(zi,2),4) :: B
969  real, dimension(0:size(zi,1)+1,0:size(zi,2)+1) :: mp
970  integer, dimension(0:size(zi,1)+1,0:size(zi,2)+1) :: nm
971  integer :: i,j,k,n
972  integer :: ni,nj
973  real :: Isum, bsum
974 
975  ni=size(zi,1) ; nj=size(zi,2)
976 
977 
978  mp(:,:) = fill_boundaries(zi,cyclic_x,tripolar_n)
979 
980  b(:,:,:) = 0.0
981  nm(:,:) = fill_boundaries(bad,cyclic_x,tripolar_n)
982 
983  do j=1,nj
984  do i=1,ni
985  if (fill(i,j) == 1) then
986  b(i,j,1)=1-nm(i+1,j);b(i,j,2)=1-nm(i-1,j)
987  b(i,j,3)=1-nm(i,j+1);b(i,j,4)=1-nm(i,j-1)
988  endif
989  enddo
990  enddo
991 
992  do n=1,niter
993  do j=1,nj
994  do i=1,ni
995  if (fill(i,j) == 1) then
996  bsum = real(b(i,j,1)+b(i,j,2)+b(i,j,3)+b(i,j,4))
997  isum = 1.0/bsum
998  res(i,j)=isum*(b(i,j,1)*mp(i+1,j)+b(i,j,2)*mp(i-1,j)+&
999  b(i,j,3)*mp(i,j+1)+b(i,j,4)*mp(i,j-1)) - mp(i,j)
1000  endif
1001  enddo
1002  enddo
1003  res(:,:)=res(:,:)*sor
1004 
1005  do j=1,nj
1006  do i=1,ni
1007  mp(i,j)=mp(i,j)+res(i,j)
1008  enddo
1009  enddo
1010 
1011  zi(:,:)=mp(1:ni,1:nj)
1012  mp = fill_boundaries(zi,cyclic_x,tripolar_n)
1013  enddo
1014 
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3