SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_atmos_dyn_tinteg_tracer Module Reference

module Atmosphere / Dynamics Temporal integration More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_tracer_setup (ATMOS_DYN_Tinteg_tracer_TYPE)
 Register. More...
 

Variables

procedure(tinteg), pointer, public atmos_dyn_tinteg_tracer => NULL()
 

Detailed Description

module Atmosphere / Dynamics Temporal integration

Description
Temporal integration scheme selecter for dynamical tracer advection
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_tracer_setup()

subroutine, public scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup ( character(len=*), intent(in)  ATMOS_DYN_Tinteg_tracer_TYPE)

Register.

Definition at line 91 of file scale_atmos_dyn_tinteg_tracer.F90.

References atmos_dyn_tinteg_tracer, scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler(), scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler_setup(), scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3(), scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_setup(), and scale_prc::prc_abort().

Referenced by scale_atmos_dyn::atmos_dyn_setup().

91 
92  use scale_precision
94  use scale_index
95  use scale_prc, only: &
96  prc_abort
103  implicit none
104  character(len=*), intent(in) :: atmos_dyn_tinteg_tracer_type
105  !---------------------------------------------------------------------------
106 
107  select case( atmos_dyn_tinteg_tracer_type )
108  case( 'EULER' )
110  atmos_dyn_tinteg_tracer_type )
111  atmos_dyn_tinteg_tracer => atmos_dyn_tinteg_tracer_euler
112  case( 'RK3WS2002' )
114  atmos_dyn_tinteg_tracer_type )
115  atmos_dyn_tinteg_tracer => atmos_dyn_tinteg_tracer_rk3
116  case( 'OFF', 'NONE' )
117  ! do nothing
118  case default
119  log_error("ATMOS_DYN_Tinteg_tracer_setup",*) 'ATMOS_DYN_TINTEG_TRACER_TYPE is invalid: ', atmos_dyn_tinteg_tracer_type
120  call prc_abort
121  end select
122 
123  return
subroutine, public atmos_dyn_tinteg_tracer_euler(QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, BND_W, BND_E, BND_S, BND_N, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
EULER.
subroutine, public atmos_dyn_tinteg_tracer_euler_setup(tinteg_type)
Setup.
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module PRECISION
subroutine, public atmos_dyn_tinteg_tracer_rk3_setup(tinteg_type)
Setup.
subroutine, public atmos_dyn_tinteg_tracer_rk3(QTRC, QTRC0, RHOQ_t, DENS0, DENS, mflx_hi, num_diff, GSQRT, MAPF, CDZ, RCDZ, RCDX, RCDY, BND_W, BND_E, BND_S, BND_N, dtl, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM)
RK3.
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ atmos_dyn_tinteg_tracer

procedure(tinteg), pointer, public scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer => NULL()

Definition at line 70 of file scale_atmos_dyn_tinteg_tracer.F90.

Referenced by atmos_dyn_tinteg_tracer_setup(), and scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve().

70  procedure(tinteg), pointer :: atmos_dyn_tinteg_tracer => null()