SCALE-RM
scale_atmos_dyn_tinteg_short.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
20  use scale_index
21  use scale_tracer
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
30 
32  abstract interface
33  subroutine short( &
34  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, & ! (inout)
35  mflx_hi, tflx_hi, & ! (inout)
36  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! (in)
37  dpres0, cvtot, corioli, & ! (in)
38  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! (in)
39  flag_fct_momentum, flag_fct_t, & ! (in)
40  flag_fct_along_stream, & ! (in)
41  cdz, fdz, fdx, fdy, & ! (in)
42  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
43  phi, gsqrt, j13g, j23g, j33g, mapf, & ! (in)
44  ref_pres, ref_dens, & ! (in)
45  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
46  dt ) ! (in)
47  use scale_precision
49  use scale_index
50  real(RP), intent(inout) :: dens(ka,ia,ja)
51  real(RP), intent(inout) :: momz(ka,ia,ja)
52  real(RP), intent(inout) :: momx(ka,ia,ja)
53  real(RP), intent(inout) :: momy(ka,ia,ja)
54  real(RP), intent(inout) :: rhot(ka,ia,ja)
55  real(RP), intent(inout) :: prog(ka,ia,ja,va)
56 
57  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
58  real(RP), intent(inout) :: tflx_hi(ka,ia,ja,3)
59 
60  real(RP), intent(in) :: dens_t(ka,ia,ja)
61  real(RP), intent(in) :: momz_t(ka,ia,ja)
62  real(RP), intent(in) :: momx_t(ka,ia,ja)
63  real(RP), intent(in) :: momy_t(ka,ia,ja)
64  real(RP), intent(in) :: rhot_t(ka,ia,ja)
65 
66  real(RP), intent(in) :: dpres0(ka,ia,ja)
67  real(RP), intent(in) :: cvtot(ka,ia,ja)
68  real(RP), intent(in) :: corioli(ia,ja)
69 
70  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
71  real(RP), intent(in) :: wdamp_coef(ka)
72  real(RP), intent(in) :: divdmp_coef
73  real(RP), intent(in) :: ddiv(ka,ia,ja)
74 
75  logical, intent(in) :: flag_fct_momentum
76  logical, intent(in) :: flag_fct_t
77  logical, intent(in) :: flag_fct_along_stream
78 
79  real(RP), intent(in) :: cdz (ka)
80  real(RP), intent(in) :: fdz (ka-1)
81  real(RP), intent(in) :: fdx (ia-1)
82  real(RP), intent(in) :: fdy (ja-1)
83  real(RP), intent(in) :: rcdz(ka)
84  real(RP), intent(in) :: rcdx(ia)
85  real(RP), intent(in) :: rcdy(ja)
86  real(RP), intent(in) :: rfdz(ka-1)
87  real(RP), intent(in) :: rfdx(ia-1)
88  real(RP), intent(in) :: rfdy(ja-1)
89 
90  real(RP), intent(in) :: phi (ka,ia,ja)
91  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
92  real(RP), intent(in) :: j13g (ka,ia,ja,7)
93  real(RP), intent(in) :: j23g (ka,ia,ja,7)
94  real(RP), intent(in) :: j33g
95  real(RP), intent(in) :: mapf (ia,ja,2,4)
96 
97  real(RP), intent(in) :: ref_pres(ka,ia,ja)
98  real(RP), intent(in) :: ref_dens(ka,ia,ja)
99 
100  logical, intent(in) :: bnd_w
101  logical, intent(in) :: bnd_e
102  logical, intent(in) :: bnd_s
103  logical, intent(in) :: bnd_n
104 
105  real(RP), intent(in) :: dt
106  end subroutine short
107  end interface
108  procedure(short), pointer :: atmos_dyn_tinteg_short => null()
109  public :: atmos_dyn_tinteg_short
110 
111  !-----------------------------------------------------------------------------
112  !
113  !++ Public parameters & variables
114  !
115  !-----------------------------------------------------------------------------
116  !
117  !++ Private procedure
118  !
119  !-----------------------------------------------------------------------------
120  !
121  !++ Private parameters & variables
122  !
123  !-----------------------------------------------------------------------------
124 contains
125  !-----------------------------------------------------------------------------
127  subroutine atmos_dyn_tinteg_short_setup( &
128  ATMOS_DYN_Tinteg_short_TYPE )
130  use scale_precision
132  use scale_index
133  use scale_prc, only: &
134  prc_abort
141  implicit none
142 
143  character(len=*), intent(in) :: ATMOS_DYN_Tinteg_short_TYPE
144  !---------------------------------------------------------------------------
145 
146  select case( atmos_dyn_tinteg_short_type )
147  case( 'RK3', 'RK3WS2002' )
149  atmos_dyn_tinteg_short_type )
151  case( 'RK4' )
153  atmos_dyn_tinteg_short_type )
155  case( 'OFF', 'NONE' )
156  ! do nothing
157  case default
158  log_error("ATMOS_DYN_Tinteg_short_setup",*) 'ATMOS_DYN_TINTEG_SHORT_TYPE is invalid: ', atmos_dyn_tinteg_short_type
159  call prc_abort
160  end select
161 
162  return
163  end subroutine atmos_dyn_tinteg_short_setup
164 
integer, public va
Definition: scale_index.F90:35
integer, public ia
of whole cells: x, local, with HALO
subroutine, public atmos_dyn_tinteg_short_rk4(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, mflx_hi, tflx_hi, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, DPRES0, CVtot, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_pres, REF_dens, BND_W, BND_E, BND_S, BND_N, dt)
RK3.
subroutine, public atmos_dyn_tinteg_short_rk4_setup(tinteg_type)
Setup.
integer, public ja
of whole cells: y, local, with HALO
module Atmosphere / Dynamics Temporal integration
module TRACER
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
subroutine, public atmos_dyn_tinteg_short_rk3(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, mflx_hi, tflx_hi, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, DPRES0, CVtot, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_pres, REF_dens, BND_W, BND_E, BND_S, BND_N, dt)
RK3.
subroutine, public atmos_dyn_tinteg_short_rk3_setup(tinteg_type)
Setup.
module profiler
Definition: scale_prof.F90:11
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
subroutine, public atmos_dyn_tinteg_short_setup(ATMOS_DYN_Tinteg_short_TYPE)
Register.
procedure(short), pointer, public atmos_dyn_tinteg_short
module STDIO
Definition: scale_io.F90:10