SCALE-RM
scale_atmos_dyn_tinteg_large_euler.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_index
26  use scale_tracer
27 
28 #ifdef DEBUG
29  use scale_debug, only: &
30  check
31  use scale_const, only: &
32  undef => const_undef, &
33  iundef => const_undef2
34 #endif
35  !-----------------------------------------------------------------------------
36  implicit none
37  private
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public procedure
41  !
44 
45  !-----------------------------------------------------------------------------
46  !
47  !++ Public parameters & variables
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private procedure
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private parameters & variables
56  !
57  !-----------------------------------------------------------------------------
58 contains
59 
60  !-----------------------------------------------------------------------------
63  tinteg_type )
64  use scale_process, only: &
66  use scale_const, only: &
67  undef => const_undef
68  implicit none
69 
70  character(len=*) :: tinteg_type
71 
72  integer :: iv
73  !---------------------------------------------------------------------------
74 
75  if ( tinteg_type /= 'EULER' ) then
76  write(*,*) 'xxx TINTEG_LARGE_TYPE is not EULER. Check!'
77  call prc_mpistop
78  end if
79 
80  return
82 
83  !-----------------------------------------------------------------------------
85  subroutine atmos_dyn_tinteg_large_euler( &
86  DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, &
87  DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, &
88  mflx_hi, tflx_hi, &
89  num_diff, num_diff_q, &
90  DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
91  CORIOLI, &
92  CDZ, CDX, CDY, FDZ, FDX, FDY, &
93  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
94  PHI, GSQRT, &
95  J13G, J23G, J33G, MAPF, &
96  AQ_CV, &
97  REF_dens, REF_pott, REF_qv, REF_pres, &
98  BND_W, BND_E, BND_S, BND_N, &
99  ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, &
100  DAMP_DENS, DAMP_VELZ, DAMP_VELX, &
101  DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
102  DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, &
103  DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
104  divdmp_coef, &
105  FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, &
106  FLAG_FCT_ALONG_STREAM, &
107  USE_AVERAGE, &
108  DTL, DTS )
109  use scale_const, only: &
110  rdry => const_rdry, &
111  rvap => const_rvap, &
112  cvdry => const_cvdry
113  use scale_comm, only: &
114  comm_vars8, &
115  comm_wait
116  use scale_gridtrans, only: &
117  i_xyz, &
118  i_xyw, &
119  i_uyz, &
120  i_xvz, &
121  i_xy, &
122  i_uy, &
123  i_xv
124  use scale_atmos_dyn_common, only: &
128  use scale_atmos_dyn_fvm_flux_ud1, only: &
132  use scale_atmos_dyn_fvm_flux, only: &
136  use scale_atmos_dyn_tstep_large, only: &
138  use scale_atmos_boundary, only: &
139  bnd_qa, &
140  bnd_smoother_fact => atmos_boundary_smoother_fact
141  use scale_history, only: &
142 #ifdef hist_tend
143  hist_in, &
144 #endif
145  hist_switch
146  implicit none
147 
148  real(RP), intent(inout) :: DENS(ka,ia,ja)
149  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
150  real(RP), intent(inout) :: MOMX(ka,ia,ja)
151  real(RP), intent(inout) :: MOMY(ka,ia,ja)
152  real(RP), intent(inout) :: RHOT(ka,ia,ja)
153  real(RP), intent(inout) :: QTRC(ka,ia,ja,qa)
154  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
155 
156  real(RP), intent(inout) :: DENS_av(ka,ia,ja)
157  real(RP), intent(inout) :: MOMZ_av(ka,ia,ja)
158  real(RP), intent(inout) :: MOMX_av(ka,ia,ja)
159  real(RP), intent(inout) :: MOMY_av(ka,ia,ja)
160  real(RP), intent(inout) :: RHOT_av(ka,ia,ja)
161  real(RP), intent(inout) :: QTRC_av(ka,ia,ja,qa)
162 
163  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
164  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
165 
166  real(RP), intent(out) :: num_diff (ka,ia,ja,5,3)
167  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
168 
169  real(RP), intent(in) :: DENS_tp(ka,ia,ja)
170  real(RP), intent(in) :: MOMZ_tp(ka,ia,ja)
171  real(RP), intent(in) :: MOMX_tp(ka,ia,ja)
172  real(RP), intent(in) :: MOMY_tp(ka,ia,ja)
173  real(RP), intent(in) :: RHOT_tp(ka,ia,ja)
174  real(RP), intent(in) :: RHOQ_tp(ka,ia,ja,qa)
175 
176  real(RP), intent(in) :: CORIOLI(ia,ja)
177 
178  real(RP), intent(in) :: CDZ (ka)
179  real(RP), intent(in) :: CDX (ia)
180  real(RP), intent(in) :: CDY (ja)
181  real(RP), intent(in) :: FDZ (ka-1)
182  real(RP), intent(in) :: FDX (ia-1)
183  real(RP), intent(in) :: FDY (ja-1)
184  real(RP), intent(in) :: RCDZ(ka)
185  real(RP), intent(in) :: RCDX(ia)
186  real(RP), intent(in) :: RCDY(ja)
187  real(RP), intent(in) :: RFDZ(ka-1)
188  real(RP), intent(in) :: RFDX(ia-1)
189  real(RP), intent(in) :: RFDY(ja-1)
190 
191  real(RP), intent(in) :: PHI (ka,ia,ja)
192  real(RP), intent(in) :: GSQRT(ka,ia,ja,7)
193  real(RP), intent(in) :: J13G (ka,ia,ja,7)
194  real(RP), intent(in) :: J23G (ka,ia,ja,7)
195  real(RP), intent(in) :: J33G
196  real(RP), intent(in) :: MAPF (ia,ja,2,4)
197 
198  real(RP), intent(in) :: AQ_CV(qqa)
199 
200  real(RP), intent(in) :: REF_dens(ka,ia,ja)
201  real(RP), intent(in) :: REF_pott(ka,ia,ja)
202  real(RP), intent(in) :: REF_qv (ka,ia,ja)
203  real(RP), intent(in) :: REF_pres(ka,ia,ja)
204 
205  logical, intent(in) :: BND_W
206  logical, intent(in) :: BND_E
207  logical, intent(in) :: BND_S
208  logical, intent(in) :: BND_N
209 
210  real(RP), intent(in) :: ND_COEF
211  real(RP), intent(in) :: ND_COEF_Q
212  integer, intent(in) :: ND_ORDER
213  real(RP), intent(in) :: ND_SFC_FACT
214  logical, intent(in) :: ND_USE_RS
215 
216  real(RP), intent(in) :: DAMP_DENS(ka,ia,ja)
217  real(RP), intent(in) :: DAMP_VELZ(ka,ia,ja)
218  real(RP), intent(in) :: DAMP_VELX(ka,ia,ja)
219  real(RP), intent(in) :: DAMP_VELY(ka,ia,ja)
220  real(RP), intent(in) :: DAMP_POTT(ka,ia,ja)
221  real(RP), intent(in) :: DAMP_QTRC(ka,ia,ja,bnd_qa)
222 
223  real(RP), intent(in) :: DAMP_alpha_DENS(ka,ia,ja)
224  real(RP), intent(in) :: DAMP_alpha_VELZ(ka,ia,ja)
225  real(RP), intent(in) :: DAMP_alpha_VELX(ka,ia,ja)
226  real(RP), intent(in) :: DAMP_alpha_VELY(ka,ia,ja)
227  real(RP), intent(in) :: DAMP_alpha_POTT(ka,ia,ja)
228  real(RP), intent(in) :: DAMP_alpha_QTRC(ka,ia,ja,bnd_qa)
229 
230  real(RP), intent(in) :: divdmp_coef
231 
232  logical, intent(in) :: FLAG_FCT_MOMENTUM
233  logical, intent(in) :: FLAG_FCT_T
234  logical, intent(in) :: FLAG_FCT_TRACER
235  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
236 
237  logical, intent(in) :: USE_AVERAGE
238 
239  real(DP), intent(in) :: DTL
240  real(DP), intent(in) :: DTS
241 
242  call atmos_dyn_tstep_large( &
243  dens, momz, momx, momy, rhot, qtrc, prog, & ! (inout)
244  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! (inout)
245  mflx_hi, tflx_hi, & ! (out)
246  num_diff, num_diff_q, & ! (out, work)
247  qtrc, & ! (in)
248  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! (in)
249  corioli, & ! (in)
250  cdz, cdx, cdy, fdz, fdx, fdy, & ! (in)
251  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
252  phi, gsqrt, & ! (in)
253  j13g, j23g, j33g, mapf, & ! (in)
254  aq_cv, & ! (in)
255  ref_dens, ref_pott, ref_qv, ref_pres, & ! (in)
256  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
257  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, & ! (in)
258  damp_dens, damp_velz, damp_velx, & ! (in)
259  damp_vely, damp_pott, damp_qtrc, & ! (in)
260  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! (in)
261  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! (in)
262  divdmp_coef, & ! (in)
263  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! (in)
264  flag_fct_along_stream, & ! (in)
265  use_average, & ! (in)
266  dtl, dts, .true. ) ! (in)
267 
268  return
269  end subroutine atmos_dyn_tinteg_large_euler
270 
module DEBUG
Definition: scale_debug.F90:13
integer, public i_xvz
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_dyn_numfilter_coef_q(num_diff_q, DENS, QTRC, CDZ, CDX, CDY, dt, REF_qv, iq, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
integer, public va
Definition: scale_index.F90:38
subroutine, public atmos_dyn_numfilter_coef(num_diff, DENS, MOMZ, MOMX, MOMY, RHOT, CDZ, CDX, CDY, FDZ, FDX, FDY, DT, REF_dens, REF_pott, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
integer, public i_xy
module Atmosphere / Dynamical scheme
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of x whole cells (local, with HALO)
module GRIDTRANS
integer, public ka
of z whole cells (local, with HALO)
integer, public i_uy
integer, public i_xyw
integer, public qqa
module COMMUNICATION
Definition: scale_comm.F90:23
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module Atmosphere / Dynamics common
module PROCESS
subroutine, public atmos_dyn_tinteg_large_euler_setup(tinteg_type)
Setup.
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
integer, public i_xv
module profiler
Definition: scale_prof.F90:10
module scale_atmos_dyn_fvm_flux
integer, public i_uyz
module PRECISION
subroutine, public atmos_dyn_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, dt, flag_vect)
Flux Correction Transport Limiter.
module ATMOSPHERE / Boundary treatment
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
module HISTORY
integer, public i_xyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz
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_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, DTL, DTS)
RK3.
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz
procedure(large), pointer, public atmos_dyn_tstep_large
module scale_atmos_dyn_fvm_flux_ud1
real(rp), public atmos_boundary_smoother_fact
integer, public ja
of y whole cells (local, with HALO)