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