SCALE-RM
scale_atmos_dyn_tinteg_large.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  subroutine large( &
34  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
35  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
36  mflx_hi, tflx_hi, &
37  num_diff, num_diff_q, &
38  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
39  CORIOLI, &
40  CDZ, CDX, CDY, FDZ, FDX, FDY, &
41  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
42  PHI, GSQRT, &
43  J13G, J23G, J33G, MAPF, &
44  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
45  REF_dens, REF_pott, REF_qv, REF_pres, &
46  BND_W, BND_E, BND_S, BND_N, &
47  ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, &
48  BND_QA, BND_SMOOTHER_FACT, &
49  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
50  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
51  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
52  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
53  wdamp_coef, &
54  divdmp_coef, &
55  FLAG_TRACER_SPLIT_TEND, &
56  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
57  FLAG_FCT_ALONG_STREAM, &
58  USE_AVERAGE, &
59  I_QV, &
60  DTL, DTS )
61  use scale_precision
63  use scale_index
64  use scale_tracer
65  real(RP), intent(inout) :: dens(ka,ia,ja)
66  real(RP), intent(inout) :: momz(ka,ia,ja)
67  real(RP), intent(inout) :: momx(ka,ia,ja)
68  real(RP), intent(inout) :: momy(ka,ia,ja)
69  real(RP), intent(inout) :: rhot(ka,ia,ja)
70  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
71  real(RP), intent(inout) :: prog(ka,ia,ja,va)
72 
73  real(RP), intent(inout) :: dens_av(ka,ia,ja)
74  real(RP), intent(inout) :: momz_av(ka,ia,ja)
75  real(RP), intent(inout) :: momx_av(ka,ia,ja)
76  real(RP), intent(inout) :: momy_av(ka,ia,ja)
77  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
78  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
79 
80  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
81  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
82  real(RP), intent(out) :: num_diff(ka,ia,ja,5,3)
83  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
84 
85  real(RP), intent(in) :: dens_tp(ka,ia,ja)
86  real(RP), intent(in) :: momz_tp(ka,ia,ja)
87  real(RP), intent(in) :: momx_tp(ka,ia,ja)
88  real(RP), intent(in) :: momy_tp(ka,ia,ja)
89  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
90  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
91 
92  real(RP), intent(in) :: corioli(ia,ja)
93 
94  real(RP), intent(in) :: cdz (ka)
95  real(RP), intent(in) :: cdx (ia)
96  real(RP), intent(in) :: cdy (ja)
97  real(RP), intent(in) :: fdz (ka-1)
98  real(RP), intent(in) :: fdx (ia-1)
99  real(RP), intent(in) :: fdy (ja-1)
100  real(RP), intent(in) :: rcdz(ka)
101  real(RP), intent(in) :: rcdx(ia)
102  real(RP), intent(in) :: rcdy(ja)
103  real(RP), intent(in) :: rfdz(ka-1)
104  real(RP), intent(in) :: rfdx(ia-1)
105  real(RP), intent(in) :: rfdy(ja-1)
106 
107  real(RP), intent(in) :: phi (ka,ia,ja)
108  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
109  real(RP), intent(in) :: j13g (ka,ia,ja,7)
110  real(RP), intent(in) :: j23g (ka,ia,ja,7)
111  real(RP), intent(in) :: j33g
112  real(RP), intent(in) :: mapf (ia,ja,2,4)
113 
114  real(RP), intent(in) :: aq_r (qa)
115  real(RP), intent(in) :: aq_cv (qa)
116  real(RP), intent(in) :: aq_cp (qa)
117  real(RP), intent(in) :: aq_mass(qa)
118 
119  real(RP), intent(in) :: ref_dens(ka,ia,ja)
120  real(RP), intent(in) :: ref_pott(ka,ia,ja)
121  real(RP), intent(in) :: ref_qv (ka,ia,ja)
122  real(RP), intent(in) :: ref_pres(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) :: nd_coef
130  real(RP), intent(in) :: nd_coef_q
131  integer, intent(in) :: nd_order
132  real(RP), intent(in) :: nd_sfc_fact
133  logical, intent(in) :: nd_use_rs
134 
135  integer, intent(in) :: bnd_qa
136  real(RP), intent(in) :: bnd_smoother_fact
137 
138  real(RP), intent(in) :: damp_dens(ka,ia,ja)
139  real(RP), intent(in) :: damp_velz(ka,ia,ja)
140  real(RP), intent(in) :: damp_velx(ka,ia,ja)
141  real(RP), intent(in) :: damp_vely(ka,ia,ja)
142  real(RP), intent(in) :: damp_pott(ka,ia,ja)
143  real(RP), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
144  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
145  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
146  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
147  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
148  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
149  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
150 
151  real(RP), intent(in) :: wdamp_coef(ka)
152  real(RP), intent(in) :: divdmp_coef
153  logical, intent(in) :: flag_tracer_split_tend
154  logical, intent(in) :: flag_fct_momentum
155  logical, intent(in) :: flag_fct_t
156  logical, intent(in) :: flag_fct_tracer
157  logical, intent(in) :: flag_fct_along_stream
158  logical, intent(in) :: use_average
159 
160  integer, intent(in) :: i_qv
161 
162  real(DP), intent(in) :: dtl
163  real(DP), intent(in) :: dts
164  end subroutine large
165  end interface
166  procedure(large), pointer :: atmos_dyn_tinteg_large => null()
167  public :: atmos_dyn_tinteg_large
168 
169  !-----------------------------------------------------------------------------
170  !
171  !++ Public parameters & variables
172  !
173  !-----------------------------------------------------------------------------
174  !
175  !++ Private procedure
176  !
177  !-----------------------------------------------------------------------------
178  !
179  !++ Private parameters & variables
180  !
181  !-----------------------------------------------------------------------------
182 contains
183  !-----------------------------------------------------------------------------
185  subroutine atmos_dyn_tinteg_large_setup( &
186  ATMOS_DYN_Tinteg_large_TYPE )
189  use scale_index
190  use scale_prc, only: &
191  prc_abort
198  implicit none
199 
200  character(len=*), intent(in) :: ATMOS_DYN_Tinteg_large_TYPE
201  !---------------------------------------------------------------------------
202 
203  select case( atmos_dyn_tinteg_large_type )
204  case( 'EULER' )
206  atmos_dyn_tinteg_large_type )
208  case( 'RK3' )
210  atmos_dyn_tinteg_large_type )
212  case( 'OFF', 'NONE' )
213  ! do nothing
214  case default
215  log_error("ATMOS_DYN_Tinteg_large_setup",*) 'ATMOS_DYN_TINTEG_LARGE_TYPE is invalid: ', atmos_dyn_tinteg_large_type
216  call prc_abort
217  end select
218 
219  return
220  end subroutine atmos_dyn_tinteg_large_setup
221 
integer, public va
Definition: scale_index.F90:35
subroutine, public atmos_dyn_tinteg_large_rk3(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, mflx_hi, tflx_hi, num_diff, num_diff_q, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_R, AQ_CV, AQ_CP, AQ_MASS, REF_dens, REF_pott, REF_qv, REF_pres, BND_W, BND_E, BND_S, BND_N, ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_SMOOTHER_FACT, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, wdamp_coef, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTL, DTS)
RK3.
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
integer, public ia
of whole cells: x, local, with HALO
module Atmosphere / Dynamics Temporal integration
integer, public qa
integer, public ja
of whole cells: y, local, with HALO
procedure(large), pointer, public atmos_dyn_tinteg_large
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
subroutine, public atmos_dyn_tinteg_large_rk3_setup(tinteg_type)
Setup.
module PROCESS
Definition: scale_prc.F90:11
subroutine, public atmos_dyn_tinteg_large_euler_setup(tinteg_type)
Setup.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module profiler
Definition: scale_prof.F90:11
subroutine, public atmos_dyn_tinteg_large_euler(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, mflx_hi, tflx_hi, num_diff, num_diff_q, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_R, AQ_CV, AQ_CP, AQ_MASS, REF_dens, REF_pott, REF_qv, REF_pres, BND_W, BND_E, BND_S, BND_N, ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_SMOOTHER_FACT, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, wdamp_coef, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTL, DTS)
RK3.
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
module STDIO
Definition: scale_io.F90:10