SCALE-RM
Functions/Subroutines
scale_ocean_phy_file Module Reference

module OCEAN / Physics File More...

Functions/Subroutines

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

Detailed Description

module OCEAN / Physics File

Description
ocean physics module, external file input
Author
Team SCALE

Function/Subroutine Documentation

◆ ocean_phy_file_setup()

subroutine, public scale_ocean_phy_file::ocean_phy_file_setup ( character(len=*), intent(in)  OCEAN_TYPE)

Setup.

Definition at line 50 of file scale_ocean_phy_file.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, 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().

50  use scale_process, only: &
52  use scale_landuse, only: &
54  implicit none
55 
56  character(len=*), intent(in) :: ocean_type
57 
58  integer :: i, j
59  integer :: ierr
60  !---------------------------------------------------------------------------
61 
62  if( io_l ) write(io_fid_log,*)
63  if( io_l ) write(io_fid_log,*) '++++++ Module[FILE] / Categ[OCEAN PHY] / Origin[SCALElib]'
64 
65  if( ocean_type /= 'FILE' ) then
66  write(*,*) 'xxx wrong OCEAN_TYPE. Check!'
67  call prc_mpistop
68  end if
69 
70  ! judge to run slab ocean model
71  allocate( is_ocn(ia,ja) )
72 
73  do j = js, je
74  do i = is, ie
75  if( landuse_fact_ocean(i,j) > 0.0_rp ) then
76  is_ocn(i,j) = .true.
77  else
78  is_ocn(i,j) = .false.
79  end if
80  end do
81  end do
82 
83  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.
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
integer, public ie
end point of inner domain: x, local
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_file()

subroutine, public scale_ocean_phy_file::ocean_phy_file ( 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 95 of file scale_ocean_phy_file.F90.

References scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_process::prc_mpistop(), and scale_time::time_nowdaysec.

Referenced by scale_ocean_phy::ocean_phy_setup().

96  use scale_time, only: &
97  nowdaysec => time_nowdaysec
98  use scale_process, only: &
100  use scale_external_input, only: &
101  extin_update
102  implicit none
103 
104  real(RP), intent(out) :: ocean_temp_t (ia,ja)
105  real(RP), intent(in) :: ocean_temp (ia,ja)
106  real(RP), intent(in) :: ocean_sflx_wh (ia,ja)
107  real(RP), intent(in) :: ocean_sflx_prec(ia,ja)
108  real(RP), intent(in) :: ocean_sflx_evap(ia,ja)
109  real(DP), intent(in) :: dt
110 
111  real(RP) :: ocean_temp_new(ia,ja)
112 
113  logical :: error
114 
115  integer :: i, j
116  !---------------------------------------------------------------------------
117 
118  if( io_l ) write(io_fid_log,*) '*** Ocean step: File'
119 
120  call extin_update( &
121  ocean_temp_new, & ! (out)
122  'OCEAN_TEMP', & ! (in)
123  'XY', & ! (in)
124  nowdaysec, & ! (in)
125  error ) ! (out)
126  if ( error ) then
127  write(*,*) 'xxx Failed to read "OCEAN_TEMP" data from file!'
128  write(*,*) 'xxx set namelist EXTITEM, for example'
129  write(*,*) ' &EXTITEM'
130  write(*,*) ' basename = "filename",'
131  write(*,*) ' varname = "OCEAN_TEMP",'
132  write(*,*) ' /'
133  call prc_mpistop
134  end if
135 
136  do j = js, je
137  do i = is, ie
138  if( is_ocn(i,j) ) then
139  ocean_temp_t(i,j) = ( ocean_temp_new(i,j) - ocean_temp(i,j) ) / dt
140  else
141  ocean_temp_t(i,j) = 0.0_rp
142  endif
143  enddo
144  enddo
145 
146  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(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
module grid index
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
module EXTERNAL INPUT
integer, public ie
end point of inner domain: x, local
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: