SCALE-RM
scale_atmos_dyn_tstep_short.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_index
22  use scale_tracer
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
31 
32  abstract interface
33 
34  subroutine short_setup
35  end subroutine short_setup
36 
38  subroutine short( DENS_new, MOMZ_new, MOMX_new, MOMY_new, RHOT_new, & ! (out)
39  prog_new, & ! (out)
40  mflx_hi, tflx_hi, & ! (out)
41  dens0, momz0, momx0, momy0, rhot0, & ! (in)
42  dens, momz, momx, momy, rhot, & ! (in)
43  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! (in)
44  prog0, prog, & ! (in)
45  dpres0, rt2p, corioli, & ! (in)
46  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! (in)
47  flag_fct_momentum, flag_fct_t, & ! (in)
48  flag_fct_along_stream, & ! (in)
49  cdz, fdz, fdx, fdy, & ! (in)
50  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
51  phi, gsqrt, j13g, j23g, j33g, mapf, & ! (in)
52  ref_dens, ref_rhot, & ! (in)
53  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
54  dtrk, last ) ! (in)
55  use scale_precision
57  use scale_index
58  real(RP), intent(out) :: dens_new(ka,ia,ja) ! prognostic variables
59  real(RP), intent(out) :: momz_new(ka,ia,ja) !
60  real(RP), intent(out) :: momx_new(ka,ia,ja) !
61  real(RP), intent(out) :: momy_new(ka,ia,ja) !
62  real(RP), intent(out) :: rhot_new(ka,ia,ja) !
63  real(RP), intent(out) :: prog_new(ka,ia,ja,va) !
64 
65  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3) ! mass flux
66  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3) ! internal energy flux
67 
68  real(RP), intent(in),target :: dens0(ka,ia,ja) ! prognostic variables at previous dynamical time step
69  real(RP), intent(in),target :: momz0(ka,ia,ja) !
70  real(RP), intent(in),target :: momx0(ka,ia,ja) !
71  real(RP), intent(in),target :: momy0(ka,ia,ja) !
72  real(RP), intent(in),target :: rhot0(ka,ia,ja) !
73 
74  real(RP), intent(in) :: dens(ka,ia,ja) ! prognostic variables at previous RK step
75  real(RP), intent(in) :: momz(ka,ia,ja) !
76  real(RP), intent(in) :: momx(ka,ia,ja) !
77  real(RP), intent(in) :: momy(ka,ia,ja) !
78  real(RP), intent(in) :: rhot(ka,ia,ja) !
79 
80  real(RP), intent(in) :: dens_t(ka,ia,ja) ! tendency
81  real(RP), intent(in) :: momz_t(ka,ia,ja) !
82  real(RP), intent(in) :: momx_t(ka,ia,ja) !
83  real(RP), intent(in) :: momy_t(ka,ia,ja) !
84  real(RP), intent(in) :: rhot_t(ka,ia,ja) !
85 
86  real(RP), intent(in) :: prog0(ka,ia,ja,va)
87  real(RP), intent(in) :: prog (ka,ia,ja,va)
88 
89  real(RP), intent(in) :: dpres0 (ka,ia,ja)
90  real(RP), intent(in) :: rt2p (ka,ia,ja)
91  real(RP), intent(in) :: corioli (1, ia,ja)
92  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
93  real(RP), intent(in) :: wdamp_coef(ka)
94  real(RP), intent(in) :: divdmp_coef
95  real(RP), intent(in) :: ddiv(ka,ia,ja)
96 
97  logical, intent(in) :: flag_fct_momentum
98  logical, intent(in) :: flag_fct_t
99  logical, intent(in) :: flag_fct_along_stream
100 
101  real(RP), intent(in) :: cdz (ka)
102  real(RP), intent(in) :: fdz (ka-1)
103  real(RP), intent(in) :: fdx (ia-1)
104  real(RP), intent(in) :: fdy (ja-1)
105  real(RP), intent(in) :: rcdz(ka)
106  real(RP), intent(in) :: rcdx(ia)
107  real(RP), intent(in) :: rcdy(ja)
108  real(RP), intent(in) :: rfdz(ka-1)
109  real(RP), intent(in) :: rfdx(ia-1)
110  real(RP), intent(in) :: rfdy(ja-1)
111 
112  real(RP), intent(in) :: phi (ka,ia,ja)
113  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
114  real(RP), intent(in) :: j13g (ka,ia,ja,7)
115  real(RP), intent(in) :: j23g (ka,ia,ja,7)
116  real(RP), intent(in) :: j33g
117  real(RP), intent(in) :: mapf (ia,ja,2,4)
118  real(RP), intent(in) :: ref_dens(ka,ia,ja)
119  real(RP), intent(in) :: ref_rhot(ka,ia,ja)
120 
121  logical, intent(in) :: bnd_w
122  logical, intent(in) :: bnd_e
123  logical, intent(in) :: bnd_s
124  logical, intent(in) :: bnd_n
125 
126  real(RP), intent(in) :: dtrk
127  logical, intent(in) :: last
128  end subroutine short
129 
130  end interface
131 
132  procedure(short_setup), pointer :: atmos_dyn_tstep_short_setup => null()
134  procedure(short), pointer :: atmos_dyn_tstep_short => null()
135  public :: atmos_dyn_tstep_short
136 
137  !-----------------------------------------------------------------------------
138  !
139  !++ Public parameters & variables
140  !
141  !-----------------------------------------------------------------------------
142  !
143  !++ Private procedure
144  !
145  !-----------------------------------------------------------------------------
146  !
147  !++ Private parameters & variables
148  !
149  !-----------------------------------------------------------------------------
150 contains
151  !-----------------------------------------------------------------------------
153  subroutine atmos_dyn_tstep_short_regist( &
154  ATMOS_DYN_TYPE, &
155  VA_out, &
156  VAR_NAME, VAR_DESC, VAR_UNIT )
159  use scale_index
160  use scale_prc, only: &
161  prc_abort
174  implicit none
175 
176  character(len=*), intent(in) :: ATMOS_DYN_TYPE
177  integer, intent(out) :: VA_out
178  character(len=H_SHORT), intent(out) :: VAR_NAME(:)
179  character(len=H_MID), intent(out) :: VAR_DESC(:)
180  character(len=H_SHORT), intent(out) :: VAR_UNIT(:)
181  !---------------------------------------------------------------------------
182 
183  select case( atmos_dyn_type )
184  case( 'FVM-HEVE', 'HEVE' )
185 
186  call atmos_dyn_tstep_short_fvm_heve_regist( atmos_dyn_type, & ! [IN]
187  va_out, & ! [OUT]
188  var_name, var_desc, var_unit ) ! [OUT]
189 
192 
193  case( 'FVM-HEVI', 'HEVI' )
194 
195  call atmos_dyn_tstep_short_fvm_hevi_regist( atmos_dyn_type, & ! [IN]
196  va_out, & ! [OUT]
197  var_name, var_desc, var_unit ) ! [OUT]
198 
201 
202  case( 'FVM-HIVI', 'HIVI' )
203 
204  log_error("ATMOS_DYN_Tstep_short_regist",*) 'HIVI is tentatively disabled'
205  call prc_abort
206 
207  call atmos_dyn_tstep_short_fvm_hivi_regist( atmos_dyn_type, & ! [IN]
208  va_out, & ! [OUT]
209  var_name, var_desc, var_unit ) ! [OUT]
210 
213 
214  case( 'OFF', 'NONE' )
215 
216  va_out = 0
217  var_name(:) = ""
218  var_desc(:) = ""
219  var_unit(:) = ""
220 
221  case default
222  log_error("ATMOS_DYN_Tstep_short_regist",*) 'ATMOS_DYN_TYPE is invalid: ', atmos_dyn_type
223  call prc_abort
224  end select
225 
226  return
227  end subroutine atmos_dyn_tstep_short_regist
228 
integer, public va
Definition: scale_index.F90:35
subroutine, public atmos_dyn_tstep_short_fvm_hevi_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
module Atmosphere / Dynamical scheme
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
subroutine, public atmos_dyn_tstep_short_fvm_heve(DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, PROG_RK, mflx_hi, tflx_hi, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, DENS, MOMZ, MOMX, MOMY, RHOT, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, PROG0, PROG, DPRES0, RT2P, 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_dens, REF_rhot, BND_W, BND_E, BND_S, BND_N, dtrk, last)
subroutine, public atmos_dyn_tstep_short_fvm_hivi_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
subroutine, public atmos_dyn_tstep_short_fvm_hivi(DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, PROG_RK, mflx_hi, tflx_hi, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, DENS, MOMZ, MOMX, MOMY, RHOT, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, PROG0, PROG, DPRES0, RT2P, 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_dens, REF_rhot, BND_W, BND_E, BND_S, BND_N, dtrk, last)
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
module PROCESS
Definition: scale_prc.F90:11
procedure(short), pointer, public atmos_dyn_tstep_short
subroutine, public atmos_dyn_tstep_short_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module profiler
Definition: scale_prof.F90:11
subroutine, public atmos_dyn_tstep_short_fvm_heve_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
subroutine, public atmos_dyn_tstep_short_fvm_hevi(DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, PROG_RK, mflx_hi, tflx_hi, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, DENS, MOMZ, MOMX, MOMY, RHOT, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, PROG0, PROG, DPRES0, RT2P, 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_dens, REF_rhot, BND_W, BND_E, BND_S, BND_N, dtrk, last)
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
module STDIO
Definition: scale_io.F90:10
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup