SCALE-RM
scale_atmos_dyn_tstep_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 
34  subroutine large( &
35  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
36  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
37  mflx_hi, tflx_hi, &
38  num_diff, num_diff_q, &
39  QTRC0, &
40  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
41  CORIOLI, &
42  CDZ, CDX, CDY, FDZ, FDX, FDY, &
43  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
44  PHI, GSQRT, &
45  J13G, J23G, J33G, MAPF, &
46  AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
47  REF_dens, REF_pott, REF_qv, REF_pres, &
48  BND_W, BND_E, BND_S, BND_N, &
49  ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, &
50  BND_QA, BND_SMOOTHER_FACT, &
51  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
52  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
53  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
54  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
55  wdamp_coef, &
56  divdmp_coef, &
57  FLAG_TRACER_SPLIT_TEND, &
58  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
59  FLAG_FCT_ALONG_STREAM, &
60  USE_AVERAGE, &
61  I_QV, &
62  DTL, DTS, Llast )
63  use scale_precision
65  use scale_index
66  use scale_tracer
67  real(RP), intent(inout) :: dens(ka,ia,ja)
68  real(RP), intent(inout) :: momz(ka,ia,ja)
69  real(RP), intent(inout) :: momx(ka,ia,ja)
70  real(RP), intent(inout) :: momy(ka,ia,ja)
71  real(RP), intent(inout) :: rhot(ka,ia,ja)
72  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
73  real(RP), intent(inout) :: prog(ka,ia,ja,va)
74 
75  real(RP), intent(inout) :: dens_av(ka,ia,ja)
76  real(RP), intent(inout) :: momz_av(ka,ia,ja)
77  real(RP), intent(inout) :: momx_av(ka,ia,ja)
78  real(RP), intent(inout) :: momy_av(ka,ia,ja)
79  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
80  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
81 
82  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
83  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
84 
85  real(RP), intent(out) :: num_diff(ka,ia,ja,5,3)
86  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
87 
88  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa)
89 
90  real(RP), intent(in) :: dens_tp(ka,ia,ja)
91  real(RP), intent(in) :: momz_tp(ka,ia,ja)
92  real(RP), intent(in) :: momx_tp(ka,ia,ja)
93  real(RP), intent(in) :: momy_tp(ka,ia,ja)
94  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
95  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
96 
97  real(RP), intent(in) :: corioli(ia,ja)
98 
99  real(RP), intent(in) :: cdz (ka)
100  real(RP), intent(in) :: cdx (ia)
101  real(RP), intent(in) :: cdy (ja)
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 
119  real(RP), intent(in) :: aq_r (qa)
120  real(RP), intent(in) :: aq_cv (qa)
121  real(RP), intent(in) :: aq_cp (qa)
122  real(RP), intent(in) :: aq_mass(qa)
123 
124  real(RP), intent(in) :: ref_dens(ka,ia,ja)
125  real(RP), intent(in) :: ref_pott(ka,ia,ja)
126  real(RP), intent(in) :: ref_qv (ka,ia,ja)
127  real(RP), intent(in) :: ref_pres(ka,ia,ja)
128 
129  logical, intent(in) :: bnd_w
130  logical, intent(in) :: bnd_e
131  logical, intent(in) :: bnd_s
132  logical, intent(in) :: bnd_n
133 
134  real(RP), intent(in) :: nd_coef
135  real(RP), intent(in) :: nd_coef_q
136  integer, intent(in) :: nd_order
137  real(RP), intent(in) :: nd_sfc_fact
138  logical, intent(in) :: nd_use_rs
139 
140  integer, intent(in) :: bnd_qa
141  real(RP), intent(in) :: bnd_smoother_fact
142 
143  real(RP), intent(in) :: damp_dens(ka,ia,ja)
144  real(RP), intent(in) :: damp_velz(ka,ia,ja)
145  real(RP), intent(in) :: damp_velx(ka,ia,ja)
146  real(RP), intent(in) :: damp_vely(ka,ia,ja)
147  real(RP), intent(in) :: damp_pott(ka,ia,ja)
148  real(RP), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
149 
150  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
151  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
152  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
153  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
154  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
155  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
156 
157  real(RP), intent(in) :: wdamp_coef(ka)
158  real(RP), intent(in) :: divdmp_coef
159 
160  logical, intent(in) :: flag_tracer_split_tend
161  logical, intent(in) :: flag_fct_momentum
162  logical, intent(in) :: flag_fct_t
163  logical, intent(in) :: flag_fct_tracer
164  logical, intent(in) :: flag_fct_along_stream
165 
166  logical, intent(in) :: use_average
167 
168  integer, intent(in) :: i_qv
169 
170  real(DP), intent(in) :: dtl
171  real(DP), intent(in) :: dts
172 
173  logical , intent(in) :: llast
174  end subroutine large
175 
176  end interface
177 
178  procedure(large), pointer :: atmos_dyn_tstep_large => null()
179  public :: atmos_dyn_tstep_large
180 
181  !-----------------------------------------------------------------------------
182  !
183  !++ Public parameters & variables
184  !
185  !-----------------------------------------------------------------------------
186  !
187  !++ Private procedure
188  !
189  !-----------------------------------------------------------------------------
190  !
191  !++ Private parameters & variables
192  !
193  !-----------------------------------------------------------------------------
194 contains
195  !-----------------------------------------------------------------------------
197  subroutine atmos_dyn_tstep_large_setup( &
198  Tstep_large_type, &
199  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
200  mflx_hi )
203  use scale_index
204  use scale_prc, only: &
205  prc_abort
209  implicit none
210  character(len=*), intent(in) :: Tstep_large_type
211  real(RP), intent(inout) :: DENS(ka,ia,ja)
212  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
213  real(RP), intent(inout) :: MOMX(ka,ia,ja)
214  real(RP), intent(inout) :: MOMY(ka,ia,ja)
215  real(RP), intent(inout) :: RHOT(ka,ia,ja)
216  real(RP), intent(inout) :: QTRC(ka,ia,ja,qa)
217  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
218  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
219 
220  !---------------------------------------------------------------------------
221 
222  select case( tstep_large_type )
223  case( 'FVM-HEVE' )
225  dens, momz, momx, momy, rhot, qtrc, prog, &
226  mflx_hi )
228  case default
229  log_error("ATMOS_DYN_Tstep_large_setup",*) 'ATMOS_DYN_Tstep_large_type is invalid: ', tstep_large_type
230  call prc_abort
231  end select
232 
233  return
234  end subroutine atmos_dyn_tstep_large_setup
235 
subroutine, public atmos_dyn_tstep_large_fvm_heve_setup(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Setup.
integer, public va
Definition: scale_index.F90:35
subroutine, public atmos_dyn_tstep_large_fvm_heve(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, QTRC0, 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, DTLS, DTSS, Llast)
Dynamical Process.
integer, public ia
of whole cells: x, local, with HALO
integer, public qa
integer, public ja
of whole cells: y, local, with HALO
module Atmosphere / Dynamical scheme
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
module profiler
Definition: scale_prof.F90:11
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Register.
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
module STDIO
Definition: scale_io.F90:10
procedure(large), pointer, public atmos_dyn_tstep_large