MOM6
mom_ale_sponge Module Reference

Detailed Description

This module contains the routines used to apply sponge layers when using the ALE mode.

Applying sponges requires the following:

  1. initialize_ALE_sponge
  2. set_up_ALE_sponge_field (tracers) and set_up_ALE_sponge_vel_field (vel)
  3. apply_ALE_sponge
  4. init_ALE_sponge_diags (not being used for now)
  5. ALE_sponge_end (not being used for now)

Data Types

type  ale_sponge_cs
 ALE sponge control structure. More...
 
interface  initialize_ale_sponge
 Ddetermine the number of points which are within sponges in this computational domain. More...
 
type  p2d
 A structure for creating arrays of pointers to 2D arrays with extra gridding information. More...
 
type  p3d
 A structure for creating arrays of pointers to 3D arrays with extra gridding information. More...
 
interface  set_up_ale_sponge_field
 Store the reference profile at h points for a variable. More...
 
interface  set_up_ale_sponge_vel_field
 This subroutine stores the reference profile at u and v points for a vector. More...
 

Functions/Subroutines

subroutine initialize_ale_sponge_fixed (Iresttime, G, param_file, CS, data_h, nz_data)
 This subroutine determines the number of points which are within sponges in this computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. It also stores the target interface heights. More...
 
integer function, public get_ale_sponge_nz_data (CS)
 Return the number of layers in the data with a fixed ALE sponge, or 0 if there are no sponge columns on this PE. More...
 
subroutine, public get_ale_sponge_thicknesses (G, data_h, sponge_mask, CS)
 Return the thicknesses used for the data with a fixed ALE sponge. More...
 
subroutine initialize_ale_sponge_varying (Iresttime, G, param_file, CS)
 This subroutine determines the number of points which are to be restoref in the computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. More...
 
subroutine, public init_ale_sponge_diags (Time, G, diag, CS)
 Initialize diagnostics for the ALE_sponge module. More...
 
subroutine set_up_ale_sponge_field_fixed (sp_val, G, f_ptr, CS)
 This subroutine stores the reference profile at h points for the variable whose address is given by f_ptr. More...
 
subroutine set_up_ale_sponge_field_varying (filename, fieldname, Time, G, GV, US, f_ptr, CS)
 This subroutine stores the reference profile at h points for the variable whose address is given by filename and fieldname. More...
 
subroutine set_up_ale_sponge_vel_field_fixed (u_val, v_val, G, u_ptr, v_ptr, CS)
 This subroutine stores the reference profile at u and v points for the variable whose address is given by u_ptr and v_ptr. More...
 
subroutine set_up_ale_sponge_vel_field_varying (filename_u, fieldname_u, filename_v, fieldname_v, Time, G, US, CS, u_ptr, v_ptr)
 This subroutine stores the reference profile at uand v points for the variable whose address is given by u_ptr and v_ptr. More...
 
subroutine, public apply_ale_sponge (h, dt, G, GV, US, CS, Time)
 This subroutine applies damping to the layers thicknesses, temp, salt and a variety of tracers for every column where there is damping. More...
 
subroutine, public rotate_ale_sponge (sponge_in, G_in, sponge, G, turns, param_file)
 Rotate the ALE sponge fields from the input to the model index map. More...
 
subroutine, public update_ale_sponge_field (sponge, p_old, G, GV, p_new)
 Scan the ALE sponge variables and replace a prescribed pointer to a new value. More...
 
subroutine, public ale_sponge_end (CS)
 This subroutine deallocates any memory associated with the ALE_sponge module. More...
 

Function/Subroutine Documentation

◆ ale_sponge_end()

subroutine, public mom_ale_sponge::ale_sponge_end ( type(ale_sponge_cs), pointer  CS)

This subroutine deallocates any memory associated with the ALE_sponge module.

Parameters
csA pointer to the control structure that is set by a previous call to initialize_ALE_sponge.

Definition at line 1168 of file MOM_ALE_sponge.F90.

1169  type(ALE_sponge_CS), pointer :: CS !< A pointer to the control structure that is
1170  !! set by a previous call to initialize_ALE_sponge.
1171 
1172  integer :: m
1173 
1174  if (.not.associated(cs)) return
1175 
1176  if (associated(cs%col_i)) deallocate(cs%col_i)
1177  if (associated(cs%col_i_u)) deallocate(cs%col_i_u)
1178  if (associated(cs%col_i_v)) deallocate(cs%col_i_v)
1179  if (associated(cs%col_j)) deallocate(cs%col_j)
1180  if (associated(cs%col_j_u)) deallocate(cs%col_j_u)
1181  if (associated(cs%col_j_v)) deallocate(cs%col_j_v)
1182 
1183  if (associated(cs%Iresttime_col)) deallocate(cs%Iresttime_col)
1184  if (associated(cs%Iresttime_col_u)) deallocate(cs%Iresttime_col_u)
1185  if (associated(cs%Iresttime_col_v)) deallocate(cs%Iresttime_col_v)
1186 
1187  do m=1,cs%fldno
1188  if (associated(cs%Ref_val(m)%p)) deallocate(cs%Ref_val(m)%p)
1189  enddo
1190 
1191  deallocate(cs)
1192 

◆ apply_ale_sponge()

subroutine, public mom_ale_sponge::apply_ale_sponge ( real, dimension(szi_(g),szj_(g),szk_(gv)), intent(inout)  h,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(ale_sponge_cs), pointer  CS,
type(time_type), intent(in), optional  Time 
)

This subroutine applies damping to the layers thicknesses, temp, salt and a variety of tracers for every column where there is damping.

Parameters
[in,out]gThe ocean's grid structure (in).
[in]gvocean vertical grid structure
[in]usA dimensional unit scaling type
[in,out]hLayer thickness [H ~> m or kg m-2] (in)
[in]dtThe amount of time covered by this call [T ~> s].
csA pointer to the control structure for this module that is set by a previous call to initialize_ALE_sponge (in).
[in]timeThe current model date

Definition at line 787 of file MOM_ALE_sponge.F90.

788  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure (in).
789  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
790  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
791  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
792  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2] (in)
793  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s].
794  type(ALE_sponge_CS), pointer :: CS !< A pointer to the control structure for this module
795  !! that is set by a previous call to initialize_ALE_sponge (in).
796  type(time_type), optional, intent(in) :: Time !< The current model date
797 
798  real :: damp ! The timestep times the local damping coefficient [nondim].
799  real :: I1pdamp ! I1pdamp is 1/(1 + damp). [nondim].
800  real :: m_to_Z ! A unit conversion factor from m to Z.
801  real, allocatable, dimension(:) :: tmp_val2 ! data values on the original grid
802  real, dimension(SZK_(G)) :: tmp_val1 ! data values remapped to model grid
803  real :: hu(SZIB_(G), SZJ_(G), SZK_(G)) ! A temporary array for h at u pts
804  real :: hv(SZI_(G), SZJB_(G), SZK_(G)) ! A temporary array for h at v pts
805  real, allocatable, dimension(:,:,:) :: sp_val ! A temporary array for fields
806  real, allocatable, dimension(:,:,:) :: mask_z ! A temporary array for field mask at h pts
807  real, dimension(:), allocatable :: hsrc ! Source thicknesses [Z ~> m].
808  ! Local variables for ALE remapping
809  real, dimension(:), allocatable :: tmpT1d
810  integer :: c, m, nkmb, i, j, k, is, ie, js, je, nz, nz_data
811  integer :: col, total_sponge_cols
812  real, allocatable, dimension(:), target :: z_in, z_edges_in
813  real :: missing_value
814  real :: h_neglect, h_neglect_edge
815  real :: zTopOfCell, zBottomOfCell ! Heights [Z ~> m].
816  integer :: nPoints
817 
818  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
819  if (.not.associated(cs)) return
820 
821  if (.not.cs%remap_answers_2018) then
822  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
823  elseif (gv%Boussinesq) then
824  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
825  else
826  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
827  endif
828 
829  if (cs%time_varying_sponges) then
830  if (.not. present(time)) &
831  call mom_error(fatal,"apply_ALE_sponge: No time information provided")
832  do m=1,cs%fldno
833  nz_data = cs%Ref_val(m)%nz_data
834  allocate(sp_val(g%isd:g%ied,g%jsd:g%jed,1:nz_data))
835  allocate(mask_z(g%isd:g%ied,g%jsd:g%jed,1:nz_data))
836  sp_val(:,:,:) = 0.0
837  mask_z(:,:,:) = 0.0
838  call horiz_interp_and_extrap_tracer(cs%Ref_val(m)%id, time, 1.0, g, sp_val, mask_z, z_in, &
839  z_edges_in, missing_value, .true., .false., .false., &
840  spongeongrid=cs%SpongeDataOngrid, m_to_z=us%m_to_Z, &
841  answers_2018=cs%hor_regrid_answers_2018)
842  allocate( hsrc(nz_data) )
843  allocate( tmpt1d(nz_data) )
844  do c=1,cs%num_col
845  i = cs%col_i(c) ; j = cs%col_j(c)
846  cs%Ref_val(m)%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
847  ! Build the source grid
848  ztopofcell = 0. ; zbottomofcell = 0. ; npoints = 0; hsrc(:) = 0.0; tmpt1d(:) = -99.9
849  do k=1,nz_data
850  if (mask_z(cs%col_i(c),cs%col_j(c),k) == 1.0) then
851  zbottomofcell = -min( z_edges_in(k+1), g%bathyT(cs%col_i(c),cs%col_j(c)) )
852  tmpt1d(k) = sp_val(cs%col_i(c),cs%col_j(c),k)
853  elseif (k>1) then
854  zbottomofcell = -g%bathyT(cs%col_i(c),cs%col_j(c))
855  tmpt1d(k) = tmpt1d(k-1)
856  else ! This next block should only ever be reached over land
857  tmpt1d(k) = -99.9
858  endif
859  hsrc(k) = ztopofcell - zbottomofcell
860  if (hsrc(k)>0.) npoints = npoints + 1
861  ztopofcell = zbottomofcell ! Bottom becomes top for next value of k
862  enddo
863  ! In case data is deeper than model
864  hsrc(nz_data) = hsrc(nz_data) + ( ztopofcell + g%bathyT(cs%col_i(c),cs%col_j(c)) )
865  cs%Ref_val(m)%h(1:nz_data,c) = gv%Z_to_H*hsrc(1:nz_data)
866  cs%Ref_val(m)%p(1:nz_data,c) = tmpt1d(1:nz_data)
867  do k=2,nz_data
868  ! if (mask_z(i,j,k)==0.) &
869  if (cs%Ref_val(m)%h(k,c) <= 0.001*gv%m_to_H) &
870  ! some confusion here about why the masks are not correct returning from horiz_interp
871  ! reverting to using a minimum thickness criteria
872  cs%Ref_val(m)%p(k,c) = cs%Ref_val(m)%p(k-1,c)
873  enddo
874  enddo
875  deallocate(sp_val, mask_z, hsrc, tmpt1d)
876  enddo
877  else
878  nz_data = cs%nz_data
879  endif
880 
881  allocate(tmp_val2(nz_data))
882  do m=1,cs%fldno
883  do c=1,cs%num_col
884 ! c is an index for the next 3 lines but a multiplier for the rest of the loop
885 ! Therefore we use c as per C code and increment the index where necessary.
886  i = cs%col_i(c) ; j = cs%col_j(c)
887  damp = dt * cs%Iresttime_col(c)
888  i1pdamp = 1.0 / (1.0 + damp)
889  tmp_val2(1:nz_data) = cs%Ref_val(m)%p(1:nz_data,c)
890  if (cs%time_varying_sponges) then
891  call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_val(m)%h(1:nz_data,c), tmp_val2, &
892  cs%nz, h(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
893  else
894  call remapping_core_h(cs%remap_cs,nz_data, cs%Ref_h%p(1:nz_data,c), tmp_val2, &
895  cs%nz, h(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
896  endif
897  !Backward Euler method
898  cs%var(m)%p(i,j,1:cs%nz) = i1pdamp * (cs%var(m)%p(i,j,1:cs%nz) + tmp_val1 * damp)
899  enddo
900  enddo
901 
902  ! for debugging
903  !c=CS%num_col
904  !do m=1,CS%fldno
905  ! write(*,*) 'APPLY SPONGE,m,CS%Ref_h(:,c),h(i,j,:),tmp_val2,tmp_val1',&
906  ! m,CS%Ref_h(:,c),h(i,j,:),tmp_val2,tmp_val1
907  !enddo
908 
909  if (cs%sponge_uv) then
910  ! u points
911  do j=cs%jsc,cs%jec; do i=cs%iscB,cs%iecB; do k=1,nz
912  hu(i,j,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
913  enddo ; enddo ; enddo
914  if (cs%time_varying_sponges) then
915  if (.not. present(time)) &
916  call mom_error(fatal,"apply_ALE_sponge: No time information provided")
917 
918  nz_data = cs%Ref_val_u%nz_data
919  allocate(sp_val(g%isdB:g%iedB,g%jsd:g%jed,1:nz_data))
920  allocate(mask_z(g%isdB:g%iedB,g%jsd:g%jed,1:nz_data))
921  ! Interpolate from the external horizontal grid and in time
922  call horiz_interp_and_extrap_tracer(cs%Ref_val_u%id, time, 1.0, g, sp_val, mask_z, z_in, &
923  z_edges_in, missing_value, .true., .false., .false., &
924  m_to_z=us%m_to_Z, answers_2018=cs%hor_regrid_answers_2018)
925 ! call pass_var(sp_val,G%Domain)
926 ! call pass_var(mask_z,G%Domain)
927  do c=1,cs%num_col
928  ! c is an index for the next 3 lines but a multiplier for the rest of the loop
929  ! Therefore we use c as per C code and increment the index where necessary.
930  i = cs%col_i(c) ; j = cs%col_j(c)
931  cs%Ref_val_u%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
932  enddo
933 
934  deallocate (sp_val, mask_z)
935 
936  nz_data = cs%Ref_val_v%nz_data
937  allocate(sp_val(g%isd:g%ied,g%jsdB:g%jedB,1:nz_data))
938  allocate(mask_z(g%isd:g%ied,g%jsdB:g%jedB,1:nz_data))
939  ! Interpolate from the external horizontal grid and in time
940  call horiz_interp_and_extrap_tracer(cs%Ref_val_v%id, time, 1.0, g, sp_val, mask_z, z_in, &
941  z_edges_in, missing_value, .true., .false., .false., &
942  m_to_z=us%m_to_Z, answers_2018=cs%hor_regrid_answers_2018)
943 
944 ! call pass_var(sp_val,G%Domain)
945 ! call pass_var(mask_z,G%Domain)
946 
947  do c=1,cs%num_col
948  ! c is an index for the next 3 lines but a multiplier for the rest of the loop
949  ! Therefore we use c as per C code and increment the index where necessary.
950  i = cs%col_i(c) ; j = cs%col_j(c)
951  cs%Ref_val_v%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
952  enddo
953 
954  deallocate (sp_val, mask_z)
955 
956  else
957  nz_data = cs%nz_data
958  endif
959 
960  do c=1,cs%num_col_u
961  i = cs%col_i_u(c) ; j = cs%col_j_u(c)
962  damp = dt * cs%Iresttime_col_u(c)
963  i1pdamp = 1.0 / (1.0 + damp)
964  if (cs%time_varying_sponges) nz_data = cs%Ref_val(m)%nz_data
965  tmp_val2(1:nz_data) = cs%Ref_val_u%p(1:nz_data,c)
966  if (cs%time_varying_sponges) then
967  call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_val_u%h(:,c), tmp_val2, &
968  cs%nz, hu(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
969  else
970  call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_hu%p(:,c), tmp_val2, &
971  cs%nz, hu(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
972  endif
973  !Backward Euler method
974  cs%var_u%p(i,j,:) = i1pdamp * (cs%var_u%p(i,j,:) + tmp_val1 * damp)
975  enddo
976 
977  ! v points
978  do j=cs%jscB,cs%jecB; do i=cs%isc,cs%iec; do k=1,nz
979  hv(i,j,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
980  enddo ; enddo ; enddo
981 
982  do c=1,cs%num_col_v
983  i = cs%col_i_v(c) ; j = cs%col_j_v(c)
984  damp = dt * cs%Iresttime_col_v(c)
985  i1pdamp = 1.0 / (1.0 + damp)
986  tmp_val2(1:nz_data) = cs%Ref_val_v%p(1:nz_data,c)
987  if (cs%time_varying_sponges) then
988  call remapping_core_h(cs%remap_cs, cs%nz_data, cs%Ref_val_v%h(:,c), tmp_val2, &
989  cs%nz, hv(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
990  else
991  call remapping_core_h(cs%remap_cs, cs%nz_data, cs%Ref_hv%p(:,c), tmp_val2, &
992  cs%nz, hv(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
993  endif
994  !Backward Euler method
995  cs%var_v%p(i,j,:) = i1pdamp * (cs%var_v%p(i,j,:) + tmp_val1 * damp)
996  enddo
997 
998  endif
999 
1000  deallocate(tmp_val2)
1001 

◆ get_ale_sponge_nz_data()

integer function, public mom_ale_sponge::get_ale_sponge_nz_data ( type(ale_sponge_cs), pointer  CS)

Return the number of layers in the data with a fixed ALE sponge, or 0 if there are no sponge columns on this PE.

Parameters
csA pointer that is set to point to the control structure for the ALE_sponge module.
Returns
The number of layers in the fixed sponge data.

Definition at line 346 of file MOM_ALE_sponge.F90.

347  type(ALE_sponge_CS), pointer :: CS !< A pointer that is set to point to the control
348  !! structure for the ALE_sponge module.
349  integer :: get_ALE_sponge_nz_data !< The number of layers in the fixed sponge data.
350 
351  if (associated(cs)) then
352  get_ale_sponge_nz_data = cs%nz_data
353  else
354  get_ale_sponge_nz_data = 0
355  endif

◆ get_ale_sponge_thicknesses()

subroutine, public mom_ale_sponge::get_ale_sponge_thicknesses ( type(ocean_grid_type), intent(in)  G,
real, dimension(:,:,:), intent(inout), allocatable  data_h,
logical, dimension(szi_(g),szj_(g)), intent(out)  sponge_mask,
type(ale_sponge_cs), pointer  CS 
)

Return the thicknesses used for the data with a fixed ALE sponge.

Parameters
[in]gThe ocean's grid structure (in).
[in,out]data_hThe thicknesses of the sponge input layers [H ~> m or kg m-2].
[out]sponge_maskA logical mask that is true where
csA pointer that is set to point to the control structure for the ALE_sponge module.

Definition at line 359 of file MOM_ALE_sponge.F90.

360  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure (in).
361  real, allocatable, dimension(:,:,:), &
362  intent(inout) :: data_h !< The thicknesses of the sponge input layers [H ~> m or kg m-2].
363  logical, dimension(SZI_(G),SZJ_(G)), &
364  intent(out) :: sponge_mask !< A logical mask that is true where
365  !! sponges are being applied.
366  type(ALE_sponge_CS), pointer :: CS !< A pointer that is set to point to the control
367  !! structure for the ALE_sponge module.
368  integer :: c, i, j, k
369 
370  if (allocated(data_h)) call mom_error(fatal, &
371  "get_ALE_sponge_thicknesses called with an allocated data_h.")
372 
373  if (.not.associated(cs)) then
374  ! There are no sponge points on this PE.
375  allocate(data_h(g%isd:g%ied,g%jsd:g%jed,1)) ; data_h(:,:,:) = -1.0
376  sponge_mask(:,:) = .false.
377  return
378  endif
379 
380  allocate(data_h(g%isd:g%ied,g%jsd:g%jed,cs%nz_data)) ; data_h(:,:,:) = -1.0
381  sponge_mask(:,:) = .false.
382 
383  do c=1,cs%num_col
384  i = cs%col_i(c) ; j = cs%col_j(c)
385  sponge_mask(i,j) = .true.
386  do k=1,cs%nz_data
387  data_h(i,j,k) = cs%Ref_h%p(k,c)
388  enddo
389  enddo
390 

◆ init_ale_sponge_diags()

subroutine, public mom_ale_sponge::init_ale_sponge_diags ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(diag_ctrl), intent(inout), target  diag,
type(ale_sponge_cs), pointer  CS 
)

Initialize diagnostics for the ALE_sponge module.

Parameters
[in]timeThe current model time
[in]gThe ocean's grid structure
[in,out]diagA structure that is used to regulate diagnostic output.
csALE sponge control structure

Definition at line 548 of file MOM_ALE_sponge.F90.

549  type(time_type), target, intent(in) :: Time !< The current model time
550  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
551  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic
552  !! output.
553  type(ALE_sponge_CS), pointer :: CS !< ALE sponge control structure
554 
555  if (.not.associated(cs)) return
556 
557  cs%diag => diag
558 

◆ initialize_ale_sponge_fixed()

subroutine mom_ale_sponge::initialize_ale_sponge_fixed ( real, dimension(szi_(g),szj_(g)), intent(in)  Iresttime,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(ale_sponge_cs), pointer  CS,
real, dimension(szi_(g),szj_(g),nz_data), intent(in)  data_h,
integer, intent(in)  nz_data 
)

This subroutine determines the number of points which are within sponges in this computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. It also stores the target interface heights.

Parameters
[in]gThe ocean's grid structure.
[in]nz_dataThe total number of sponge input layers.
[in]iresttimeThe inverse of the restoring time [T-1 ~> s-1].
[in]param_fileA structure indicating the open file to parse for model parameter values.
csA pointer that is set to point to the control structure for this module (in/out).
[in]data_hThe thicknesses of the sponge input layers [H ~> m or kg m-2].

Definition at line 148 of file MOM_ALE_sponge.F90.

149 
150  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
151  integer, intent(in) :: nz_data !< The total number of sponge input layers.
152  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: Iresttime !< The inverse of the restoring time [T-1 ~> s-1].
153  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
154  !! to parse for model parameter values.
155  type(ALE_sponge_CS), pointer :: CS !< A pointer that is set to point to the control
156  !! structure for this module (in/out).
157  real, dimension(SZI_(G),SZJ_(G),nz_data), intent(in) :: data_h !< The thicknesses of the sponge
158  !! input layers [H ~> m or kg m-2].
159 
160 
161 ! This include declares and sets the variable "version".
162 #include "version_variable.h"
163  character(len=40) :: mdl = "MOM_sponge" ! This module's name.
164  logical :: use_sponge
165  real, allocatable, dimension(:,:,:) :: data_hu !< thickness at u points [H ~> m or kg m-2]
166  real, allocatable, dimension(:,:,:) :: data_hv !< thickness at v points [H ~> m or kg m-2]
167  real, allocatable, dimension(:,:) :: Iresttime_u !< inverse of the restoring time at u points [T-1 ~> s-1]
168  real, allocatable, dimension(:,:) :: Iresttime_v !< inverse of the restoring time at v points [T-1 ~> s-1]
169  logical :: bndExtrapolation = .true. ! If true, extrapolate boundaries
170  logical :: default_2018_answers
171  integer :: i, j, k, col, total_sponge_cols, total_sponge_cols_u, total_sponge_cols_v
172  character(len=10) :: remapScheme
173  if (associated(cs)) then
174  call mom_error(warning, "initialize_ALE_sponge_fixed called with an associated "// &
175  "control structure.")
176  return
177  endif
178 
179 ! Set default, read and log parameters
180  call log_version(param_file, mdl, version, "")
181  call get_param(param_file, mdl, "SPONGE", use_sponge, &
182  "If true, sponges may be applied anywhere in the domain. "//&
183  "The exact location and properties of those sponges are "//&
184  "specified from MOM_initialization.F90.", default=.false.)
185 
186  if (.not.use_sponge) return
187 
188  allocate(cs)
189 
190  call get_param(param_file, mdl, "SPONGE_UV", cs%sponge_uv, &
191  "Apply sponges in u and v, in addition to tracers.", &
192  default=.false.)
193 
194  call get_param(param_file, mdl, "REMAPPING_SCHEME", remapscheme, &
195  "This sets the reconstruction scheme used "//&
196  " for vertical remapping for all variables.", &
197  default="PLM", do_not_log=.true.)
198 
199  call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", bndextrapolation, &
200  "When defined, a proper high-order reconstruction "//&
201  "scheme is used within boundary cells rather "//&
202  "than PCM. E.g., if PPM is used for remapping, a "//&
203  "PPM reconstruction will also be used within boundary cells.", &
204  default=.false., do_not_log=.true.)
205  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
206  "This sets the default value for the various _2018_ANSWERS parameters.", &
207  default=.false.)
208  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
209  "If true, use the order of arithmetic and expressions that recover the "//&
210  "answers from the end of 2018. Otherwise, use updated and more robust "//&
211  "forms of the same expressions.", default=default_2018_answers)
212  call get_param(param_file, mdl, "HOR_REGRID_2018_ANSWERS", cs%hor_regrid_answers_2018, &
213  "If true, use the order of arithmetic for horizonal regridding that recovers "//&
214  "the answers from the end of 2018. Otherwise, use rotationally symmetric "//&
215  "forms of the same expressions.", default=default_2018_answers)
216 
217  cs%time_varying_sponges = .false.
218  cs%nz = g%ke
219  cs%isc = g%isc ; cs%iec = g%iec ; cs%jsc = g%jsc ; cs%jec = g%jec
220  cs%isd = g%isd ; cs%ied = g%ied ; cs%jsd = g%jsd ; cs%jed = g%jed
221  cs%iscB = g%iscB ; cs%iecB = g%iecB; cs%jscB = g%jscB ; cs%jecB = g%jecB
222 
223  ! number of columns to be restored
224  cs%num_col = 0 ; cs%fldno = 0
225  do j=g%jsc,g%jec ; do i=g%isc,g%iec
226  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
227  cs%num_col = cs%num_col + 1
228  enddo ; enddo
229 
230  if (cs%num_col > 0) then
231  allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
232  allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
233  allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
234  ! pass indices, restoring time to the CS structure
235  col = 1
236  do j=g%jsc,g%jec ; do i=g%isc,g%iec
237  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) then
238  cs%col_i(col) = i ; cs%col_j(col) = j
239  cs%Iresttime_col(col) = iresttime(i,j)
240  col = col +1
241  endif
242  enddo ; enddo
243  ! same for total number of arbritary layers and correspondent data
244  cs%nz_data = nz_data
245  allocate(cs%Ref_h%p(cs%nz_data,cs%num_col))
246  do col=1,cs%num_col ; do k=1,cs%nz_data
247  cs%Ref_h%p(k,col) = data_h(cs%col_i(col),cs%col_j(col),k)
248  enddo ; enddo
249  endif
250 
251  total_sponge_cols = cs%num_col
252  call sum_across_pes(total_sponge_cols)
253 
254 ! Call the constructor for remapping control structure
255  call initialize_remapping(cs%remap_cs, remapscheme, boundary_extrapolation=bndextrapolation, &
256  answers_2018=cs%remap_answers_2018)
257 
258  call log_param(param_file, mdl, "!Total sponge columns at h points", total_sponge_cols, &
259  "The total number of columns where sponges are applied at h points.", like_default=.true.)
260 
261  if (cs%sponge_uv) then
262 
263  allocate(data_hu(g%isdB:g%iedB,g%jsd:g%jed,nz_data)) ; data_hu(:,:,:) = 0.0
264  allocate(data_hv(g%isd:g%ied,g%jsdB:g%jedB,nz_data)) ; data_hv(:,:,:) = 0.0
265  allocate(iresttime_u(g%isdB:g%iedB,g%jsd:g%jed)) ; iresttime_u(:,:) = 0.0
266  allocate(iresttime_v(g%isd:g%ied,g%jsdB:g%jedB)) ; iresttime_v(:,:) = 0.0
267 
268  ! u points
269  cs%num_col_u = 0 ; !CS%fldno_u = 0
270  do j=cs%jsc,cs%jec ; do i=cs%iscB,cs%iecB
271  data_hu(i,j,:) = 0.5 * (data_h(i,j,:) + data_h(i+1,j,:))
272  iresttime_u(i,j) = 0.5 * (iresttime(i,j) + iresttime(i+1,j))
273  if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) &
274  cs%num_col_u = cs%num_col_u + 1
275  enddo ; enddo
276 
277  if (cs%num_col_u > 0) then
278 
279  allocate(cs%Iresttime_col_u(cs%num_col_u)) ; cs%Iresttime_col_u(:) = 0.0
280  allocate(cs%col_i_u(cs%num_col_u)) ; cs%col_i_u(:) = 0
281  allocate(cs%col_j_u(cs%num_col_u)) ; cs%col_j_u(:) = 0
282 
283  ! pass indices, restoring time to the CS structure
284  col = 1
285  do j=cs%jsc,cs%jec ; do i=cs%iscB,cs%iecB
286  if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) then
287  cs%col_i_u(col) = i ; cs%col_j_u(col) = j
288  cs%Iresttime_col_u(col) = iresttime_u(i,j)
289  col = col + 1
290  endif
291  enddo ; enddo
292 
293  ! same for total number of arbritary layers and correspondent data
294 
295  allocate(cs%Ref_hu%p(cs%nz_data,cs%num_col_u))
296  do col=1,cs%num_col_u ; do k=1,cs%nz_data
297  cs%Ref_hu%p(k,col) = data_hu(cs%col_i_u(col),cs%col_j_u(col),k)
298  enddo ; enddo
299  endif
300  total_sponge_cols_u = cs%num_col_u
301  call sum_across_pes(total_sponge_cols_u)
302  call log_param(param_file, mdl, "!Total sponge columns at u points", total_sponge_cols_u, &
303  "The total number of columns where sponges are applied at u points.", like_default=.true.)
304 
305  ! v points
306  cs%num_col_v = 0 ; !CS%fldno_v = 0
307  do j=cs%jscB,cs%jecB; do i=cs%isc,cs%iec
308  data_hv(i,j,:) = 0.5 * (data_h(i,j,:) + data_h(i,j+1,:))
309  iresttime_v(i,j) = 0.5 * (iresttime(i,j) + iresttime(i,j+1))
310  if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) &
311  cs%num_col_v = cs%num_col_v + 1
312  enddo ; enddo
313 
314  if (cs%num_col_v > 0) then
315 
316  allocate(cs%Iresttime_col_v(cs%num_col_v)) ; cs%Iresttime_col_v = 0.0
317  allocate(cs%col_i_v(cs%num_col_v)) ; cs%col_i_v = 0
318  allocate(cs%col_j_v(cs%num_col_v)) ; cs%col_j_v = 0
319 
320  ! pass indices, restoring time to the CS structure
321  col = 1
322  do j=cs%jscB,cs%jecB ; do i=cs%isc,cs%iec
323  if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) then
324  cs%col_i_v(col) = i ; cs%col_j_v(col) = j
325  cs%Iresttime_col_v(col) = iresttime_v(i,j)
326  col = col + 1
327  endif
328  enddo ; enddo
329 
330  ! same for total number of arbritary layers and correspondent data
331  allocate(cs%Ref_hv%p(cs%nz_data,cs%num_col_v))
332  do col=1,cs%num_col_v ; do k=1,cs%nz_data
333  cs%Ref_hv%p(k,col) = data_hv(cs%col_i_v(col),cs%col_j_v(col),k)
334  enddo ; enddo
335  endif
336  total_sponge_cols_v = cs%num_col_v
337  call sum_across_pes(total_sponge_cols_v)
338  call log_param(param_file, mdl, "!Total sponge columns at v points", total_sponge_cols_v, &
339  "The total number of columns where sponges are applied at v points.", like_default=.true.)
340  endif
341 

◆ initialize_ale_sponge_varying()

subroutine mom_ale_sponge::initialize_ale_sponge_varying ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  Iresttime,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(ale_sponge_cs), pointer  CS 
)
private

This subroutine determines the number of points which are to be restoref in the computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges.

Parameters
[in]gThe ocean's grid structure.
[in]iresttimeThe inverse of the restoring time [T-1 ~> s-1].
[in]param_fileA structure indicating the open file to parse for model parameter values.
csA pointer that is set to point to the control structure for this module (in/out).

Definition at line 396 of file MOM_ALE_sponge.F90.

397 
398  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
399  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: Iresttime !< The inverse of the restoring time [T-1 ~> s-1].
400  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to parse
401  !! for model parameter values.
402  type(ALE_sponge_CS), pointer :: CS !< A pointer that is set to point to the control
403  !! structure for this module (in/out).
404 
405 ! This include declares and sets the variable "version".
406 #include "version_variable.h"
407  character(len=40) :: mdl = "MOM_sponge" ! This module's name.
408  logical :: use_sponge
409  real, allocatable, dimension(:,:) :: Iresttime_u !< inverse of the restoring time at u points [T-1 ~> s-1]
410  real, allocatable, dimension(:,:) :: Iresttime_v !< inverse of the restoring time at v points [T-1 ~> s-1]
411  logical :: bndExtrapolation = .true. ! If true, extrapolate boundaries
412  logical :: default_2018_answers
413  logical :: spongeDataOngrid = .false.
414  integer :: i, j, k, col, total_sponge_cols, total_sponge_cols_u, total_sponge_cols_v
415  character(len=10) :: remapScheme
416 
417  if (associated(cs)) then
418  call mom_error(warning, "initialize_ALE_sponge_varying called with an associated "// &
419  "control structure.")
420  return
421  endif
422 ! Set default, read and log parameters
423  call log_version(param_file, mdl, version, "")
424  call get_param(param_file, mdl, "SPONGE", use_sponge, &
425  "If true, sponges may be applied anywhere in the domain. "//&
426  "The exact location and properties of those sponges are "//&
427  "specified from MOM_initialization.F90.", default=.false.)
428  if (.not.use_sponge) return
429  allocate(cs)
430  call get_param(param_file, mdl, "SPONGE_UV", cs%sponge_uv, &
431  "Apply sponges in u and v, in addition to tracers.", &
432  default=.false.)
433  call get_param(param_file, mdl, "REMAPPING_SCHEME", remapscheme, &
434  "This sets the reconstruction scheme used "//&
435  " for vertical remapping for all variables.", &
436  default="PLM", do_not_log=.true.)
437  call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", bndextrapolation, &
438  "When defined, a proper high-order reconstruction "//&
439  "scheme is used within boundary cells rather "//&
440  "than PCM. E.g., if PPM is used for remapping, a "//&
441  "PPM reconstruction will also be used within boundary cells.", &
442  default=.false., do_not_log=.true.)
443  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
444  "This sets the default value for the various _2018_ANSWERS parameters.", &
445  default=.false.)
446  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
447  "If true, use the order of arithmetic and expressions that recover the "//&
448  "answers from the end of 2018. Otherwise, use updated and more robust "//&
449  "forms of the same expressions.", default=default_2018_answers)
450  call get_param(param_file, mdl, "SPONGE_DATA_ONGRID", cs%spongeDataOngrid, &
451  "When defined, the incoming sponge data are "//&
452  "assumed to be on the model grid " , &
453  default=.false.)
454  cs%time_varying_sponges = .true.
455  cs%nz = g%ke
456  cs%isc = g%isc ; cs%iec = g%iec ; cs%jsc = g%jsc ; cs%jec = g%jec
457  cs%isd = g%isd ; cs%ied = g%ied ; cs%jsd = g%jsd ; cs%jed = g%jed
458  cs%iscB = g%iscB ; cs%iecB = g%iecB; cs%jscB = g%jscB ; cs%jecB = g%jecB
459 
460  ! number of columns to be restored
461  cs%num_col = 0 ; cs%fldno = 0
462  do j=g%jsc,g%jec ; do i=g%isc,g%iec
463  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
464  cs%num_col = cs%num_col + 1
465  enddo ; enddo
466  if (cs%num_col > 0) then
467  allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
468  allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
469  allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
470  ! pass indices, restoring time to the CS structure
471  col = 1
472  do j=g%jsc,g%jec ; do i=g%isc,g%iec
473  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) then
474  cs%col_i(col) = i ; cs%col_j(col) = j
475  cs%Iresttime_col(col) = iresttime(i,j)
476  col = col +1
477  endif
478  enddo ; enddo
479  endif
480  total_sponge_cols = cs%num_col
481  call sum_across_pes(total_sponge_cols)
482 
483 ! Call the constructor for remapping control structure
484  call initialize_remapping(cs%remap_cs, remapscheme, boundary_extrapolation=bndextrapolation, &
485  answers_2018=cs%remap_answers_2018)
486  call log_param(param_file, mdl, "!Total sponge columns at h points", total_sponge_cols, &
487  "The total number of columns where sponges are applied at h points.", like_default=.true.)
488  if (cs%sponge_uv) then
489  allocate(iresttime_u(g%isdB:g%iedB,g%jsd:g%jed)) ; iresttime_u(:,:) = 0.0
490  allocate(iresttime_v(g%isd:g%ied,g%jsdB:g%jedB)) ; iresttime_v(:,:) = 0.0
491  ! u points
492  cs%num_col_u = 0 ; !CS%fldno_u = 0
493  do j=cs%jsc,cs%jec; do i=cs%iscB,cs%iecB
494  iresttime_u(i,j) = 0.5 * (iresttime(i,j) + iresttime(i+1,j))
495  if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) &
496  cs%num_col_u = cs%num_col_u + 1
497  enddo ; enddo
498  if (cs%num_col_u > 0) then
499  allocate(cs%Iresttime_col_u(cs%num_col_u)) ; cs%Iresttime_col_u = 0.0
500  allocate(cs%col_i_u(cs%num_col_u)) ; cs%col_i_u = 0
501  allocate(cs%col_j_u(cs%num_col_u)) ; cs%col_j_u = 0
502  ! pass indices, restoring time to the CS structure
503  col = 1
504  do j=cs%jsc,cs%jec ; do i=cs%iscB,cs%iecB
505  if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) then
506  cs%col_i_u(col) = i ; cs%col_j_u(col) = j
507  cs%Iresttime_col_u(col) = iresttime_u(i,j)
508  col = col +1
509  endif
510  enddo ; enddo
511  ! same for total number of arbritary layers and correspondent data
512  endif
513  total_sponge_cols_u = cs%num_col_u
514  call sum_across_pes(total_sponge_cols_u)
515  call log_param(param_file, mdl, "!Total sponge columns at u points", total_sponge_cols_u, &
516  "The total number of columns where sponges are applied at u points.", like_default=.true.)
517  ! v points
518  cs%num_col_v = 0 ; !CS%fldno_v = 0
519  do j=cs%jscB,cs%jecB; do i=cs%isc,cs%iec
520  iresttime_v(i,j) = 0.5 * (iresttime(i,j) + iresttime(i,j+1))
521  if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) &
522  cs%num_col_v = cs%num_col_v + 1
523  enddo ; enddo
524  if (cs%num_col_v > 0) then
525  allocate(cs%Iresttime_col_v(cs%num_col_v)) ; cs%Iresttime_col_v = 0.0
526  allocate(cs%col_i_v(cs%num_col_v)) ; cs%col_i_v = 0
527  allocate(cs%col_j_v(cs%num_col_v)) ; cs%col_j_v = 0
528  ! pass indices, restoring time to the CS structure
529  col = 1
530  do j=cs%jscB,cs%jecB ; do i=cs%isc,cs%iec
531  if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) then
532  cs%col_i_v(col) = i ; cs%col_j_v(col) = j
533  cs%Iresttime_col_v(col) = iresttime_v(i,j)
534  col = col +1
535  endif
536  enddo ; enddo
537  endif
538  total_sponge_cols_v = cs%num_col_v
539  call sum_across_pes(total_sponge_cols_v)
540  call log_param(param_file, mdl, "!Total sponge columns at v points", total_sponge_cols_v, &
541  "The total number of columns where sponges are applied at v points.", like_default=.true.)
542  endif
543 

◆ rotate_ale_sponge()

subroutine, public mom_ale_sponge::rotate_ale_sponge ( type(ale_sponge_cs), intent(in)  sponge_in,
type(ocean_grid_type), intent(in)  G_in,
type(ale_sponge_cs), pointer  sponge,
type(ocean_grid_type), intent(in)  G,
integer, intent(in)  turns,
type(param_file_type), intent(in)  param_file 
)

Rotate the ALE sponge fields from the input to the model index map.

Parameters
[in]sponge_inThe control structure for this module with the original grid rotation
[in]g_inThe ocean's grid structure with the original rotation.
spongeA pointer to the control that will be set up with the new grid rotation
[in]gThe ocean's grid structure with the new rotation.
[in]turnsThe number of 90-degree turns between grids
[in]param_fileA structure indicating the open file to parse for model parameter values.

Definition at line 1005 of file MOM_ALE_sponge.F90.

1006  type(ALE_sponge_CS), intent(in) :: sponge_in !< The control structure for this module with the
1007  !! original grid rotation
1008  type(ocean_grid_type), intent(in) :: G_in !< The ocean's grid structure with the original rotation.
1009  type(ALE_sponge_CS), pointer :: sponge !< A pointer to the control that will be set up with
1010  !! the new grid rotation
1011  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure with the new rotation.
1012  integer, intent(in) :: turns !< The number of 90-degree turns between grids
1013  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
1014  !! to parse for model parameter values.
1015 
1016  ! First part: Index construction
1017  ! 1. Reconstruct Iresttime(:,:) from sponge_in
1018  ! 2. rotate Iresttime(:,:)
1019  ! 3. Call initialize_ALE_sponge using new grid and rotated Iresttime(:,:)
1020  ! All the index adjustment should follow from the Iresttime rotation
1021 
1022  real, dimension(:,:), allocatable :: Iresttime_in, Iresttime
1023  real, dimension(:,:,:), allocatable :: data_h_in, data_h
1024  real, dimension(:,:,:), allocatable :: sp_val_in, sp_val
1025  real, dimension(:,:,:), pointer :: sp_ptr => null()
1026  integer :: c, c_i, c_j
1027  integer :: k, nz_data
1028  integer :: n
1029  logical :: fixed_sponge
1030 
1031  fixed_sponge = .not. sponge_in%time_varying_sponges
1032  ! NOTE: nz_data is only conditionally set when fixed_sponge is true.
1033 
1034  allocate(iresttime_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed))
1035  allocate(iresttime(g%isd:g%ied, g%jsd:g%jed))
1036  iresttime_in(:,:) = 0.0
1037 
1038  if (fixed_sponge) then
1039  nz_data = sponge_in%nz_data
1040  allocate(data_h_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed, nz_data))
1041  allocate(data_h(g%isd:g%ied, g%jsd:g%jed, nz_data))
1042  data_h_in(:,:,:) = 0.
1043  endif
1044 
1045  ! Re-populate the 2D Iresttime and data_h arrays on the original grid
1046  do c=1,sponge_in%num_col
1047  c_i = sponge_in%col_i(c)
1048  c_j = sponge_in%col_j(c)
1049  iresttime_in(c_i, c_j) = sponge_in%Iresttime_col(c)
1050  if (fixed_sponge) then
1051  do k = 1, nz_data
1052  data_h(c_i, c_j, k) = sponge_in%Ref_h%p(k,c)
1053  enddo
1054  endif
1055  enddo
1056 
1057  call rotate_array(iresttime_in, turns, iresttime)
1058  if (fixed_sponge) then
1059  call rotate_array(data_h_in, turns, data_h)
1060  call initialize_ale_sponge_fixed(iresttime, g, param_file, sponge, &
1061  data_h, nz_data)
1062  else
1063  call initialize_ale_sponge_varying(iresttime, g, param_file, sponge)
1064  endif
1065 
1066  deallocate(iresttime_in)
1067  deallocate(iresttime)
1068  if (fixed_sponge) then
1069  deallocate(data_h_in)
1070  deallocate(data_h)
1071  endif
1072 
1073  ! Second part: Provide rotated fields for which relaxation is applied
1074 
1075  sponge%fldno = sponge_in%fldno
1076 
1077  if (fixed_sponge) then
1078  allocate(sp_val_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed, nz_data))
1079  allocate(sp_val(g%isd:g%ied, g%jsd:g%jed, nz_data))
1080  endif
1081 
1082  do n=1,sponge_in%fldno
1083  ! Assume that tracers are pointers and are remapped in other functions(?)
1084  sp_ptr => sponge_in%var(n)%p
1085  if (fixed_sponge) then
1086  sp_val_in(:,:,:) = 0.0
1087  do c = 1, sponge_in%num_col
1088  c_i = sponge_in%col_i(c)
1089  c_j = sponge_in%col_j(c)
1090  do k = 1, nz_data
1091  sp_val_in(c_i, c_j, k) = sponge_in%Ref_val(n)%p(k,c)
1092  enddo
1093  enddo
1094 
1095  call rotate_array(sp_val_in, turns, sp_val)
1096 
1097  ! NOTE: This points sp_val with the unrotated field. See note below.
1098  call set_up_ale_sponge_field(sp_val, g, sp_ptr, sponge)
1099 
1100  deallocate(sp_val_in)
1101  else
1102  ! We don't want to repeat FMS init in set_up_ALE_sponge_field_varying()
1103  ! (time_interp_external_init, init_external_field, etc), so we manually
1104  ! do a portion of this function below.
1105  sponge%Ref_val(n)%id = sponge_in%Ref_val(n)%id
1106  sponge%Ref_val(n)%num_tlevs = sponge_in%Ref_val(n)%num_tlevs
1107 
1108  nz_data = sponge_in%Ref_val(n)%nz_data
1109  sponge%Ref_val(n)%nz_data = nz_data
1110 
1111  allocate(sponge%Ref_val(n)%p(nz_data, sponge_in%num_col))
1112  allocate(sponge%Ref_val(n)%h(nz_data, sponge_in%num_col))
1113  sponge%Ref_val(n)%p(:,:) = 0.0
1114  sponge%Ref_val(n)%h(:,:) = 0.0
1115 
1116  ! TODO: There is currently no way to associate a generic field pointer to
1117  ! its rotated equivalent without introducing a new data structure which
1118  ! explicitly tracks the pairing.
1119  !
1120  ! As a temporary fix, we store the pointer to the unrotated field in
1121  ! the rotated sponge, and use this reference to replace the pointer
1122  ! to the rotated field update_ALE_sponge field.
1123  !
1124  ! This makes a lot of unverifiable assumptions, and should not be
1125  ! considered the final solution.
1126  sponge%var(n)%p => sp_ptr
1127  endif
1128  enddo
1129 
1130  ! TODO: var_u and var_v sponge dampling is not yet supported.
1131  if (associated(sponge_in%var_u%p) .or. associated(sponge_in%var_v%p)) &
1132  call mom_error(fatal, "Rotation of ALE sponge velocities is not yet " &
1133  // "implemented.")
1134 
1135  ! Transfer any existing diag_CS reference pointer
1136  sponge%diag => sponge_in%diag
1137 
1138  ! NOTE: initialize_ALE_sponge_* resolves remap_cs

◆ set_up_ale_sponge_field_fixed()

subroutine mom_ale_sponge::set_up_ale_sponge_field_fixed ( real, dimension( g %isd: g %ied, g %jsd: g %jed,cs%nz_data), intent(in)  sp_val,
type(ocean_grid_type), intent(in)  G,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in), target  f_ptr,
type(ale_sponge_cs), pointer  CS 
)
private

This subroutine stores the reference profile at h points for the variable whose address is given by f_ptr.

Parameters
[in]gGrid structure
csALE sponge control structure (in/out).
[in]sp_valField to be used in the sponge, it has arbitrary number of layers.
[in]f_ptrPointer to the field to be damped

Definition at line 563 of file MOM_ALE_sponge.F90.

564  type(ocean_grid_type), intent(in) :: G !< Grid structure
565  type(ALE_sponge_CS), pointer :: CS !< ALE sponge control structure (in/out).
566  real, dimension(SZI_(G),SZJ_(G),CS%nz_data), &
567  intent(in) :: sp_val !< Field to be used in the sponge, it has arbitrary number of layers.
568  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
569  target, intent(in) :: f_ptr !< Pointer to the field to be damped
570 
571  integer :: j, k, col
572  character(len=256) :: mesg ! String for error messages
573 
574  if (.not.associated(cs)) return
575 
576  cs%fldno = cs%fldno + 1
577 
578  if (cs%fldno > max_fields_) then
579  write(mesg,'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
580  &the number of fields to be damped in the call to &
581  &initialize_ALE_sponge." )') cs%fldno
582  call mom_error(fatal,"set_up_ALE_sponge_field: "//mesg)
583  endif
584 
585  ! stores the reference profile
586  allocate(cs%Ref_val(cs%fldno)%p(cs%nz_data,cs%num_col))
587  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
588  do col=1,cs%num_col
589  do k=1,cs%nz_data
590  cs%Ref_val(cs%fldno)%p(k,col) = sp_val(cs%col_i(col),cs%col_j(col),k)
591  enddo
592  enddo
593 
594  cs%var(cs%fldno)%p => f_ptr
595 

◆ set_up_ale_sponge_field_varying()

subroutine mom_ale_sponge::set_up_ale_sponge_field_varying ( character(len=*), intent(in)  filename,
character(len=*), intent(in)  fieldname,
type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in), target  f_ptr,
type(ale_sponge_cs), pointer  CS 
)
private

This subroutine stores the reference profile at h points for the variable whose address is given by filename and fieldname.

Parameters
[in]filenameThe name of the file with the time varying field data
[in]fieldnameThe name of the field in the file with the time varying field data
[in]timeThe current model time
[in]gGrid structure (in).
[in]gvocean vertical grid structure
[in]usA dimensional unit scaling type
[in]f_ptrPointer to the field to be damped (in).
csSponge control structure (in/out).

Definition at line 600 of file MOM_ALE_sponge.F90.

601  character(len=*), intent(in) :: filename !< The name of the file with the
602  !! time varying field data
603  character(len=*), intent(in) :: fieldname !< The name of the field in the file
604  !! with the time varying field data
605  type(time_type), intent(in) :: Time !< The current model time
606  type(ocean_grid_type), intent(in) :: G !< Grid structure (in).
607  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
608  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
609  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
610  target, intent(in) :: f_ptr !< Pointer to the field to be damped (in).
611  type(ALE_sponge_CS), pointer :: CS !< Sponge control structure (in/out).
612 
613  ! Local variables
614  real, allocatable, dimension(:,:,:) :: sp_val !< Field to be used in the sponge
615  real, allocatable, dimension(:,:,:) :: mask_z !< Field mask for the sponge data
616  real, allocatable, dimension(:), target :: z_in, z_edges_in ! Heights [Z ~> m].
617  real :: missing_value
618  integer :: j, k, col
619  integer :: isd,ied,jsd,jed
620  integer :: nPoints
621  integer, dimension(4) :: fld_sz
622  integer :: nz_data !< the number of vertical levels in this input field
623  character(len=256) :: mesg ! String for error messages
624  ! Local variables for ALE remapping
625  real, dimension(:), allocatable :: tmpT1d
626  real :: zTopOfCell, zBottomOfCell ! Heights [Z ~> m].
627  type(remapping_CS) :: remapCS ! Remapping parameters and work arrays
628 
629  if (.not.associated(cs)) return
630  ! initialize time interpolator module
631  call time_interp_external_init()
632  isd = g%isd; ied = g%ied; jsd = g%jsd; jed = g%jed
633  cs%fldno = cs%fldno + 1
634  if (cs%fldno > max_fields_) then
635  write(mesg,'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
636  &the number of fields to be damped in the call to &
637  &initialize_ALE_sponge." )') cs%fldno
638  call mom_error(fatal,"set_up_ALE_sponge_field: "//mesg)
639  endif
640  ! get a unique time interp id for this field. If sponge data is ongrid, then setup
641  ! to only read on the computational domain
642  if (cs%spongeDataOngrid) then
643  cs%Ref_val(cs%fldno)%id = init_external_field(filename, fieldname,domain=g%Domain%mpp_domain)
644  else
645  cs%Ref_val(cs%fldno)%id = init_external_field(filename, fieldname)
646  endif
647  fld_sz(1:4)=-1
648  fld_sz = get_external_field_size(cs%Ref_val(cs%fldno)%id)
649  nz_data = fld_sz(3)
650  cs%Ref_val(cs%fldno)%nz_data = nz_data !< individual sponge fields may reside on a different vertical grid
651  cs%Ref_val(cs%fldno)%num_tlevs = fld_sz(4)
652  ! initializes the target profile array for this field
653  ! for all columns which will be masked
654  allocate(cs%Ref_val(cs%fldno)%p(nz_data,cs%num_col))
655  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
656  allocate( cs%Ref_val(cs%fldno)%h(nz_data,cs%num_col) )
657  cs%Ref_val(cs%fldno)%h(:,:) = 0.0
658  cs%var(cs%fldno)%p => f_ptr
659 
660 

◆ set_up_ale_sponge_vel_field_fixed()

subroutine mom_ale_sponge::set_up_ale_sponge_vel_field_fixed ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed,cs%nz_data), intent(in)  u_val,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb,cs%nz_data), intent(in)  v_val,
type(ocean_grid_type), intent(in)  G,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in), target  u_ptr,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in), target  v_ptr,
type(ale_sponge_cs), pointer  CS 
)
private

This subroutine stores the reference profile at u and v points for the variable whose address is given by u_ptr and v_ptr.

Parameters
[in]gGrid structure (in).
csSponge structure (in/out).
[in]u_valu field to be used in the sponge, it has arbritary number of layers.
[in]v_valv field to be used in the sponge, it has arbritary number of layers.
[in]u_ptru pointer to the field to be damped
[in]v_ptrv pointer to the field to be damped

Definition at line 665 of file MOM_ALE_sponge.F90.

666  type(ocean_grid_type), intent(in) :: G !< Grid structure (in).
667  type(ALE_sponge_CS), pointer :: CS !< Sponge structure (in/out).
668  real, dimension(SZIB_(G),SZJ_(G),CS%nz_data), &
669  intent(in) :: u_val !< u field to be used in the sponge, it has arbritary number of layers.
670  real, dimension(SZI_(G),SZJB_(G),CS%nz_data), &
671  intent(in) :: v_val !< v field to be used in the sponge, it has arbritary number of layers.
672  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), target, intent(in) :: u_ptr !< u pointer to the field to be damped
673  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), target, intent(in) :: v_ptr !< v pointer to the field to be damped
674 
675  integer :: j, k, col
676  character(len=256) :: mesg ! String for error messages
677 
678  if (.not.associated(cs)) return
679 
680  ! stores the reference profile
681  allocate(cs%Ref_val_u%p(cs%nz_data,cs%num_col_u))
682  cs%Ref_val_u%p(:,:) = 0.0
683  do col=1,cs%num_col_u
684  do k=1,cs%nz_data
685  cs%Ref_val_u%p(k,col) = u_val(cs%col_i_u(col),cs%col_j_u(col),k)
686  enddo
687  enddo
688  cs%var_u%p => u_ptr
689  allocate(cs%Ref_val_v%p(cs%nz_data,cs%num_col_v))
690  cs%Ref_val_v%p(:,:) = 0.0
691  do col=1,cs%num_col_v
692  do k=1,cs%nz_data
693  cs%Ref_val_v%p(k,col) = v_val(cs%col_i_v(col),cs%col_j_v(col),k)
694  enddo
695  enddo
696  cs%var_v%p => v_ptr
697 

◆ set_up_ale_sponge_vel_field_varying()

subroutine mom_ale_sponge::set_up_ale_sponge_vel_field_varying ( character(len=*), intent(in)  filename_u,
character(len=*), intent(in)  fieldname_u,
character(len=*), intent(in)  filename_v,
character(len=*), intent(in)  fieldname_v,
type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(ale_sponge_cs), pointer  CS,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in), target  u_ptr,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in), target  v_ptr 
)
private

This subroutine stores the reference profile at uand v points for the variable whose address is given by u_ptr and v_ptr.

Parameters
[in]filename_uFile name for u field
[in]fieldname_uName of u variable in file
[in]filename_vFile name for v field
[in]fieldname_vName of v variable in file
[in]timeModel time
[in,out]gOcean grid (in)
[in]usA dimensional unit scaling type
csSponge structure (in/out).
[in]u_ptru pointer to the field to be damped (in).
[in]v_ptrv pointer to the field to be damped (in).

Definition at line 702 of file MOM_ALE_sponge.F90.

704  character(len=*), intent(in) :: filename_u !< File name for u field
705  character(len=*), intent(in) :: fieldname_u !< Name of u variable in file
706  character(len=*), intent(in) :: filename_v !< File name for v field
707  character(len=*), intent(in) :: fieldname_v !< Name of v variable in file
708  type(time_type), intent(in) :: Time !< Model time
709  type(ocean_grid_type), intent(inout) :: G !< Ocean grid (in)
710  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
711  type(ALE_sponge_CS), pointer :: CS !< Sponge structure (in/out).
712  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), target, intent(in) :: u_ptr !< u pointer to the field to be damped (in).
713  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), target, intent(in) :: v_ptr !< v pointer to the field to be damped (in).
714  ! Local variables
715  real, allocatable, dimension(:,:,:) :: u_val !< U field to be used in the sponge.
716  real, allocatable, dimension(:,:,:) :: mask_u !< U field mask for the sponge data.
717  real, allocatable, dimension(:,:,:) :: v_val !< V field to be used in the sponge.
718  real, allocatable, dimension(:,:,:) :: mask_v !< V field mask for the sponge data.
719 
720  real, allocatable, dimension(:), target :: z_in, z_edges_in
721  real :: missing_value
722 
723  integer :: j, k, col
724  integer :: isd, ied, jsd, jed
725  integer :: isdB, iedB, jsdB, jedB
726  integer, dimension(4) :: fld_sz
727  character(len=256) :: mesg ! String for error messages
728 
729  if (.not.associated(cs)) return
730 
731  isd = g%isd; ied = g%ied; jsd = g%jsd; jed = g%jed
732  isdb = g%isdB; iedb = g%iedB; jsdb = g%jsdB; jedb = g%jedB
733  ! get a unique id for this field which will allow us to return an array
734  ! containing time-interpolated values from an external file corresponding
735  ! to the current model date.
736  cs%Ref_val_u%id = init_external_field(filename_u, fieldname_u)
737  fld_sz(1:4)=-1
738  fld_sz = get_external_field_size(cs%Ref_val_u%id)
739  cs%Ref_val_u%nz_data = fld_sz(3)
740  cs%Ref_val_u%num_tlevs = fld_sz(4)
741  cs%Ref_val_v%id = init_external_field(filename_v, fieldname_v)
742  fld_sz(1:4)=-1
743  fld_sz = get_external_field_size(cs%Ref_val_v%id)
744  cs%Ref_val_v%nz_data = fld_sz(3)
745  cs%Ref_val_v%num_tlevs = fld_sz(4)
746  allocate( u_val(isdb:iedb,jsd:jed, fld_sz(3)) )
747  allocate( mask_u(isdb:iedb,jsd:jed, fld_sz(3)) )
748  allocate( v_val(isd:ied,jsdb:jedb, fld_sz(3)) )
749  allocate( mask_v(isd:ied,jsdb:jedb, fld_sz(3)) )
750  ! Interpolate external file data to the model grid
751  ! I am hard-wiring this call to assume that the input grid is zonally re-entrant
752  ! In the future, this should be generalized using an interface to return the
753  ! modulo attribute of the zonal axis (mjh).
754  call horiz_interp_and_extrap_tracer(cs%Ref_val_u%id, time, 1.0, g, u_val, mask_u, z_in, z_edges_in, &
755  missing_value, .true., .false., .false., m_to_z=us%m_to_Z, &
756  answers_2018=cs%hor_regrid_answers_2018)
757  !!! TODO: add a velocity interface! (mjh)
758  ! Interpolate external file data to the model grid
759  ! I am hard-wiring this call to assume that the input grid is zonally re-entrant
760  ! In the future, this should be generalized using an interface to return the
761  ! modulo attribute of the zonal axis (mjh).
762  call horiz_interp_and_extrap_tracer(cs%Ref_val_v%id, time, 1.0, g, v_val, mask_v, z_in, z_edges_in, &
763  missing_value, .true., .false., .false., m_to_z=us%m_to_Z, &
764  answers_2018=cs%hor_regrid_answers_2018)
765  ! stores the reference profile
766  allocate(cs%Ref_val_u%p(fld_sz(3),cs%num_col_u))
767  cs%Ref_val_u%p(:,:) = 0.0
768  do col=1,cs%num_col_u
769  do k=1,fld_sz(3)
770  cs%Ref_val_u%p(k,col) = u_val(cs%col_i_u(col),cs%col_j_u(col),k)
771  enddo
772  enddo
773  cs%var_u%p => u_ptr
774  allocate(cs%Ref_val_v%p(fld_sz(3),cs%num_col_v))
775  cs%Ref_val_v%p(:,:) = 0.0
776  do col=1,cs%num_col_v
777  do k=1,fld_sz(3)
778  cs%Ref_val_v%p(k,col) = v_val(cs%col_i_v(col),cs%col_j_v(col),k)
779  enddo
780  enddo
781  cs%var_v%p => v_ptr
782 

◆ update_ale_sponge_field()

subroutine, public mom_ale_sponge::update_ale_sponge_field ( type(ale_sponge_cs), pointer  sponge,
real, dimension(:,:,:), intent(in), target  p_old,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in), target  p_new 
)

Scan the ALE sponge variables and replace a prescribed pointer to a new value.

Parameters
spongeA pointer to the control structure for this module that is set by a previous call to initialize_ALE_sponge.
[in]p_oldThe previous array of target values
[in]gThe updated ocean grid structure
[in]gvocean vertical grid structure
[in]p_newThe new array of target values

Definition at line 1146 of file MOM_ALE_sponge.F90.

1147  type(ALE_sponge_CS), pointer :: sponge !< A pointer to the control structure for this module
1148  !! that is set by a previous call to initialize_ALE_sponge.
1149  real, dimension(:,:,:), &
1150  target, intent(in) :: p_old !< The previous array of target values
1151  type(ocean_grid_type), intent(in) :: G !< The updated ocean grid structure
1152  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
1153  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1154  target, intent(in) :: p_new !< The new array of target values
1155 
1156  integer :: n
1157 
1158  do n=1,sponge%fldno
1159  if (associated(sponge%var(n)%p, p_old)) sponge%var(n)%p => p_new
1160  enddo
1161