SCALE-RM
scale_ocean_phy_slab.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
18  use scale_prof
19  use scale_debug
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: ocean_phy_slab_setup
29  public :: ocean_phy_slab
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  !-----------------------------------------------------------------------------
36  !
37  !++ Private procedure
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private parameters & variables
42  !
43  real(RP), private :: ocean_phy_slab_depth = 10.0_rp
44  logical, private :: ocean_phy_slab_fixedsst = .false.
45  real(RP), private :: ocean_phy_slab_heatcapacity
46 
47  logical, allocatable, private :: is_ocn(:,:)
48 
49  !-----------------------------------------------------------------------------
50 contains
51  !-----------------------------------------------------------------------------
53  subroutine ocean_phy_slab_setup( OCEAN_TYPE )
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
116  end subroutine ocean_phy_slab_setup
117 
118  !-----------------------------------------------------------------------------
120  subroutine ocean_phy_slab( &
121  OCEAN_TEMP_t, &
122  OCEAN_TEMP, &
123  OCEAN_SFLX_WH, &
124  OCEAN_SFLX_prec, &
125  OCEAN_SFLX_evap, &
126  dt )
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
161  end subroutine ocean_phy_slab
162 
163 end module scale_ocean_phy_slab
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
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
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
module STDIO
Definition: scale_stdio.F90:12
subroutine, public ocean_phy_slab_setup(OCEAN_TYPE)
Setup.
module grid index
module OCEAN / Physics Slab model
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
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
subroutine, public ocean_phy_slab(OCEAN_TEMP_t, OCEAN_TEMP, OCEAN_SFLX_WH, OCEAN_SFLX_prec, OCEAN_SFLX_evap, dt)
Slab ocean model.
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)