SCALE-RM
Functions/Subroutines
scale_ocean_phy_slab Module Reference

module OCEAN / Physics Slab model More...

Functions/Subroutines

subroutine, public ocean_phy_slab_setup (OCEAN_TYPE)
 Setup. More...
 
subroutine, public ocean_phy_slab (OCEAN_TEMP_t, OCEAN_TEMP, OCEAN_SFLX_WH, OCEAN_SFLX_prec, OCEAN_SFLX_evap, dt)
 Slab ocean model. More...
 

Detailed Description

module OCEAN / Physics Slab model

Description
ocean physics module, slab model
Author
Team SCALE
NAMELIST
  • PARAM_OCEAN_PHY_SLAB
    nametypedefault valuecomment
    OCEAN_PHY_SLAB_DEPTH real(RP) 10.0_RP water depth of slab ocean [m]

History Output
No history output

Function/Subroutine Documentation

◆ ocean_phy_slab_setup()

subroutine, public scale_ocean_phy_slab::ocean_phy_slab_setup ( character(len=*), intent(in)  OCEAN_TYPE)

Setup.

Definition at line 54 of file scale_ocean_phy_slab.F90.

References scale_const::const_cl, scale_const::const_dwatr, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_landuse::landuse_fact_ocean, and scale_process::prc_mpistop().

Referenced by scale_ocean_phy::ocean_phy_setup().

54  use scale_process, only: &
56  use scale_const, only: &
57  dwatr => const_dwatr, &
58  cl => const_cl
59  use scale_landuse, only: &
61  implicit none
62 
63  character(len=*), intent(in) :: ocean_type
64 
65  namelist / param_ocean_phy_slab / &
66  ocean_phy_slab_depth
67 
68  integer :: i, j
69  integer :: ierr
70  !---------------------------------------------------------------------------
71 
72  if( io_l ) write(io_fid_log,*)
73  if( io_l ) write(io_fid_log,*) '++++++ Module[SLAB] / Categ[OCEAN PHY] / Origin[SCALElib]'
74 
75  !--- read namelist
76  rewind(io_fid_conf)
77  read(io_fid_conf,nml=param_ocean_phy_slab,iostat=ierr)
78  if( ierr < 0 ) then !--- missing
79  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
80  elseif( ierr > 0 ) then !--- fatal error
81  write(*,*) 'xxx Not appropriate names in namelist PARAM_OCEAN_PHY_SLAB. Check!'
82  call prc_mpistop
83  endif
84  if( io_lnml ) write(io_fid_log,nml=param_ocean_phy_slab)
85 
86  if( ocean_type == 'CONST' ) then
87  ocean_phy_slab_fixedsst = .true.
88  else if( ocean_type == 'SLAB' ) then
89  ocean_phy_slab_fixedsst = .false.
90  else
91  write(*,*) 'xxx wrong OCEAN_TYPE. Check!'
92  call prc_mpistop
93  end if
94 
95  ocean_phy_slab_heatcapacity = dwatr * cl * ocean_phy_slab_depth
96 
97  if( io_l ) write(io_fid_log,*)
98  if( io_l ) write(io_fid_log,*) '*** Prevent SST change? : ', ocean_phy_slab_fixedsst
99  if( io_l ) write(io_fid_log,*) '*** Slab ocean depth [m] : ', ocean_phy_slab_depth
100  if( io_l ) write(io_fid_log,*) '*** Ocean heat capacity [J/K/m2] : ', ocean_phy_slab_heatcapacity
101 
102  ! judge to run slab ocean model
103  allocate( is_ocn(ia,ja) )
104 
105  do j = js, je
106  do i = is, ie
107  if( landuse_fact_ocean(i,j) > 0.0_rp ) then
108  is_ocn(i,j) = .true.
109  else
110  is_ocn(i,j) = .false.
111  end if
112  end do
113  end do
114 
115  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:87
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
integer, public ia
of x whole cells (local, with HALO)
module LANDUSE
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
integer, public js
start point of inner domain: y, local
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ie
end point of inner domain: x, local
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ocean_phy_slab()

subroutine, public scale_ocean_phy_slab::ocean_phy_slab ( real(rp), dimension (ia,ja), intent(out)  OCEAN_TEMP_t,
real(rp), dimension (ia,ja), intent(in)  OCEAN_TEMP,
real(rp), dimension (ia,ja), intent(in)  OCEAN_SFLX_WH,
real(rp), dimension(ia,ja), intent(in)  OCEAN_SFLX_prec,
real(rp), dimension(ia,ja), intent(in)  OCEAN_SFLX_evap,
real(dp), intent(in)  dt 
)

Slab ocean model.

Definition at line 127 of file scale_ocean_phy_slab.F90.

References scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, and scale_grid_index::js.

Referenced by scale_ocean_phy::ocean_phy_setup().

127  use scale_grid_index
128  implicit none
129 
130  real(RP), intent(out) :: ocean_temp_t (ia,ja)
131  real(RP), intent(in) :: ocean_temp (ia,ja)
132  real(RP), intent(in) :: ocean_sflx_wh (ia,ja)
133  real(RP), intent(in) :: ocean_sflx_prec(ia,ja)
134  real(RP), intent(in) :: ocean_sflx_evap(ia,ja)
135  real(DP), intent(in) :: dt
136 
137  integer :: i, j
138  !---------------------------------------------------------------------------
139 
140  if( io_l ) write(io_fid_log,*) '*** Ocean step: Slab'
141 
142  if( ocean_phy_slab_fixedsst )then
143  do j = js, je
144  do i = is, ie
145  ocean_temp_t(i,j) = 0.0_rp
146  enddo
147  enddo
148  else
149  do j = js, je
150  do i = is, ie
151  if( is_ocn(i,j) ) then
152  ocean_temp_t(i,j) = - ocean_sflx_wh(i,j) / ocean_phy_slab_heatcapacity
153  else
154  ocean_temp_t(i,j) = 0.0_rp
155  endif
156  enddo
157  enddo
158  end if
159 
160  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
module grid index
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
integer, public ie
end point of inner domain: x, local
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function: