SCALE-RM
scale_atmos_dyn_tstep_short.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
14 #include "inc_openmp.h"
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  use scale_index
25  use scale_tracer
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
34 
35  abstract interface
36 
37  subroutine short_setup
38  end subroutine short_setup
39 
41  subroutine short( DENS_new, MOMZ_new, MOMX_new, MOMY_new, RHOT_new, & ! (out)
42  prog_new, & ! (out)
43  mflx_hi, tflx_hi, & ! (out)
44  dens0, momz0, momx0, momy0, rhot0, & ! (in)
45  dens, momz, momx, momy, rhot, & ! (in)
46  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! (in)
47  prog0, prog, & ! (in)
48  dpres0, rt2p, corioli, & ! (in)
49  num_diff, divdmp_coef, ddiv, & ! (in)
50  flag_fct_momentum, flag_fct_t, & ! (in)
51  flag_fct_along_stream, & ! (in)
52  cdz, fdz, fdx, fdy, & ! (in)
53  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
54  phi, gsqrt, j13g, j23g, j33g, mapf, & ! (in)
55  ref_dens, ref_rhot, & ! (in)
56  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
57  dtrk, dt ) ! (in)
58  use scale_precision
60  use scale_index
61  real(RP), intent(out) :: dens_new(ka,ia,ja) ! prognostic variables
62  real(RP), intent(out) :: momz_new(ka,ia,ja) !
63  real(RP), intent(out) :: momx_new(ka,ia,ja) !
64  real(RP), intent(out) :: momy_new(ka,ia,ja) !
65  real(RP), intent(out) :: rhot_new(ka,ia,ja) !
66  real(RP), intent(out) :: prog_new(ka,ia,ja,va) !
67 
68  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3) ! mass flux
69  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3) ! internal energy flux
70 
71  real(RP), intent(in),target :: dens0(ka,ia,ja) ! prognostic variables at previous dynamical time step
72  real(RP), intent(in),target :: momz0(ka,ia,ja) !
73  real(RP), intent(in),target :: momx0(ka,ia,ja) !
74  real(RP), intent(in),target :: momy0(ka,ia,ja) !
75  real(RP), intent(in),target :: rhot0(ka,ia,ja) !
76 
77  real(RP), intent(in) :: dens(ka,ia,ja) ! prognostic variables at previous RK step
78  real(RP), intent(in) :: momz(ka,ia,ja) !
79  real(RP), intent(in) :: momx(ka,ia,ja) !
80  real(RP), intent(in) :: momy(ka,ia,ja) !
81  real(RP), intent(in) :: rhot(ka,ia,ja) !
82 
83  real(RP), intent(in) :: dens_t(ka,ia,ja) ! tendency
84  real(RP), intent(in) :: momz_t(ka,ia,ja) !
85  real(RP), intent(in) :: momx_t(ka,ia,ja) !
86  real(RP), intent(in) :: momy_t(ka,ia,ja) !
87  real(RP), intent(in) :: rhot_t(ka,ia,ja) !
88 
89  real(RP), intent(in) :: prog0(ka,ia,ja,va)
90  real(RP), intent(in) :: prog (ka,ia,ja,va)
91 
92  real(RP), intent(in) :: dpres0 (ka,ia,ja)
93  real(RP), intent(in) :: rt2p (ka,ia,ja)
94  real(RP), intent(in) :: corioli (1, ia,ja)
95  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
96  real(RP), intent(in) :: divdmp_coef
97  real(RP), intent(in) :: ddiv(ka,ia,ja)
98 
99  logical, intent(in) :: flag_fct_momentum
100  logical, intent(in) :: flag_fct_t
101  logical, intent(in) :: flag_fct_along_stream
102 
103  real(RP), intent(in) :: cdz (ka)
104  real(RP), intent(in) :: fdz (ka-1)
105  real(RP), intent(in) :: fdx (ia-1)
106  real(RP), intent(in) :: fdy (ja-1)
107  real(RP), intent(in) :: rcdz(ka)
108  real(RP), intent(in) :: rcdx(ia)
109  real(RP), intent(in) :: rcdy(ja)
110  real(RP), intent(in) :: rfdz(ka-1)
111  real(RP), intent(in) :: rfdx(ia-1)
112  real(RP), intent(in) :: rfdy(ja-1)
113 
114  real(RP), intent(in) :: phi (ka,ia,ja)
115  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
116  real(RP), intent(in) :: j13g (ka,ia,ja,7)
117  real(RP), intent(in) :: j23g (ka,ia,ja,7)
118  real(RP), intent(in) :: j33g
119  real(RP), intent(in) :: mapf (ia,ja,2,4)
120  real(RP), intent(in) :: ref_dens(ka,ia,ja)
121  real(RP), intent(in) :: ref_rhot(ka,ia,ja)
122 
123  logical, intent(in) :: bnd_w
124  logical, intent(in) :: bnd_e
125  logical, intent(in) :: bnd_s
126  logical, intent(in) :: bnd_n
127 
128  real(RP), intent(in) :: dtrk
129  real(RP), intent(in) :: dt
130  end subroutine short
131 
132  end interface
133 
134  procedure(short_setup), pointer :: atmos_dyn_tstep_short_setup => null()
136  procedure(short), pointer :: atmos_dyn_tstep_short => null()
137  public :: atmos_dyn_tstep_short
138 
139  !-----------------------------------------------------------------------------
140  !
141  !++ Public parameters & variables
142  !
143  !-----------------------------------------------------------------------------
144  !
145  !++ Private procedure
146  !
147  !-----------------------------------------------------------------------------
148  !
149  !++ Private parameters & variables
150  !
151  !-----------------------------------------------------------------------------
152 contains
153  !-----------------------------------------------------------------------------
155  subroutine atmos_dyn_tstep_short_regist( &
156  ATMOS_DYN_TYPE, &
157  VA_out, &
158  VAR_NAME, VAR_DESC, VAR_UNIT )
160  use scale_grid_index
161  use scale_index
162  use scale_process, only: &
164 #define EXTM(pre, name, post) pre ## name ## post
165 #define NAME(pre, name, post) EXTM(pre, name, post)
166 #ifdef DYNAMICS
167  use name(scale_atmos_dyn_tstep_, dynamics,), only: &
168  name(atmos_dyn_rk_tstep_short_, dynamics, _regist), &
169  name(atmos_dyn_rk_tstep_short_, dynamics, _setup), &
170  name(atmos_dyn_rk_tstep_short_, dynamics,)
171 #else
184 #endif
185  implicit none
186 
187  character(len=*), intent(in) :: ATMOS_DYN_TYPE
188  integer, intent(out) :: VA_out
189  character(len=H_SHORT), intent(out) :: VAR_NAME(:)
190  character(len=H_MID), intent(out) :: VAR_DESC(:)
191  character(len=H_SHORT), intent(out) :: VAR_UNIT(:)
192  !---------------------------------------------------------------------------
193 
194 #ifdef DYNAMICS
195  name(atmos_dyn_tstep_, dynamics, _regist)( &
196  atmos_dyn_type )
197  atmos_dyn_tstep_short_ => name(atmos_dyn_tstep_short_, dynamics,)
198  atmos_dyn_tstep_short_setup => name(atmos_dyn_tstep_short_, dynamics, _setup)
199 #else
200  select case ( atmos_dyn_type )
201  case ( 'FVM-HEVE', 'HEVE' )
202 
203  call atmos_dyn_tstep_short_fvm_heve_regist( atmos_dyn_type, & ! [IN]
204  va_out, & ! [OUT]
205  var_name, var_desc, var_unit ) ! [OUT]
206 
209 
210  case ( 'FVM-HEVI', 'HEVI' )
211 
212  call atmos_dyn_tstep_short_fvm_hevi_regist( atmos_dyn_type, & ! [IN]
213  va_out, & ! [OUT]
214  var_name, var_desc, var_unit ) ! [OUT]
215 
218 
219  case ( 'FVM-HIVI', 'HIVI' )
220 
221  if( io_l ) write(io_fid_log,*) 'xxx HIVI is temtatively disabled'
222  call prc_mpistop
223 
224  call atmos_dyn_tstep_short_fvm_hivi_regist( atmos_dyn_type, & ! [IN]
225  va_out, & ! [OUT]
226  var_name, var_desc, var_unit ) ! [OUT]
227 
230 
231  case ( 'OFF', 'NONE' )
232 
233  va_out = 0
234  var_name(:) = ""
235  var_desc(:) = ""
236  var_unit(:) = ""
237 
238  case default
239  write(*,*) 'xxx ATMOS_DYN_TYPE is invalid: ', atmos_dyn_type
240  call prc_mpistop
241  end select
242 #endif
243 
244  return
245  end subroutine atmos_dyn_tstep_short_regist
246 
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
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
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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, 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, dt)
module STDIO
Definition: scale_stdio.F90:12
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, 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, dt)
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, 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, dt)
subroutine, public atmos_dyn_tstep_short_fvm_hivi_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
procedure(short), pointer, public atmos_dyn_tstep_short
module PROCESS
subroutine, public atmos_dyn_tstep_short_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
module profiler
Definition: scale_prof.F90:10
subroutine, public atmos_dyn_tstep_short_fvm_heve_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
module PRECISION
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup
integer, public ja
of y whole cells (local, with HALO)