SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tstep_short_fvm_hevi Module Reference

module Atmosphere / Dynamics RK More...

Functions/Subroutines

subroutine, public atmos_dyn_tstep_short_fvm_hevi_regist (ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
 Register. More...
 
subroutine, public atmos_dyn_tstep_short_fvm_hevi_setup
 Setup. More...
 
subroutine, public atmos_dyn_tstep_short_fvm_hevi (DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, PROG_RK, mflx_hi, tflx_hi, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, DENS, MOMZ, MOMX, MOMY, RHOT, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, PROG0, PROG, DPRES0, RT2P, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_dens, REF_rhot, BND_W, BND_E, BND_S, BND_N, dtrk, last)
 
subroutine solve_direct (C, F1, F2, F3)
 

Detailed Description

module Atmosphere / Dynamics RK

Description
HEVI FVM scheme for Atmospheric dynamical process
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
DENS_t_advch tendency of density (horiz. advection) (w/ HIST_TEND) kg/m3/s advch_t
DENS_t_advcv tendency of density (vert. advection) (w/ HIST_TEND) kg/m3/s advcv_t
MOMX_t_advch tendency of momentum x (horiz. advection) (w/ HIST_TEND) kg/m2/s2 advch_t
MOMX_t_advcv tendency of momentum x (vert. advection) (w/ HIST_TEND) kg/m2/s2 advcv_t
MOMX_t_cf tendency of momentum x (coliolis force) (w/ HIST_TEND) kg/m2/s2 cf_t
MOMX_t_ddiv tendency of momentum x (divergence damping) (w/ HIST_TEND) kg/m2/s2 ddiv_t
MOMX_t_pg tendency of momentum x (pressure gradient) (w/ HIST_TEND) kg/m2/s2 pg_t
MOMY_t_advch tendency of momentum y (horiz. advection) (w/ HIST_TEND) kg/m2/s2 advch_t
MOMY_t_advcv tendency of momentum y (vert. advection) (w/ HIST_TEND) kg/m2/s2 advcv_t
MOMY_t_cf tendency of momentum y (coliolis force) (w/ HIST_TEND) kg/m2/s2 cf_t
MOMY_t_ddiv tendency of momentum y (divergence damping) (w/ HIST_TEND) kg/m2/s2 ddiv_t
MOMY_t_pg tendency of momentum y (pressure gradient) (w/ HIST_TEND) kg/m2/s2 pg_t
MOMZ_t_advch tendency of momentum z (horiz. advection) (w/ HIST_TEND) kg/m2/s2 advch_t
MOMZ_t_advcv tendency of momentum z (vert. advection) (w/ HIST_TEND) kg/m2/s2 advcv_t
MOMZ_t_ddiv tendency of momentum z (divergence damping) (w/ HIST_TEND) kg/m2/s2 ddiv_t
MOMZ_t_pg tendency of momentum z (pressure gradient) (w/ HIST_TEND) kg/m2/s2 pg_t
MOMZ_t_wdamp tendency of momentum z (Rayleigh damping) (w/ HIST_TEND) kg/m2/s2 wdmp_t
RHOT_t_advch tendency of rho*theta (horiz. advection) (w/ HIST_TEND) K kg/m3/s advch_t
RHOT_t_advcv tendency of rho*theta (vert. advection) (w/ HIST_TEND) K kg/m3/s advcv_t

Function/Subroutine Documentation

◆ atmos_dyn_tstep_short_fvm_hevi_regist()

subroutine, public scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi_regist ( character(len=*), intent(in)  ATMOS_DYN_TYPE,
integer, intent(out)  VA_out,
character(len=h_short), dimension(:), intent(out)  VAR_NAME,
character(len=h_mid), dimension(:), intent(out)  VAR_DESC,
character(len=h_short), dimension(:), intent(out)  VAR_UNIT 
)

Register.

Parameters
[out]va_outnumber of prognostic variables
[out]var_namename of the variables
[out]var_descdesc. of the variables
[out]var_unitunit of the variables

Definition at line 86 of file scale_atmos_dyn_tstep_short_fvm_hevi.F90.

References scale_prc::prc_abort().

Referenced by scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_regist().

86  use scale_prc, only: &
87  prc_abort
88  implicit none
89 
90  character(len=*), intent(in) :: atmos_dyn_type
91  integer, intent(out) :: va_out
92  character(len=H_SHORT), intent(out) :: var_name(:)
93  character(len=H_MID), intent(out) :: var_desc(:)
94  character(len=H_SHORT), intent(out) :: var_unit(:)
95  !---------------------------------------------------------------------------
96 
97  if ( atmos_dyn_type /= 'FVM-HEVI' .AND. atmos_dyn_type /= 'HEVI' ) then
98  log_error("ATMOS_DYN_Tstep_short_fvm_hevi_regist",*) 'ATMOS_DYN_TYPE is not FVM-HEVI. Check!'
99  call prc_abort
100  endif
101 
102  va_out = va_fvm_hevi
103  var_name(:) = ""
104  var_desc(:) = ""
105  var_unit(:) = ""
106 
107  log_newline
108  log_info("ATMOS_DYN_Tstep_short_fvm_hevi_regist",*) 'Register additional prognostic variables (HEVI)'
109  if ( va_out < 1 ) then
110  log_info_cont(*) '=> nothing.'
111  endif
112 
113  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tstep_short_fvm_hevi_setup()

subroutine, public scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi_setup ( )

Setup.

Definition at line 119 of file scale_atmos_dyn_tstep_short_fvm_hevi.F90.

Referenced by scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_regist().

119  implicit none
120  !---------------------------------------------------------------------------
121 
122  log_info("ATMOS_DYN_Tstep_short_fvm_hevi_setup",*) 'HEVI Setup'
123 
124  return
Here is the caller graph for this function:

◆ atmos_dyn_tstep_short_fvm_hevi()

subroutine, public scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi ( real(rp), dimension(ka,ia,ja), intent(out)  DENS_RK,
real(rp), dimension(ka,ia,ja), intent(out)  MOMZ_RK,
real(rp), dimension(ka,ia,ja), intent(out)  MOMX_RK,
real(rp), dimension(ka,ia,ja), intent(out)  MOMY_RK,
real(rp), dimension(ka,ia,ja), intent(out)  RHOT_RK,
real(rp), dimension(ka,ia,ja,va), intent(out)  PROG_RK,
real(rp), dimension(ka,ia,ja,3), intent(inout)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(out)  tflx_hi,
real(rp), dimension(ka,ia,ja), intent(in), target  DENS0,
real(rp), dimension(ka,ia,ja), intent(in), target  MOMZ0,
real(rp), dimension(ka,ia,ja), intent(in), target  MOMX0,
real(rp), dimension(ka,ia,ja), intent(in), target  MOMY0,
real(rp), dimension(ka,ia,ja), intent(in), target  RHOT0,
real(rp), dimension(ka,ia,ja), intent(in)  DENS,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT,
real(rp), dimension(ka,ia,ja), intent(in)  DENS_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY_t,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT_t,
real(rp), dimension(ka,ia,ja,va), intent(in)  PROG0,
real(rp), dimension (ka,ia,ja,va), intent(in)  PROG,
real(rp), dimension(ka,ia,ja), intent(in)  DPRES0,
real(rp), dimension(ka,ia,ja), intent(in)  RT2P,
real(rp), dimension(1,ia,ja), intent(in)  CORIOLI,
real(rp), dimension(ka,ia,ja,5,3), intent(in)  num_diff,
real(rp), dimension(ka), intent(in)  wdamp_coef,
real(rp), intent(in)  divdmp_coef,
real(rp), dimension(ka,ia,ja), intent(in)  DDIV,
logical, intent(in)  FLAG_FCT_MOMENTUM,
logical, intent(in)  FLAG_FCT_T,
logical, intent(in)  FLAG_FCT_ALONG_STREAM,
real(rp), dimension(ka), intent(in)  CDZ,
real(rp), dimension(ka-1), intent(in)  FDZ,
real(rp), dimension(ia-1), intent(in)  FDX,
real(rp), dimension(ja-1), intent(in)  FDY,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
real(rp), dimension(ka-1), intent(in)  RFDZ,
real(rp), dimension(ia-1), intent(in)  RFDX,
real(rp), dimension(ja-1), intent(in)  RFDY,
real(rp), dimension (ka,ia,ja), intent(in)  PHI,
real(rp), dimension (ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ka,ia,ja,7), intent(in)  J13G,
real(rp), dimension (ka,ia,ja,7), intent(in)  J23G,
real(rp), intent(in)  J33G,
real(rp), dimension (ia,ja,2,4), intent(in)  MAPF,
real(rp), dimension(ka,ia,ja), intent(in)  REF_dens,
real(rp), dimension(ka,ia,ja), intent(in)  REF_rhot,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
real(rp), intent(in)  dtrk,
logical, intent(in)  last 
)
Parameters
[in]phigeopotential
[in]gsqrtvertical metrics {G}^1/2
[in]j13g(1,3) element of Jacobian matrix
[in]j23g(2,3) element of Jacobian matrix
[in]j33g(3,3) element of Jacobian matrix
[in]mapfmap factor
[in]ref_densreference density

Definition at line 146 of file scale_atmos_dyn_tstep_short_fvm_hevi.F90.

References scale_atmos_dyn_common::atmos_dyn_fct(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_valuew_z, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_uyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_xvz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_xyw, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_uyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_xvz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_xyw, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_uyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xvz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyw, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_uyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xvz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyw, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_uyz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xvz, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyw, scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyz, scale_debug::check(), scale_index::i_dens, scale_index::i_momx, scale_index::i_momy, scale_index::i_momz, scale_index::i_rhot, scale_atmos_grid_cartesc_index::i_uv, scale_atmos_grid_cartesc_index::i_uvz, scale_atmos_grid_cartesc_index::i_uy, scale_atmos_grid_cartesc_index::i_uyw, scale_atmos_grid_cartesc_index::i_uyz, scale_atmos_grid_cartesc_index::i_xv, scale_atmos_grid_cartesc_index::i_xvw, scale_atmos_grid_cartesc_index::i_xvz, scale_atmos_grid_cartesc_index::i_xy, scale_atmos_grid_cartesc_index::i_xyw, scale_atmos_grid_cartesc_index::i_xyz, scale_atmos_grid_cartesc_index::iblock, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieh, scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::jblock, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeh, scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, solve_direct(), scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_regist().

147  use scale_const, only: &
148 #ifdef dry
149  rdry => const_rdry, &
150  cvdry => const_cvdry, &
151  cpdry => const_cpdry, &
152 #endif
153  eps => const_eps, &
154  grav => const_grav, &
155  p00 => const_pre00
156  use scale_atmos_dyn_common, only: &
158  use scale_atmos_dyn_fvm_flux, only: &
178 #ifdef HIST_TEND
179  use scale_file_history, only: &
180  file_history_in
181 #endif
182  implicit none
183 
184  real(RP), intent(out) :: dens_rk(ka,ia,ja) ! prognostic variables
185  real(RP), intent(out) :: momz_rk(ka,ia,ja) !
186  real(RP), intent(out) :: momx_rk(ka,ia,ja) !
187  real(RP), intent(out) :: momy_rk(ka,ia,ja) !
188  real(RP), intent(out) :: rhot_rk(ka,ia,ja) !
189 
190  real(RP), intent(out) :: prog_rk(ka,ia,ja,va) !
191 
192  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3) ! rho * vel(x,y,z)
193  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3) ! rho * theta * vel(x,y,z)
194 
195  real(RP), intent(in),target :: dens0(ka,ia,ja) ! prognostic variables
196  real(RP), intent(in),target :: momz0(ka,ia,ja) ! at previous dynamical time step
197  real(RP), intent(in),target :: momx0(ka,ia,ja) !
198  real(RP), intent(in),target :: momy0(ka,ia,ja) !
199  real(RP), intent(in),target :: rhot0(ka,ia,ja) !
200 
201  real(RP), intent(in) :: dens(ka,ia,ja) ! prognostic variables
202  real(RP), intent(in) :: momz(ka,ia,ja) ! at previous RK step
203  real(RP), intent(in) :: momx(ka,ia,ja) !
204  real(RP), intent(in) :: momy(ka,ia,ja) !
205  real(RP), intent(in) :: rhot(ka,ia,ja) !
206 
207  real(RP), intent(in) :: dens_t(ka,ia,ja)
208  real(RP), intent(in) :: momz_t(ka,ia,ja)
209  real(RP), intent(in) :: momx_t(ka,ia,ja)
210  real(RP), intent(in) :: momy_t(ka,ia,ja)
211  real(RP), intent(in) :: rhot_t(ka,ia,ja)
212 
213  real(RP), intent(in) :: prog0(ka,ia,ja,va)
214  real(RP), intent(in) :: prog (ka,ia,ja,va)
215 
216  real(RP), intent(in) :: dpres0(ka,ia,ja)
217  real(RP), intent(in) :: rt2p(ka,ia,ja)
218  real(RP), intent(in) :: corioli(1,ia,ja)
219  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
220  real(RP), intent(in) :: wdamp_coef(ka)
221  real(RP), intent(in) :: divdmp_coef
222  real(RP), intent(in) :: ddiv(ka,ia,ja)
223 
224  logical, intent(in) :: flag_fct_momentum
225  logical, intent(in) :: flag_fct_t
226  logical, intent(in) :: flag_fct_along_stream
227 
228  real(RP), intent(in) :: cdz(ka)
229  real(RP), intent(in) :: fdz(ka-1)
230  real(RP), intent(in) :: fdx(ia-1)
231  real(RP), intent(in) :: fdy(ja-1)
232  real(RP), intent(in) :: rcdz(ka)
233  real(RP), intent(in) :: rcdx(ia)
234  real(RP), intent(in) :: rcdy(ja)
235  real(RP), intent(in) :: rfdz(ka-1)
236  real(RP), intent(in) :: rfdx(ia-1)
237  real(RP), intent(in) :: rfdy(ja-1)
238 
239  real(RP), intent(in) :: phi (ka,ia,ja)
240  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
241  real(RP), intent(in) :: j13g (ka,ia,ja,7)
242  real(RP), intent(in) :: j23g (ka,ia,ja,7)
243  real(RP), intent(in) :: j33g
244  real(RP), intent(in) :: mapf (ia,ja,2,4)
245  real(RP), intent(in) :: ref_dens(ka,ia,ja)
246  real(RP), intent(in) :: ref_rhot(ka,ia,ja)
247 
248  logical, intent(in) :: bnd_w
249  logical, intent(in) :: bnd_e
250  logical, intent(in) :: bnd_s
251  logical, intent(in) :: bnd_n
252 
253  real(RP), intent(in) :: dtrk
254  logical, intent(in) :: last
255 
256 
257  ! diagnostic variables (work space)
258  real(RP) :: pott(ka,ia,ja) ! potential temperature [K]
259  real(RP) :: dpres(ka,ia,ja) ! pressure deviation from reference pressure
260 
261  real(RP) :: qflx_hi (ka,ia,ja,3)
262  real(RP) :: qflx_j13(ka,ia,ja)
263  real(RP) :: qflx_j23(ka,ia,ja)
264 
265  real(RP) :: advch ! horizontal advection
266  real(RP) :: advcv ! vertical advection
267  real(RP) :: wdmp ! Raileight damping
268  real(RP) :: div ! divergence damping
269  real(RP) :: pg ! pressure gradient force
270  real(RP) :: cf ! colioris force
271 #ifdef HIST_TEND
272  real(RP) :: advch_t(ka,ia,ja,5)
273  real(RP) :: advcv_t(ka,ia,ja,5)
274  real(RP) :: wdmp_t(ka,ia,ja)
275  real(RP) :: ddiv_t(ka,ia,ja,3)
276  real(RP) :: pg_t(ka,ia,ja,3)
277  real(RP) :: cf_t(ka,ia,ja,2)
278  logical :: lhist
279 #endif
280 
281  ! for implicit solver
282  real(RP) :: a(ka)
283  real(RP) :: b
284  real(RP) :: sr(ka,ia,ja)
285  real(RP) :: sw(ka,ia,ja)
286  real(RP) :: st(ka,ia,ja)
287  real(RP) :: pt(ka)
288  real(RP) :: c(kmax-1)
289 
290  real(RP) :: f1(ka)
291  real(RP) :: f2(ka)
292  real(RP) :: f3(ka)
293 
294  integer :: iis, iie, jjs, jje
295  integer :: k, i, j
296  integer :: iss, iee
297 
298 #ifdef DEBUG
299  pott(:,:,:) = undef
300  dpres(:,:,:) = undef
301 
302  pt(:) = undef
303 
304  qflx_hi(:,:,:,:) = undef
305  qflx_j13(:,:,:) = undef
306  qflx_j23(:,:,:) = undef
307 #endif
308 
309 #if defined DEBUG || defined QUICKDEBUG
310  dens_rk( 1:ks-1,:,:) = undef
311  dens_rk(ke+1:ka ,:,:) = undef
312  momz_rk( 1:ks-1,:,:) = undef
313  momz_rk(ke+1:ka ,:,:) = undef
314  momx_rk( 1:ks-1,:,:) = undef
315  momx_rk(ke+1:ka ,:,:) = undef
316  momy_rk( 1:ks-1,:,:) = undef
317  momy_rk(ke+1:ka ,:,:) = undef
318  rhot_rk( 1:ks-1,:,:) = undef
319  rhot_rk(ke+1:ka ,:,:) = undef
320  prog_rk( 1:ks-1,:,:,:) = undef
321  prog_rk(ke+1:ka ,:,:,:) = undef
322 #endif
323 
324 #ifdef HIST_TEND
325  lhist = last
326 #endif
327 
328 #ifdef PROFILE_FIPP
329  call fipp_start()
330 #endif
331 
332  ifs_off = 1
333  jfs_off = 1
334  if ( bnd_w ) ifs_off = 0
335  if ( bnd_s ) jfs_off = 0
336 
337 
338  do jjs = js, je, jblock
339  jje = jjs+jblock-1
340  do iis = is, ie, iblock
341  iie = iis+iblock-1
342 
343  profile_start("hevi_pres")
344  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
345  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,DPRES0,RT2P,RHOT,REF_rhot,DPRES,DENS,PHI)
346  do j = jjs, jje+1
347  do i = iis, iie+1
348  do k = ks, ke
349 #ifdef DEBUG
350  call check( __line__, dpres0(k,i,j) )
351  call check( __line__, rt2p(k,i,j) )
352  call check( __line__, rhot(k,i,j) )
353  call check( __line__, ref_rhot(k,i,j) )
354 #endif
355  dpres(k,i,j) = dpres0(k,i,j) + rt2p(k,i,j) * ( rhot(k,i,j) - ref_rhot(k,i,j) )
356  enddo
357  dpres(ks-1,i,j) = dpres0(ks-1,i,j) - dens(ks,i,j) * ( phi(ks-1,i,j) - phi(ks+1,i,j) )
358  dpres(ke+1,i,j) = dpres0(ke+1,i,j) - dens(ke,i,j) * ( phi(ke+1,i,j) - phi(ke-1,i,j) )
359  enddo
360  enddo
361  profile_stop("hevi_pres")
362 
363  !##### continuity equation #####
364 
365  profile_start("hevi_mflx_z")
366  ! at (x, y, w)
367  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
368  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,GSQRT,I_XYW,MOMX,MOMY,J13G,J23G,mflx_hi,MAPF,I_XY,num_diff)
369  do j = jjs, jje
370  do i = iis-1, iie
371  mflx_hi(ks-1,i,j,zdir) = 0.0_rp
372  do k = ks, ke-1
373 #ifdef DEBUG
374  call check( __line__, momx(k+1,i ,j) )
375  call check( __line__, momx(k+1,i-1,j) )
376  call check( __line__, momx(k ,i ,j) )
377  call check( __line__, momx(k ,i+1,j) )
378  call check( __line__, momy(k+1,i,j) )
379  call check( __line__, momy(k+1,i,j-1) )
380  call check( __line__, momy(k ,i,j) )
381  call check( __line__, momy(k ,i,j-1) )
382 #endif
383  mflx_hi(k,i,j,zdir) = j13g(k,i,j,i_xyw) * 0.25_rp / mapf(i,j,2,i_xy) &
384  * ( momx(k+1,i,j)+momx(k+1,i-1,j) &
385  + momx(k ,i,j)+momx(k ,i-1,j) ) &
386  + j23g(k,i,j,i_xyw) * 0.25_rp / mapf(i,j,1,i_xy) &
387  * ( momy(k+1,i,j)+momy(k+1,i,j-1) &
388  + momy(k ,i,j)+momy(k ,i,j-1) ) &
389  + gsqrt(k,i,j,i_xyw) / ( mapf(i,j,1,i_xy)*mapf(i,j,2,i_xy) ) * num_diff(k,i,j,i_dens,zdir)
390  enddo
391  mflx_hi(ke ,i,j,zdir) = 0.0_rp
392  enddo
393  enddo
394 #ifdef DEBUG
395  k = iundef; i = iundef; j = iundef
396 #endif
397  profile_stop("hevi_mflx_z")
398 
399  profile_start("hevi_mflx_x")
400  iss = max(iis-1,is-ifs_off)
401  iee = min(iie,ieh)
402  ! at (u, y, z)
403  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
404  !$omp shared(JJS,JJE,iss,iee,KS,KE,GSQRT,I_UYZ,MOMX,num_diff,mflx_hi,MAPF,I_UY)
405  do j = jjs, jje
406  do i = iss, iee
407  do k = ks, ke
408 #ifdef DEBUG
409  call check( __line__, gsqrt(k,i,j,i_uyz) )
410  call check( __line__, momx(k,i,j) )
411  call check( __line__, num_diff(k,i,j,i_dens,xdir) )
412 #endif
413  mflx_hi(k,i,j,xdir) = gsqrt(k,i,j,i_uyz) / mapf(i,j,2,i_uy) &
414  * ( momx(k,i,j) + num_diff(k,i,j,i_dens,xdir) )
415  enddo
416  enddo
417  enddo
418 #ifdef DEBUG
419  k = iundef; i = iundef; j = iundef
420 #endif
421  profile_stop("hevi_mflx_x")
422 
423  ! at (x, v, z)
424  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
425  !$omp shared(JJS,JS,JFS_OFF,JJE,JEH,IIS,IIE,KS,KE,GSQRT,I_XVZ,MOMY,num_diff,mflx_hi,MAPF,I_XV)
426  do j = max(jjs-1,js-jfs_off), min(jje,jeh)
427  do i = iis, iie
428  do k = ks, ke
429 #ifdef DEBUG
430  call check( __line__, gsqrt(k,i,j,i_xvz) )
431  call check( __line__, momy(k,i,j) )
432  call check( __line__, num_diff(k,i,j,i_dens,ydir) )
433 #endif
434  mflx_hi(k,i,j,ydir) = gsqrt(k,i,j,i_xvz) / mapf(i,j,1,i_xv) &
435  * ( momy(k,i,j) + num_diff(k,i,j,i_dens,ydir) )
436  enddo
437  enddo
438  enddo
439 #ifdef DEBUG
440  k = iundef; i = iundef; j = iundef
441 #endif
442 
443  !--- update density
444  profile_start("hevi_sr")
445  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
446  !$omp private(i,j,k,advcv,advch) &
447 #ifdef HIST_TEND
448  !$omp shared(advch_t,lhist) &
449 #endif
450  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
451  !$omp shared(DENS0,Sr,mflx_hi,DENS_t) &
452  !$omp shared(RCDZ,RCDX,RCDY,MAPF,GSQRT,I_XY,I_XYZ)
453  do j = jjs, jje
454  do i = iis, iie
455  do k = ks, ke
456 #ifdef DEBUG
457  call check( __line__, dens0(k,i,j) )
458  call check( __line__, mflx_hi(k ,i ,j ,xdir) )
459  call check( __line__, mflx_hi(k ,i-1,j ,xdir) )
460  call check( __line__, mflx_hi(k ,i ,j ,ydir) )
461  call check( __line__, mflx_hi(k ,i ,j-1,ydir) )
462  call check( __line__, dens_t(k,i,j) )
463 #endif
464  advcv = - ( mflx_hi(k,i,j,zdir)-mflx_hi(k-1,i ,j, zdir) ) * rcdz(k)
465  advch = - ( ( mflx_hi(k,i,j,xdir)-mflx_hi(k ,i-1,j, xdir) ) * rcdx(i) &
466  + ( mflx_hi(k,i,j,ydir)-mflx_hi(k ,i, j-1,ydir) ) * rcdy(j) )
467  sr(k,i,j) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) + dens_t(k,i,j)
468 #ifdef HIST_TEND
469  if ( lhist ) then
470  advch_t(k,i,j,i_dens) = ( advch + advcv ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
471  end if
472 #endif
473  enddo
474  enddo
475  enddo
476 #ifdef DEBUG
477  k = iundef; i = iundef; j = iundef
478 #endif
479  profile_stop("hevi_sr")
480 
481  !##### momentum equation (z) #####
482 
483  ! at (x, y, z)
484  ! not that z-index is added by -1
485  profile_start("hevi_momz_qflxhi_z")
486  call atmos_dyn_fvm_fluxz_xyw( qflx_hi(:,:,:,zdir), & ! (out)
487  momz, momz, dens, & ! (in)
488  gsqrt(:,:,:,i_xyz), j33g, & ! (in)
489  num_diff(:,:,:,i_momz,zdir), & ! (in)
490  cdz, fdz, dtrk, &
491  iis, iie, jjs, jje ) ! (in)
492  profile_stop("hevi_momz_qflxhi_z")
493 
494  profile_start("hevi_momz_qflxj")
495  call atmos_dyn_fvm_fluxj13_xyw( qflx_j13, & ! (out)
496  momx, momz, dens, & ! (in)
497  gsqrt(:,:,:,i_xyz), j13g(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
498  cdz, &
499  iis, iie, jjs, jje ) ! (in)
500  call atmos_dyn_fvm_fluxj23_xyw( qflx_j23, & ! (out)
501  momy, momz, dens, & ! (in)
502  gsqrt(:,:,:,i_xyz), j23g(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
503  cdz, &
504  iis, iie, jjs, jje ) ! (in)
505  profile_stop("hevi_momz_qflxj")
506 
507  ! at (u, y, w)
508  profile_start("hevi_momz_qflxhi_x")
509  call atmos_dyn_fvm_fluxx_xyw( qflx_hi(:,:,:,xdir), & ! (out)
510  momx, momz, dens, & ! (in)
511  gsqrt(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
512  num_diff(:,:,:,i_momz,xdir), & ! (in)
513  cdz, & ! (in)
514  iis, iie, jjs, jje ) ! (in)
515  profile_stop("hevi_momz_qflxhi_x")
516 
517  ! at (x, v, w)
518  profile_start("hevi_momz_qflxhi_y")
519  call atmos_dyn_fvm_fluxy_xyw( qflx_hi(:,:,:,ydir), & ! (out)
520  momy, momz, dens, & ! (in)
521  gsqrt(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
522  num_diff(:,:,:,i_momz,ydir), & ! (in)
523  cdz, & ! (in)
524  iis, iie, jjs, jje ) ! (in)
525  profile_stop("hevi_momz_qflxhi_y")
526 
527  !--- update momentum(z)
528  profile_start("hevi_sw")
529  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
530  !$omp private(i,j,k,advcv,advch,cf,wdmp,div) &
531 #ifdef HIST_TEND
532  !$omp shared(lhist,advcv_t,advch_t,wdmp_t,ddiv_t) &
533 #endif
534  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
535  !$omp shared(qflx_hi,qflx_J13,qflx_J23,DDIV,MOMZ0,MOMZ_t,Sw) &
536  !$omp shared(RFDZ,RCDX,RCDY,FDZ,dtrk,wdamp_coef,divdmp_coef) &
537  !$omp shared(MAPF,GSQRT,I_XY,I_XYW)
538  do j = jjs, jje
539  do i = iis, iie
540  do k = ks, ke-1
541 #ifdef DEBUG
542  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
543  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
544  call check( __line__, qflx_j13(k ,i ,j) )
545  call check( __line__, qflx_j13(k-1,i ,j) )
546  call check( __line__, qflx_j23(k ,i ,j) )
547  call check( __line__, qflx_j23(k-1,i ,j) )
548  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
549  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
550  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
551  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
552  call check( __line__, ddiv(k ,i,j) )
553  call check( __line__, ddiv(k+1,i,j) )
554  call check( __line__, momz0(k,i,j) )
555  call check( __line__, momz_t(k,i,j) )
556 #endif
557  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
558  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
559  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rfdz(k)
560  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k,i-1,j ,xdir) ) * rcdx(i) &
561  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k,i ,j-1,ydir) ) * rcdy(j) ) &
562  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy)
563  wdmp = - wdamp_coef(k) * momz0(k,i,j)
564  div = divdmp_coef / dtrk * ( ddiv(k+1,i,j)-ddiv(k,i,j) ) * fdz(k) ! divergence damping
565  sw(k,i,j) = ( advcv + advch ) / gsqrt(k,i,j,i_xyw) &
566  + wdmp + div + momz_t(k,i,j)
567 #ifdef HIST_TEND
568  if ( lhist ) then
569  advcv_t(k,i,j,i_momz) = advcv / gsqrt(k,i,j,i_xyw)
570  advch_t(k,i,j,i_momz) = advch / gsqrt(k,i,j,i_xyw)
571  wdmp_t(k,i,j) = wdmp
572  ddiv_t(k,i,j,1) = div
573  endif
574 #endif
575  enddo
576  enddo
577  enddo
578  profile_stop("hevi_sw")
579 #ifdef DEBUG
580  k = iundef; i = iundef; j = iundef
581 #endif
582 
583 
584  !##### Thermodynamic Equation #####
585 
586  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
587  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,JHALO,IHALO,RHOT,DENS,POTT)
588  do j = jjs-jhalo, jje+jhalo
589  do i = iis-ihalo, iie+ihalo
590  do k = ks, ke
591 #ifdef DEBUG
592  call check( __line__, rhot(k,i,j) )
593  call check( __line__, dens(k,i,j) )
594 #endif
595  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
596  enddo
597  enddo
598  enddo
599 #ifdef DEBUG
600  k = iundef; i = iundef; j = iundef
601 #endif
602 
603  ! at (x, y, w)
604  call atmos_dyn_fvm_fluxz_xyz( tflx_hi(:,:,:,zdir), & ! (out)
605  mflx_hi(:,:,:,zdir), pott, gsqrt(:,:,:,i_xyw), & ! (in)
606  num_diff(:,:,:,i_rhot,zdir), & ! (in)
607  cdz, & ! (in)
608  iis, iie, jjs, jje ) ! (in)
609 
610  ! at (u, y, z)
611  call atmos_dyn_fvm_fluxx_xyz( tflx_hi(:,:,:,xdir), & ! (out)
612  mflx_hi(:,:,:,xdir), pott, gsqrt(:,:,:,i_uyz), & ! (in)
613  num_diff(:,:,:,i_rhot,xdir), & ! (in)
614  cdz, & ! (in)
615  iis, iie, jjs, jje ) ! (in)
616 
617  ! at (x, v, z)
618  call atmos_dyn_fvm_fluxy_xyz( tflx_hi(:,:,:,ydir), & ! (out)
619  mflx_hi(:,:,:,ydir), pott, gsqrt(:,:,:,i_xvz), & ! (in)
620  num_diff(:,:,:,i_rhot,ydir), & ! (in)
621  cdz, & ! (in)
622  iis, iie, jjs, jje ) ! (in)
623 
624 
625  profile_start("hevi_st")
626  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
627  !$omp private(i,j,k,advcv,advch) &
628 #ifdef HIST_TEND
629  !$omp shared(lhist,advcv_t,advch_t) &
630 #endif
631  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
632  !$omp shared(tflx_hi,RHOT_t,St,RCDZ,RCDX,RCDY) &
633  !$omp shared(MAPF,GSQRT,I_XY,I_XYZ)
634  do j = jjs, jje
635  do i = iis, iie
636  do k = ks, ke
637 #ifdef DEBUG
638  call check( __line__, tflx_hi(k ,i ,j ,zdir) )
639  call check( __line__, tflx_hi(k-1,i ,j ,zdir) )
640  call check( __line__, tflx_hi(k ,i ,j ,xdir) )
641  call check( __line__, tflx_hi(k ,i-1,j ,xdir) )
642  call check( __line__, tflx_hi(k ,i ,j ,ydir) )
643  call check( __line__, tflx_hi(k ,i ,j-1,ydir) )
644  call check( __line__, rhot_t(k,i,j) )
645 #endif
646  advcv = - ( tflx_hi(k,i,j,zdir) - tflx_hi(k-1,i ,j ,zdir) ) * rcdz(k)
647  advch = - ( ( tflx_hi(k,i,j,xdir) - tflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
648  + ( tflx_hi(k,i,j,ydir) - tflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) )
649  st(k,i,j) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) + rhot_t(k,i,j)
650 #ifdef HIST_TEND
651  if ( lhist ) then
652  advch_t(k,i,j,i_rhot) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
653  endif
654 #endif
655  enddo
656  enddo
657  enddo
658 #ifdef DEBUG
659  k = iundef; i = iundef; j = iundef
660 #endif
661  profile_stop("hevi_st")
662 
663  ! implicit solver
664 
665  profile_start("hevi_solver")
666 
667 !OCL INDEPENDENT
668 !OCL PREFETCH_SEQUENTIAL(SOFT)
669 #ifndef __GFORTRAN__
670  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
671  !$omp private(i,j,A,B,C,F1,F2,F3,PT,pg,advcv) &
672 #ifdef HIST_TEND
673  !$omp shared(lhist,pg_t,advcv_t) &
674 #endif
675 #ifdef DEBUG
676  !$omp shared(RHOT) &
677 #endif
678  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
679  !$omp shared(mflx_hi,tflx_hi,MOMZ_RK,MOMZ0) &
680  !$omp shared(DENS_RK,RHOT_RK,DENS0,RHOT0,DENS,MOMZ,POTT,DPRES) &
681  !$omp shared(GRAV,dtrk,REF_dens,Sr,Sw,St,RT2P) &
682  !$omp shared(ATMOS_DYN_FVM_flux_valueW_Z) &
683  !$omp shared(MAPF,GSQRT,J33G,I_XY,I_XYZ,I_XYW,CDZ,RCDZ,RFDZ)
684 #else
685  !$omp parallel do default(shared) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
686  !$omp private(A,B,C,F1,F2,F3,PT,pg,advcv)
687 #endif
688  do j = jjs, jje
689  do i = iis, iie
690 
691  call atmos_dyn_fvm_flux_valuew_z( pt(:), & ! (out)
692  momz(:,i,j), pott(:,i,j), gsqrt(:,i,j,i_xyz), & ! (in)
693  cdz )
694 
695  do k = ks, ke
696  a(k) = dtrk**2 * j33g * rcdz(k) * rt2p(k,i,j) * j33g / gsqrt(k,i,j,i_xyz)
697  enddo
698  b = grav * dtrk**2 * j33g / ( cdz(ks+1) + cdz(ks) )
699  f1(ks) = - ( pt(ks+1) * rfdz(ks) * a(ks+1) + b ) / gsqrt(ks,i,j,i_xyw)
700  f2(ks) = 1.0_rp + ( pt(ks ) * rfdz(ks) * ( a(ks+1)+a(ks) ) ) / gsqrt(ks,i,j,i_xyw)
701  do k = ks+1, ke-2
702  b = grav * dtrk**2 * j33g / ( cdz(k+1) + cdz(k) )
703  f1(k) = - ( pt(k+1) * rfdz(k) * a(k+1) + b ) / gsqrt(k,i,j,i_xyw)
704  f2(k) = 1.0_rp + ( pt(k ) * rfdz(k) * ( a(k+1)+a(k) ) ) / gsqrt(k,i,j,i_xyw)
705  f3(k) = - ( pt(k-1) * rfdz(k) * a(k) - b ) / gsqrt(k,i,j,i_xyw)
706  enddo
707  b = grav * dtrk**2 * j33g / ( cdz(ke) + cdz(ke-1) )
708  f2(ke-1) = 1.0_rp + ( pt(ke-1) * rfdz(ke-1) * ( a(ke)+a(ke-1) ) ) / gsqrt(ke-1,i,j,i_xyw)
709  f3(ke-1) = - ( pt(ke-2) * rfdz(ke-1) * a(ke-1) - b ) / gsqrt(ke-1,i,j,i_xyw)
710  do k = ks, ke-1
711  ! use not density at the half level but mean density between CZ(k) and C(k+1)
712  pg = - ( dpres(k+1,i,j) + rt2p(k+1,i,j)*dtrk*st(k+1,i,j) &
713  - dpres(k ,i,j) - rt2p(k ,i,j)*dtrk*st(k ,i,j) ) &
714  * rfdz(k) * j33g / gsqrt(k,i,j,i_xyw) &
715  - grav &
716  * ( f2h(k,1,i_xyz) * ( dens(k+1,i,j) - ref_dens(k+1,i,j) + sr(k+1,i,j) * dtrk ) &
717  + f2h(k,2,i_xyz) * ( dens(k ,i,j) - ref_dens(k ,i,j) + sr(k ,i,j) * dtrk ) )
718  c(k-ks+1) = momz(k,i,j) + dtrk * ( pg + sw(k,i,j) )
719 #ifdef HIST_TEND
720  if ( lhist ) pg_t(k,i,j,1) = pg
721 #endif
722  enddo
723 
724  call solve_direct( &
725  c(:), & ! (inout)
726  f1(:), f2(:), f3(:) ) ! (in)
727 
728  do k = ks, ke-1
729 #ifdef DEBUG_HEVI2HEVE
730  ! for debug (change to explicit integration)
731  c(k-ks+1) = momz(k,i,j)
732  mflx_hi(k,i,j,zdir) = mflx_hi(k,i,j,zdir) &
733  + j33g * momz(k,i,j) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
734  tflx_hi(k,i,j,zdir) = tflx_hi(k,i,j,zdir) &
735  + j33g * momz(k,i,j) * pt(k) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
736  ! use not density at the half level but mean density between CZ(k) and C(k+1)
737  momz_rk(k,i,j) = momz0(k,i,j) &
738  + dtrk*( &
739  - j33g * ( dpres(k+1,i,j)-dpres(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,i_xyw) &
740  - grav * ( f2h(k,2,i_xyz)*(dens(k,i,j)-ref_dens(k,i,j))+f2h(k,1,i_xyz)*(dens(k+1,i,j)-ref_dens(k+1,i,j)) ) &
741  + sw(k,i,j) )
742 #else
743  ! z-flux
744  mflx_hi(k,i,j,zdir) = mflx_hi(k,i,j,zdir) &
745  + j33g * c(k-ks+1) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
746  tflx_hi(k,i,j,zdir) = tflx_hi(k,i,j,zdir) &
747  + j33g * c(k-ks+1) * pt(k) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
748  ! z-momentum
749  momz_rk(k,i,j) = momz0(k,i,j) &
750  + ( c(k-ks+1) - momz(k,i,j) )
751 #endif
752  enddo
753  momz_rk(ks-1,i,j) = 0.0_rp
754  momz_rk(ke ,i,j) = 0.0_rp
755 
756  ! density and rho*theta
757  advcv = - c(1) * j33g * rcdz(ks) / gsqrt(ks,i,j,i_xyz) ! C(0) = 0
758  dens_rk(ks,i,j) = dens0(ks,i,j) + dtrk * ( advcv + sr(ks,i,j) )
759 #ifdef HIST_TEND
760  if ( lhist ) advcv_t(ks,i,j,i_dens) = advcv
761 #endif
762  advcv = - c(1) * pt(ks) * j33g * rcdz(ks) / gsqrt(ks,i,j,i_xyz) ! C(0) = 0
763  rhot_rk(ks,i,j) = rhot0(ks,i,j) + dtrk * ( advcv + st(ks,i,j) )
764 #ifdef HIST_TEND
765  if ( lhist ) advcv_t(ks,i,j,i_rhot) = advcv
766 #endif
767  do k = ks+1, ke-1
768  advcv = - ( c(k-ks+1) - c(k-ks) ) &
769  * j33g * rcdz(k) / gsqrt(k,i,j,i_xyz)
770  dens_rk(k,i,j) = dens0(k,i,j) + dtrk * ( advcv + sr(k,i,j) )
771 #ifdef HIST_TEND
772  if ( lhist ) advcv_t(k,i,j,i_dens) = advcv
773 #endif
774  advcv = - ( c(k-ks+1) * pt(k) - c(k-ks) * pt(k-1) ) &
775  * j33g * rcdz(k) / gsqrt(k,i,j,i_xyz)
776  rhot_rk(k,i,j) = rhot0(k,i,j) + dtrk * ( advcv + st(k,i,j) )
777 #ifdef HIST_TEND
778  if ( lhist ) advcv_t(k,i,j,i_rhot) = advcv
779 #endif
780  enddo
781  advcv = c(ke-ks) * j33g * rcdz(ke) / gsqrt(ke,i,j,i_xyz) ! C(KE-KS+1) = 0
782  dens_rk(ke,i,j) = dens0(ke,i,j) + dtrk * ( advcv + sr(ke,i,j) )
783 #ifdef HIST_TEND
784  if ( lhist ) advcv_t(ke,i,j,i_dens) = advcv
785 #endif
786  advcv = c(ke-ks) * pt(ke-1) * j33g * rcdz(ke) / gsqrt(ke,i,j,i_xyz) ! C(KE-KS+1) = 0
787  rhot_rk(ke,i,j) = rhot0(ke,i,j) + dtrk * ( advcv + st(ke,i,j) )
788 #ifdef HIST_TEND
789  if ( lhist ) advcv_t(ke,i,j,i_rhot) = advcv
790 #endif
791 
792 #ifdef DEBUG
793  call check_equation( &
794  c(:), &
795  dens(:,i,j), momz(:,i,j), rhot(:,i,j), dpres(:,i,j), &
796  ref_dens(:,i,j), &
797  sr(:,i,j), sw(:,i,j), st(:,i,j), &
798  j33g, gsqrt(:,i,j,:), &
799  rt2p(:,i,j), &
800  dtrk, i, j )
801 #endif
802 
803  enddo
804  enddo
805 #ifdef DEBUG
806  k = iundef; i = iundef; j = iundef
807 #endif
808 
809  profile_stop("hevi_solver")
810 
811 
812  !##### momentum equation (x) #####
813 
814  profile_start("hevi_momx")
815  ! at (u, y, w)
816  call atmos_dyn_fvm_fluxz_uyz( qflx_hi(:,:,:,zdir), & ! (out)
817  momz, momx, dens, & ! (in)
818  gsqrt(:,:,:,i_uyw), j33g, & ! (in)
819  num_diff(:,:,:,i_momx,zdir), & ! (in)
820  cdz, & ! (in)
821  iis, iie, jjs, jje ) ! (in)
822  call atmos_dyn_fvm_fluxj13_uyz( qflx_j13, & ! (out)
823  momx, momx, dens, & ! (in)
824  gsqrt(:,:,:,i_uyz), j13g(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
825  cdz, & ! (in)
826  iis, iie, jjs, jje ) ! (in)
827  call atmos_dyn_fvm_fluxj23_uyz( qflx_j23, & ! (out)
828  momy, momx, dens, & ! (in)
829  gsqrt(:,:,:,i_uyz), j23g(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
830  cdz, & ! (in)
831  iis, iie, jjs, jje ) ! (in)
832 
833  ! at (x, y, z)
834  ! note that x-index is added by -1
835  call atmos_dyn_fvm_fluxx_uyz( qflx_hi(:,:,:,xdir), & ! (out)
836  momx, momx, dens, & ! (in)
837  gsqrt(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
838  num_diff(:,:,:,i_momx,xdir), & ! (in)
839  cdz, & ! (in)
840  iis, iie, jjs, jje ) ! (in)
841 
842  ! at (u, v, z)
843  call atmos_dyn_fvm_fluxy_uyz( qflx_hi(:,:,:,ydir), & ! (out)
844  momy, momx, dens, & ! (in)
845  gsqrt(:,:,:,i_uvz), mapf(:,:,1,i_uv), & ! (in)
846  num_diff(:,:,:,i_momx,ydir), & ! (in)
847  cdz, & ! (in)
848  iis, iie, jjs, jje ) ! (in)
849 
850  !--- update momentum(x)
851  iee = min(iie,ieh)
852  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
853  !$omp private(i,j,k,advch,advcv,pg,cf,div) &
854 #ifdef HIST_TEND
855  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
856 #endif
857  !$omp shared(JJS,JJE,IIS,iee,KS,KE) &
858  !$omp shared(MOMX_RK,DPRES,DENS,MOMX,MOMY,DDIV,MOMX0,MOMX_t) &
859  !$omp shared(qflx_hi,qflx_J13,qflx_J23) &
860  !$omp shared(RCDZ,RCDY,RFDX,CDZ,FDX) &
861  !$omp shared(MAPF,GSQRT,J13G,I_XY,I_UY,I_UV,I_XYZ,I_UYW,I_UYZ) &
862  !$omp shared(dtrk,CORIOLI,divdmp_coef)
863  do j = jjs, jje
864  do i = iis, iee
865  do k = ks, ke
866 #ifdef DEBUG
867  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
868  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
869  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
870  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
871  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
872  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
873  call check( __line__, dpres(k,i+1,j) )
874  call check( __line__, dpres(k,i ,j) )
875  call check( __line__, corioli(1,i ,j) )
876  call check( __line__, corioli(1,i+1,j) )
877  call check( __line__, momy(k,i ,j ) )
878  call check( __line__, momy(k,i+1,j ) )
879  call check( __line__, momy(k,i ,j-1) )
880  call check( __line__, momy(k,i+1,j-1) )
881  call check( __line__, ddiv(k,i+1,j) )
882  call check( __line__, ddiv(k,i ,j) )
883  call check( __line__, momx0(k,i,j) )
884 #endif
885  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
886  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
887  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rcdz(k)
888  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rfdx(i) &
889  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) ) &
890  * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy)
891  pg = ( ( gsqrt(k,i+1,j,i_xyz) * dpres(k,i+1,j) & ! [x,y,z]
892  - gsqrt(k,i ,j,i_xyz) * dpres(k,i ,j) & ! [x,y,z]
893  ) * rfdx(i) &
894  + ( j13g(k ,i,j,i_uyw) &
895  * 0.5_rp * ( f2h(k,1,i_uyz) * ( dpres(k+1,i+1,j)+dpres(k+1,i,j) ) &
896  + f2h(k,2,i_uyz) * ( dpres(k ,i+1,j)+dpres(k ,i,j) ) ) & ! [x,y,z->u,y,w]
897  - j13g(k-1,i,j,i_uyw) &
898  * 0.5_rp * ( f2h(k,1,i_uyz) * ( dpres(k ,i+1,j)+dpres(k ,i,j) ) &
899  + f2h(k,2,i_uyz) * ( dpres(k-1,i+1,j)+dpres(k-1,i,j) ) ) & ! [x,y,z->u,y,w]
900  ) * rcdz(k) ) &
901  * mapf(i,j,1,i_uy)
902  cf = 0.125_rp * ( corioli(1,i+1,j )+corioli(1,i,j ) ) & ! [x,y,z->u,y,z]
903  * ( momy(k,i+1,j )+momy(k,i,j ) &
904  + momy(k,i+1,j-1)+momy(k,i,j-1) ) & ! [x,v,z->u,y,z]
905  + 0.25_rp * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy) &
906  * ( momy(k,i,j) + momy(k,i,j-1) + momy(k,i+1,j) + momy(k,i+1,j-1) ) &
907  * ( ( momy(k,i,j) + momy(k,i,j-1) + momy(k,i+1,j) + momy(k,i+1,j-1) ) * 0.25_rp &
908  * ( 1.0_rp/mapf(i+1,j,2,i_xy) - 1.0_rp/mapf(i,j,2,i_xy) ) * rfdx(i) &
909  - momx(k,i,j) &
910  * ( 1.0_rp/mapf(i,j,1,i_uv) - 1.0_rp/mapf(i,j-1,1,i_uv) ) * rcdy(j) ) &
911  * 2.0_rp / ( dens(k,i+1,j) + dens(k,i,j) ) ! metric term
912  div = divdmp_coef / dtrk * ( ddiv(k,i+1,j)/mapf(i+1,j,2,i_xy) - ddiv(k,i,j)/mapf(i,j,1,i_xy) ) &
913  * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy) * fdx(i) ! divergence damping
914  momx_rk(k,i,j) = momx0(k,i,j) &
915  + dtrk * ( ( advcv + advch - pg ) / gsqrt(k,i,j,i_uyz) + cf + div + momx_t(k,i,j) )
916 #ifdef HIST_TEND
917  if ( lhist ) then
918  advcv_t(k,i,j,i_momx) = advcv / gsqrt(k,i,j,i_uyz)
919  advch_t(k,i,j,i_momx) = advch / gsqrt(k,i,j,i_uyz)
920  pg_t(k,i,j,2) = - pg / gsqrt(k,i,j,i_uyz)
921  cf_t(k,i,j,1) = cf
922  ddiv_t(k,i,j,2) = div
923  endif
924 #endif
925  enddo
926  enddo
927  enddo
928  profile_stop("hevi_momx")
929 #ifdef DEBUG
930  k = iundef; i = iundef; j = iundef
931 #endif
932 
933  !##### momentum equation (y) #####
934  profile_start("hevi_momy")
935  ! at (x, v, w)
936  call atmos_dyn_fvm_fluxz_xvz( qflx_hi(:,:,:,zdir), & ! (out)
937  momz, momy, dens, & ! (in)
938  gsqrt(:,:,:,i_xvw), j33g, & ! (in)
939  num_diff(:,:,:,i_momy,zdir), & ! (in)
940  cdz, & ! (in)
941  iis, iie, jjs, jje ) ! (in)
942  call atmos_dyn_fvm_fluxj13_xvz( qflx_j13, & ! (out)
943  momx, momy, dens, & ! (in)
944  gsqrt(:,:,:,i_xvz), j13g(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
945  cdz, & ! (in)
946  iis, iie, jjs, jje ) ! (in)
947  call atmos_dyn_fvm_fluxj23_xvz( qflx_j23, & ! (out)
948  momy, momy, dens, & ! (in)
949  gsqrt(:,:,:,i_xvz), j23g(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
950  cdz, & ! (in)
951  iis, iie, jjs, jje ) ! (in)
952 
953  ! at (u, v, z)
954  call atmos_dyn_fvm_fluxx_xvz( qflx_hi(:,:,:,xdir), & ! (out)
955  momx, momy, dens, & ! (in)
956  gsqrt(:,:,:,i_uvz), mapf(:,:,:,i_uv), & ! (in)
957  num_diff(:,:,:,i_momy,xdir), & ! (in)
958  cdz, & ! (in)
959  iis, iie, jjs, jje ) ! (in)
960 
961  ! at (x, y, z)
962  ! note that y-index is added by -1
963  call atmos_dyn_fvm_fluxy_xvz( qflx_hi(:,:,:,ydir), & ! (out)
964  momy, momy, dens, & ! (in)
965  gsqrt(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
966  num_diff(:,:,:,i_momy,ydir), & ! (in
967  cdz, & ! (in)
968  iis, iie, jjs, jje ) ! (in)
969 
970  !--- update momentum(y)
971  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
972  !$omp private(i,j,k,advch,advcv,pg,cf,div) &
973 #ifdef HIST_TEND
974  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
975 #endif
976  !$omp shared(JJS,JJE,JEH,IIS,IIE,KS,KE) &
977  !$omp shared(MOMY_RK,DPRES,DENS,MOMX,MOMY,DDIV,MOMY0,MOMY_t) &
978  !$omp shared(qflx_hi,qflx_J13,qflx_J23) &
979  !$omp shared(RCDZ,RCDX,RFDY,CDZ,FDY) &
980  !$omp shared(MAPF,GSQRT,J23G,I_XY,I_XV,I_UV,I_XYZ,I_XVW,I_XVZ) &
981  !$omp shared(dtrk,CORIOLI,divdmp_coef)
982  do j = jjs, min(jje,jeh)
983  do i = iis, iie
984  do k = ks, ke
985 #ifdef DEBUG
986  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
987  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
988  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
989  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
990  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
991  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
992  call check( __line__, dpres(k,i,j ) )
993  call check( __line__, dpres(k,i,j+1) )
994  call check( __line__, corioli(1,i,j ) )
995  call check( __line__, corioli(1,i,j+1) )
996  call check( __line__, momx(k,i ,j ) )
997  call check( __line__, momx(k,i ,j+1) )
998  call check( __line__, momx(k,i-1,j ) )
999  call check( __line__, momx(k,i-1,j+1) )
1000  call check( __line__, ddiv(k,i,j+1) )
1001  call check( __line__, ddiv(k,i,j ) )
1002  call check( __line__, momy_t(k,i,j) )
1003  call check( __line__, momy0(k,i,j) )
1004 #endif
1005  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
1006  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
1007  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rcdz(k)
1008  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
1009  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rfdy(j) ) &
1010  * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv)
1011  pg = ( ( gsqrt(k,i,j+1,i_xyz) * dpres(k,i,j+1) & ! [x,y,z]
1012  - gsqrt(k,i,j ,i_xyz) * dpres(k,i,j ) & ! [x,y,z]
1013  ) * rfdy(j) &
1014  + ( j23g(k ,i,j,i_xvw) &
1015  * 0.5_rp * ( f2h(k ,1,i_xvz) * ( dpres(k+1,i,j+1)+dpres(k+1,i,j) ) &
1016  + f2h(k ,2,i_xvz) * ( dpres(k ,i,j+1)+dpres(k ,i,j) ) ) & ! [x,y,z->x,v,w]
1017  - j23g(k-1,i,j,i_xvw) &
1018  * 0.5_rp * ( f2h(k-1,1,i_xvz) * ( dpres(k ,i,j+1)+dpres(k ,i,j) ) &
1019  + f2h(k-1,2,i_xvz) * ( dpres(k-1,i,j+1)+dpres(k-1,i,j) ) ) & ! [x,y,z->x,v,w]
1020  ) * rcdz(k) ) &
1021  * mapf(i,j,2,i_xv)
1022  cf = - 0.125_rp * ( corioli(1,i ,j+1)+corioli(1,i ,j) ) & ! [x,y,z->x,v,z]
1023  * ( momx(k,i ,j+1)+momx(k,i ,j) &
1024  + momx(k,i-1,j+1)+momx(k,i-1,j) ) & ! [u,y,z->x,v,z]
1025  - 0.25_rp * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv) &
1026  * ( momx(k,i,j) + momx(k,i-1,j) + momx(k,i,j+1) + momx(k,i-1,j+1) ) &
1027  * ( momy(k,i,j) &
1028  * ( 1.0_rp/mapf(i,j,2,i_uv) - 1.0_rp/mapf(i-1,j,2,i_uv) ) * rcdx(i) &
1029  - 0.25_rp * ( momx(k,i,j)+momx(k,i-1,j)+momx(k,i,j+1)+momx(k,i-1,j+1) ) &
1030  * ( 1.0_rp/mapf(i,j+1,1,i_xy) - 1.0_rp/mapf(i,j,1,i_xy) ) * rfdy(j) ) &
1031  * 2.0_rp / ( dens(k,i,j+1) + dens(k,i,j) ) ! metoric term
1032  div = divdmp_coef / dtrk * ( ddiv(k,i,j+1)/mapf(i,j+1,1,i_xy) - ddiv(k,i,j)/mapf(i,j,1,i_xy) ) &
1033  * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv) * fdy(j) ! divergence damping
1034  momy_rk(k,i,j) = momy0(k,i,j) &
1035  + dtrk * ( ( advcv + advch - pg ) / gsqrt(k,i,j,i_xvz) + cf + div + momy_t(k,i,j) )
1036 #ifdef HIST_TEND
1037  if ( lhist ) then
1038  advcv_t(k,i,j,i_momy) = advcv / gsqrt(k,i,j,i_xvz)
1039  advch_t(k,i,j,i_momy) = advch / gsqrt(k,i,j,i_xvz)
1040  pg_t(k,i,j,3) = - pg / gsqrt(k,i,j,i_xvz)
1041  cf_t(k,i,j,2) = cf
1042  ddiv_t(k,i,j,3) = div
1043  endif
1044 #endif
1045  enddo
1046  enddo
1047  enddo
1048  profile_stop("hevi_momy")
1049 #ifdef DEBUG
1050  k = iundef; i = iundef; j = iundef
1051 #endif
1052 
1053  enddo
1054  enddo
1055 
1056 #ifdef PROFILE_FIPP
1057  call fipp_stop()
1058 #endif
1059 
1060 #ifdef HIST_TEND
1061  if ( lhist ) then
1062  call file_history_in(advcv_t(:,:,:,i_dens), 'DENS_t_advcv', 'tendency of density (vert. advection) (w/ HIST_TEND)', 'kg/m3/s' )
1063  call file_history_in(advcv_t(:,:,:,i_momz), 'MOMZ_t_advcv', 'tendency of momentum z (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1064  call file_history_in(advcv_t(:,:,:,i_momx), 'MOMX_t_advcv', 'tendency of momentum x (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1065  call file_history_in(advcv_t(:,:,:,i_momy), 'MOMY_t_advcv', 'tendency of momentum y (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1066  call file_history_in(advcv_t(:,:,:,i_rhot), 'RHOT_t_advcv', 'tendency of rho*theta (vert. advection) (w/ HIST_TEND)', 'K kg/m3/s' )
1067 
1068  call file_history_in(advch_t(:,:,:,i_dens), 'DENS_t_advch', 'tendency of density (horiz. advection) (w/ HIST_TEND)', 'kg/m3/s' )
1069  call file_history_in(advch_t(:,:,:,i_momz), 'MOMZ_t_advch', 'tendency of momentum z (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1070  call file_history_in(advch_t(:,:,:,i_momx), 'MOMX_t_advch', 'tendency of momentum x (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1071  call file_history_in(advch_t(:,:,:,i_momy), 'MOMY_t_advch', 'tendency of momentum y (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1072  call file_history_in(advch_t(:,:,:,i_rhot), 'RHOT_t_advch', 'tendency of rho*theta (horiz. advection) (w/ HIST_TEND)', 'K kg/m3/s' )
1073 
1074  call file_history_in(pg_t(:,:,:,1), 'MOMZ_t_pg', 'tendency of momentum z (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1075  call file_history_in(pg_t(:,:,:,2), 'MOMX_t_pg', 'tendency of momentum x (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1076  call file_history_in(pg_t(:,:,:,3), 'MOMY_t_pg', 'tendency of momentum y (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1077 
1078  call file_history_in(wdmp_t(:,:,:), 'MOMZ_t_wdamp', 'tendency of momentum z (Rayleigh damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1079 
1080  call file_history_in(ddiv_t(:,:,:,1), 'MOMZ_t_ddiv', 'tendency of momentum z (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1081  call file_history_in(ddiv_t(:,:,:,2), 'MOMX_t_ddiv', 'tendency of momentum x (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1082  call file_history_in(ddiv_t(:,:,:,3), 'MOMY_t_ddiv', 'tendency of momentum y (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1083 
1084  call file_history_in(cf_t(:,:,:,1), 'MOMX_t_cf', 'tendency of momentum x (coliolis force) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1085  call file_history_in(cf_t(:,:,:,2), 'MOMY_t_cf', 'tendency of momentum y (coliolis force) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1086  endif
1087 #endif
1088 
integer, parameter, public i_rhot
Definition: scale_index.F90:32
integer, public va
Definition: scale_index.F90:35
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_xvz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_uyz
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_momx
Definition: scale_index.F90:30
integer, parameter, public i_momz
Definition: scale_index.F90:29
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_xyw
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz
integer, public iblock
block size for cache blocking: x
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_xyw
integer, public ja
of whole cells: y, local, with HALO
integer, parameter, public i_dens
Definition: scale_index.F90:28
integer, parameter, public i_momy
Definition: scale_index.F90:31
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_xyw
integer, public is
start point of inner domain: x, local
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_xyw
integer, public ie
end point of inner domain: x, local
procedure(flux_z), pointer, public atmos_dyn_fvm_fluxz_xvz
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ieh
end point of inner domain: x, local (half level)
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
integer, public jblock
block size for cache blocking: y
integer, public kmax
of computational cells: z, local
procedure(valuew), pointer, public atmos_dyn_fvm_flux_valuew_z
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_xvz
procedure(flux_z), pointer, public atmos_dyn_fvm_fluxz_uyz
module scale_atmos_dyn_fvm_flux
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.
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_uyz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_xvz
integer, public jeh
end point of inner domain: y, local (half level)
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_uyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz
module file_history
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_uyz
procedure(flux_wz), pointer, public atmos_dyn_fvm_fluxz_xyw
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_xvz
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solve_direct()

subroutine scale_atmos_dyn_tstep_short_fvm_hevi::solve_direct ( real(rp), dimension(kmax-1), intent(inout)  C,
real(rp), dimension(ka), intent(in)  F1,
real(rp), dimension(ka), intent(in)  F2,
real(rp), dimension(ka), intent(in)  F3 
)

Definition at line 1095 of file scale_atmos_dyn_tstep_short_fvm_hevi.F90.

References scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdz, scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdz, scale_const::const_eps, scale_const::const_grav, scale_atmos_grid_cartesc_index::i_xyw, scale_atmos_grid_cartesc_index::i_xyz, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, and scale_prc::prc_abort().

Referenced by atmos_dyn_tstep_short_fvm_hevi().

1095  use scale_prc, only: &
1096  prc_abort
1097  implicit none
1098  real(RP), intent(inout) :: c(kmax-1)
1099  real(RP), intent(in) :: f1(ka)
1100  real(RP), intent(in) :: f2(ka)
1101  real(RP), intent(in) :: f3(ka)
1102 
1103  real(RP) :: e(kmax-2)
1104  real(RP) :: f(kmax-2)
1105 
1106  real(RP) :: rdenom
1107 
1108  integer :: k
1109 
1110  rdenom = 1.0_rp / f2(ks)
1111  e(1) = - f1(ks) * rdenom
1112  f(1) = c(1) * rdenom
1113  do k = 2, kmax-2
1114  rdenom = 1.0_rp / ( f2(k+ks-1) + f3(k+ks-1) * e(k-1) )
1115  e(k) = - f1(k+ks-1) * rdenom
1116  f(k) = ( c(k) - f3(k+ks-1) * f(k-1) ) * rdenom
1117  enddo
1118 
1119  ! C = \rho w
1120  c(kmax-1) = ( c(kmax-1) - f3(ke-1) * f(kmax-2) ) &
1121  / ( f2(ke-1) + f3(ke-1) * e(kmax-2) ) ! C(KMAX-1) = f(KMAX-1)
1122  do k = kmax-2, 1, -1
1123  c(k) = e(k) * c(k+1) + f(k)
1124  enddo
1125 
1126  return
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public ks
start point of inner domain: z, local
integer, public kmax
of computational cells: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
Here is the call graph for this function:
Here is the caller graph for this function: