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, wdamp_coef, 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, last ) ! (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) :: wdamp_coef(KA)
97  real(RP), intent(in) :: divdmp_coef
98  real(RP), intent(in) :: DDIV(KA,IA,JA)
99 
100  logical, intent(in) :: FLAG_FCT_MOMENTUM
101  logical, intent(in) :: FLAG_FCT_T
102  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
103 
104  real(RP), intent(in) :: CDZ (KA)
105  real(RP), intent(in) :: FDZ (KA-1)
106  real(RP), intent(in) :: FDX (IA-1)
107  real(RP), intent(in) :: FDY (JA-1)
108  real(RP), intent(in) :: RCDZ(KA)
109  real(RP), intent(in) :: RCDX(IA)
110  real(RP), intent(in) :: RCDY(JA)
111  real(RP), intent(in) :: RFDZ(KA-1)
112  real(RP), intent(in) :: RFDX(IA-1)
113  real(RP), intent(in) :: RFDY(JA-1)
114 
115  real(RP), intent(in) :: PHI (KA,IA,JA)
116  real(RP), intent(in) :: GSQRT (KA,IA,JA,7)
117  real(RP), intent(in) :: J13G (KA,IA,JA,7)
118  real(RP), intent(in) :: J23G (KA,IA,JA,7)
119  real(RP), intent(in) :: J33G
120  real(RP), intent(in) :: MAPF (IA,JA,2,4)
121  real(RP), intent(in) :: REF_dens(KA,IA,JA)
122  real(RP), intent(in) :: REF_rhot(KA,IA,JA)
123 
124  logical, intent(in) :: BND_W
125  logical, intent(in) :: BND_E
126  logical, intent(in) :: BND_S
127  logical, intent(in) :: BND_N
128 
129  real(RP), intent(in) :: dtrk
130  logical, intent(in) :: last
131  end subroutine short
132 
133  end interface
134 
135  procedure(short_setup), pointer :: atmos_dyn_tstep_short_setup => null()
137  procedure(short), pointer :: atmos_dyn_tstep_short => null()
138  public :: atmos_dyn_tstep_short
139 
140  !-----------------------------------------------------------------------------
141  !
142  !++ Public parameters & variables
143  !
144  !-----------------------------------------------------------------------------
145  !
146  !++ Private procedure
147  !
148  !-----------------------------------------------------------------------------
149  !
150  !++ Private parameters & variables
151  !
152  !-----------------------------------------------------------------------------
153 contains
154  !-----------------------------------------------------------------------------
156  subroutine atmos_dyn_tstep_short_regist( &
157  ATMOS_DYN_TYPE, &
158  VA_out, &
159  VAR_NAME, VAR_DESC, VAR_UNIT )
161  use scale_grid_index
162  use scale_index
163  use scale_process, only: &
177  implicit none
178 
179  character(len=*), intent(in) :: atmos_dyn_type
180  integer, intent(out) :: va_out
181  character(len=H_SHORT), intent(out) :: var_name(:)
182  character(len=H_MID), intent(out) :: var_desc(:)
183  character(len=H_SHORT), intent(out) :: var_unit(:)
184  !---------------------------------------------------------------------------
185 
186  select case( atmos_dyn_type )
187  case( 'FVM-HEVE', 'HEVE' )
188 
189  call atmos_dyn_tstep_short_fvm_heve_regist( atmos_dyn_type, & ! [IN]
190  va_out, & ! [OUT]
191  var_name, var_desc, var_unit ) ! [OUT]
192 
195 
196  case( 'FVM-HEVI', 'HEVI' )
197 
198  call atmos_dyn_tstep_short_fvm_hevi_regist( atmos_dyn_type, & ! [IN]
199  va_out, & ! [OUT]
200  var_name, var_desc, var_unit ) ! [OUT]
201 
204 
205  case( 'FVM-HIVI', 'HIVI' )
206 
207  write(*,*) 'xxx HIVI is tentatively disabled'
208  call prc_mpistop
209 
210  call atmos_dyn_tstep_short_fvm_hivi_regist( atmos_dyn_type, & ! [IN]
211  va_out, & ! [OUT]
212  var_name, var_desc, var_unit ) ! [OUT]
213 
216 
217  case( 'OFF', 'NONE' )
218 
219  va_out = 0
220  var_name(:) = ""
221  var_desc(:) = ""
222  var_unit(:) = ""
223 
224  case default
225  write(*,*) 'xxx ATMOS_DYN_TYPE is invalid: ', atmos_dyn_type
226  call prc_mpistop
227  end select
228 
229  return
230  end subroutine atmos_dyn_tstep_short_regist
231 
subroutine, public prc_mpistop
Abort MPI.
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
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, 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 grid index
module TRACER
module Index
Definition: scale_index.F90:14
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.
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
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup