Skip to content
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
86 changes: 45 additions & 41 deletions src/parameterizations/stochastic/MOM_stochastics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -11,12 +11,12 @@ module MOM_stochastics
use MOM_diag_mediator, only : register_diag_field, diag_ctrl, time_type
use MOM_grid, only : ocean_grid_type
use MOM_verticalGrid, only : verticalGrid_type
use MOM_error_handler, only : MOM_error, FATAL, WARNING, is_root_pe
use MOM_error_handler, only : MOM_error, MOM_mesg, FATAL, WARNING, is_root_pe
use MOM_error_handler, only : callTree_enter, callTree_leave
use MOM_file_parser, only : get_param, log_version, close_param_file, param_file_type
use mpp_domains_mod, only : domain2d, mpp_get_layout, mpp_get_global_domain
use mpp_domains_mod, only : mpp_define_domains, mpp_get_compute_domain, mpp_get_data_domain
use MOM_domains, only : root_PE,num_PEs
use MOM_domains, only : root_PE, num_PEs
use MOM_coms, only : Get_PElist
use stochastic_physics, only : init_stochastic_physics_ocn, run_stochastic_physics_ocn

Expand All @@ -31,8 +31,8 @@ module MOM_stochastics
logical :: do_sppt !< If true, stochastically perturb the diabatic
logical :: pert_epbl !< If true, then randomly perturb the KE dissipation and genration terms
integer :: id_sppt_wts = -1 !< Diagnostic id for SPPT
integer :: id_epbl1_wts=-1 !< Diagnostic id for epbl generation perturbation
integer :: id_epbl2_wts=-1 !< Diagnostic id for epbl dissipation perturbation
integer :: id_epbl1_wts = -1 !< Diagnostic id for epbl generation perturbation
integer :: id_epbl2_wts = -1 !< Diagnostic id for epbl dissipation perturbation
! stochastic patterns
real, allocatable :: sppt_wts(:,:) !< Random pattern for ocean SPPT
!! tendencies with a number between 0 and 2
Expand All @@ -46,24 +46,25 @@ module MOM_stochastics

!! This subroutine initializes the stochastics physics control structure.
subroutine stochastics_init(dt, grid, GV, CS, param_file, diag, Time)
real, intent(in) :: dt !< time step [T ~> s]
type(ocean_grid_type), intent(in) :: grid !< horizontal grid information
type(verticalGrid_type), intent(in) :: GV !< vertical grid structure
type(stochastic_CS), pointer, intent(inout):: CS !< stochastic control structure
real, intent(in) :: dt !< time step [T ~> s]
type(ocean_grid_type), intent(in) :: grid !< horizontal grid information
type(verticalGrid_type), intent(in) :: GV !< vertical grid structure
type(stochastic_CS), pointer, intent(inout) :: CS !< stochastic control structure
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
type(diag_ctrl), target, intent(inout) :: diag !< structure to regulate diagnostic output
type(time_type), target :: Time !< model time
type(diag_ctrl), target, intent(inout) :: diag !< structure to regulate diagnostic output
type(time_type), target :: Time !< model time

! Local variables
integer,allocatable :: pelist(:) ! list of pes for this instance of the ocean
integer, allocatable :: pelist(:) ! list of pes for this instance of the ocean
integer :: mom_comm ! list of pes for this instance of the ocean
integer :: num_procs ! number of processors to pass to stochastic physics
integer :: iret ! return code from stochastic physics
integer :: pe_zero ! root pe
integer :: nx ! number of x-points including halo
integer :: ny ! number of x-points including halo

! This include declares and sets the variable "version".
#include "version_variable.h"
! This include declares and sets the variable "version".
# include "version_variable.h"
character(len=40) :: mdl = "ocean_stochastics_init" ! This module's name.

call callTree_enter("ocean_model_stochastic_init(), MOM_stochastics.F90")
Expand All @@ -79,7 +80,7 @@ subroutine stochastics_init(dt, grid, GV, CS, param_file, diag, Time)
! Read all relevant parameters and write them to the model log.
call log_version(param_file, mdl, version, "")

! get number of processors and PE list for stocasthci physics initialization
! get number of processors and PE list for stochastic physics initialization
call get_param(param_file, mdl, "DO_SPPT", CS%do_sppt, &
"If true, then stochastically perturb the thermodynamic "//&
"tendemcies of T,S, amd h. Amplitude and correlations are "//&
Expand All @@ -91,37 +92,40 @@ subroutine stochastics_init(dt, grid, GV, CS, param_file, diag, Time)
"controlled by the nam_stoch namelist in the UFS model only.", &
default=.false.)
if (CS%do_sppt .OR. CS%pert_epbl) then
num_procs=num_PEs()
allocate(pelist(num_procs))
call Get_PElist(pelist,commID = mom_comm)
pe_zero=root_PE()
nx = grid%ied - grid%isd + 1
ny = grid%jed - grid%jsd + 1
call init_stochastic_physics_ocn(dt,grid%geoLonT,grid%geoLatT,nx,ny,GV%ke, &
CS%pert_epbl,CS%do_sppt,pe_zero,mom_comm,iret)
if (iret/=0) then
call MOM_error(FATAL, "call to init_stochastic_physics_ocn failed")
return
endif

if (CS%do_sppt) allocate(CS%sppt_wts(grid%isd:grid%ied,grid%jsd:grid%jed))
if (CS%pert_epbl) then
allocate(CS%epbl1_wts(grid%isd:grid%ied,grid%jsd:grid%jed))
allocate(CS%epbl2_wts(grid%isd:grid%ied,grid%jsd:grid%jed))
endif
num_procs = num_PEs()
allocate(pelist(num_procs))
call Get_PElist(pelist,commID = mom_comm)
pe_zero = root_PE()
nx = grid%ied - grid%isd + 1
ny = grid%jed - grid%jsd + 1
call init_stochastic_physics_ocn(dt,grid%geoLonT,grid%geoLatT,nx,ny,GV%ke, &
CS%pert_epbl,CS%do_sppt,pe_zero,mom_comm,iret)
if (iret/=0) then
call MOM_error(FATAL, "call to init_stochastic_physics_ocn failed")
endif

if (CS%do_sppt) allocate(CS%sppt_wts(grid%isd:grid%ied,grid%jsd:grid%jed), source=0.0)
if (CS%pert_epbl) then
allocate(CS%epbl1_wts(grid%isd:grid%ied,grid%jsd:grid%jed), source=0.0)
allocate(CS%epbl2_wts(grid%isd:grid%ied,grid%jsd:grid%jed), source=0.0)
endif
endif
if (CS%do_sppt) then
CS%id_sppt_wts = register_diag_field('ocean_model', 'sppt_pattern', CS%diag%axesT1, Time, &
'random pattern for sppt', 'None')
endif
if (CS%pert_epbl) then
CS%id_epbl1_wts = register_diag_field('ocean_model', 'epbl1_wts', CS%diag%axesT1, Time, &
'random pattern for KE generation', 'None')
CS%id_epbl2_wts = register_diag_field('ocean_model', 'epbl2_wts', CS%diag%axesT1, Time, &
'random pattern for KE dissipation', 'None')
endif
CS%id_sppt_wts = register_diag_field('ocean_model', 'sppt_pattern', CS%diag%axesT1, Time, &
'random pattern for sppt', 'None')
CS%id_epbl1_wts = register_diag_field('ocean_model', 'epbl1_wts', CS%diag%axesT1, Time, &
'random pattern for KE generation', 'None')
CS%id_epbl2_wts = register_diag_field('ocean_model', 'epbl2_wts', CS%diag%axesT1, Time, &
'random pattern for KE dissipation', 'None')

if (is_root_pe()) &
write(*,'(/12x,a/)') '=== COMPLETED MOM STOCHASTIC INITIALIZATION ====='
if (CS%do_sppt .OR. CS%pert_epbl) &
call MOM_mesg(' === COMPLETED MOM STOCHASTIC INITIALIZATION =====')

call callTree_leave("ocean_model_init(")
return

end subroutine stochastics_init

!> update_ocean_model uses the forcing in Ice_ocean_boundary to advance the
Expand Down