SCALE-RM
scale_atmos_dyn_tstep_large.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 large( &
38  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
39  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
40  mflx_hi, tflx_hi, &
41  num_diff, num_diff_q, &
42  QTRC0, &
43  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
44  CORIOLI, &
45  CDZ, CDX, CDY, FDZ, FDX, FDY, &
46  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
47  PHI, GSQRT, &
48  J13G, J23G, J33G, MAPF, &
49  AQ_CV, &
50  REF_dens, REF_pott, REF_qv, REF_pres, &
51  BND_W, BND_E, BND_S, BND_N, & ! (in)
52  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, &
53  damp_dens, damp_velz, damp_velx, &
54  damp_vely, damp_pott, damp_qtrc, &
55  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, &
56  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, &
57  divdmp_coef, &
58  flag_fct_momentum, flag_fct_t, flag_fct_tracer, &
59  flag_fct_along_stream, &
60  use_average, &
61  dtl, dts, llast )
62  use scale_precision
64  use scale_index
65  use scale_tracer
66  use scale_atmos_boundary, only: &
67  bnd_qa
68  real(RP), intent(inout) :: dens(ka,ia,ja)
69  real(RP), intent(inout) :: momz(ka,ia,ja)
70  real(RP), intent(inout) :: momx(ka,ia,ja)
71  real(RP), intent(inout) :: momy(ka,ia,ja)
72  real(RP), intent(inout) :: rhot(ka,ia,ja)
73  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
74  real(RP), intent(inout) :: prog(ka,ia,ja,va)
75 
76  real(RP), intent(inout) :: dens_av(ka,ia,ja)
77  real(RP), intent(inout) :: momz_av(ka,ia,ja)
78  real(RP), intent(inout) :: momx_av(ka,ia,ja)
79  real(RP), intent(inout) :: momy_av(ka,ia,ja)
80  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
81  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
82 
83  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
84  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
85 
86  real(RP), intent(out) :: num_diff(ka,ia,ja,5,3)
87  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
88 
89  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa)
90 
91  real(RP), intent(in) :: dens_tp(ka,ia,ja)
92  real(RP), intent(in) :: momz_tp(ka,ia,ja)
93  real(RP), intent(in) :: momx_tp(ka,ia,ja)
94  real(RP), intent(in) :: momy_tp(ka,ia,ja)
95  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
96  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
97 
98  real(RP), intent(in) :: corioli(ia,ja)
99 
100  real(RP), intent(in) :: cdz (ka)
101  real(RP), intent(in) :: cdx (ia)
102  real(RP), intent(in) :: cdy (ja)
103  real(RP), intent(in) :: fdz (ka-1)
104  real(RP), intent(in) :: fdx (ia-1)
105  real(RP), intent(in) :: fdy (ja-1)
106  real(RP), intent(in) :: rcdz(ka)
107  real(RP), intent(in) :: rcdx(ia)
108  real(RP), intent(in) :: rcdy(ja)
109  real(RP), intent(in) :: rfdz(ka-1)
110  real(RP), intent(in) :: rfdx(ia-1)
111  real(RP), intent(in) :: rfdy(ja-1)
112 
113  real(RP), intent(in) :: phi (ka,ia,ja)
114  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
115  real(RP), intent(in) :: j13g (ka,ia,ja,7)
116  real(RP), intent(in) :: j23g (ka,ia,ja,7)
117  real(RP), intent(in) :: j33g
118  real(RP), intent(in) :: mapf (ia,ja,2,4)
119 
120  real(RP), intent(in) :: aq_cv(qqa)
121 
122  real(RP), intent(in) :: ref_dens(ka,ia,ja)
123  real(RP), intent(in) :: ref_pott(ka,ia,ja)
124  real(RP), intent(in) :: ref_qv (ka,ia,ja)
125  real(RP), intent(in) :: ref_pres(ka,ia,ja)
126 
127  logical, intent(in) :: bnd_w
128  logical, intent(in) :: bnd_e
129  logical, intent(in) :: bnd_s
130  logical, intent(in) :: bnd_n
131 
132  real(RP), intent(in) :: nd_coef
133  real(RP), intent(in) :: nd_coef_q
134  integer, intent(in) :: nd_order
135  real(RP), intent(in) :: nd_sfc_fact
136  logical, intent(in) :: nd_use_rs
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 
145  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
146  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
147  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
148  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
149  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
150  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
151 
152  real(RP), intent(in) :: divdmp_coef
153 
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 
159  logical, intent(in) :: use_average
160 
161  real(DP), intent(in) :: dtl
162  real(DP), intent(in) :: dts
163 
164  logical , intent(in) :: llast
165  end subroutine large
166 
167  end interface
168 
169  procedure(large), pointer :: atmos_dyn_tstep_large => null()
170  public :: atmos_dyn_tstep_large
171 
172  !-----------------------------------------------------------------------------
173  !
174  !++ Public parameters & variables
175  !
176  !-----------------------------------------------------------------------------
177  !
178  !++ Private procedure
179  !
180  !-----------------------------------------------------------------------------
181  !
182  !++ Private parameters & variables
183  !
184  !-----------------------------------------------------------------------------
185 contains
186  !-----------------------------------------------------------------------------
188  subroutine atmos_dyn_tstep_large_setup( &
189  Tstep_large_type, &
190  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
191  mflx_hi )
193  use scale_grid_index
194  use scale_index
195  use scale_process, only: &
200  implicit none
201  character(len=*), intent(in) :: Tstep_large_type
202  real(RP), intent(inout) :: DENS(ka,ia,ja)
203  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
204  real(RP), intent(inout) :: MOMX(ka,ia,ja)
205  real(RP), intent(inout) :: MOMY(ka,ia,ja)
206  real(RP), intent(inout) :: RHOT(ka,ia,ja)
207  real(RP), intent(inout) :: QTRC(ka,ia,ja,qa)
208  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
209  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
210 
211  !---------------------------------------------------------------------------
212 
213  select case ( tstep_large_type )
214  case ( 'FVM-HEVE' )
216  dens, momz, momx, momy, rhot, qtrc, prog, &
217  mflx_hi )
219  case default
220  write(*,*) 'xxx ATMOS_DYN_Tstep_large_type is invalid: ', tstep_large_type
221  call prc_mpistop
222  end select
223 
224  return
225  end subroutine atmos_dyn_tstep_large_setup
226 
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_CV, 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, 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, divdmp_coef, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, DTLS, DTSS, Llast)
Dynamical Process.
subroutine, public atmos_dyn_tstep_large_fvm_heve_setup(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Setup.
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
module Atmosphere / Dynamical scheme
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)
integer, public qqa
module PROCESS
module profiler
Definition: scale_prof.F90:10
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Register.
module PRECISION
module ATMOSPHERE / Boundary treatment
procedure(large), pointer, public atmos_dyn_tstep_large
integer, public ja
of y whole cells (local, with HALO)