SCALE-RM
mod_atmos_phy_cp_driver.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
12 !-------------------------------------------------------------------------------
13 #include "inc_openmp.h"
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23  use scale_tracer
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
33  public :: atmos_phy_cp_driver
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48 contains
49  !-----------------------------------------------------------------------------
51  subroutine atmos_phy_cp_driver_setup
52  use scale_atmos_phy_cp, only: &
54  use mod_atmos_admin, only: &
57  implicit none
58  !---------------------------------------------------------------------------
59 
60  if( io_l ) write(io_fid_log,*)
61  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_CP] / Origin[SCALE-RM]'
62 
63  if ( atmos_sw_phy_cp ) then
64 
65  ! setup library component
67 
68  else
69  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
70  endif
71 
72  return
73  end subroutine atmos_phy_cp_driver_setup
74 
75  !-----------------------------------------------------------------------------
78  use mod_atmos_admin, only: &
80  implicit none
81 
82  if ( atmos_sw_phy_cp ) then
83 
84  ! run once (only for the diagnostic value)
85  call prof_rapstart('ATM_Cumulus', 1)
86  call atmos_phy_cp_driver( update_flag = .true. )
87  call prof_rapend ('ATM_Cumulus', 1)
88 
89  endif
90 
91  return
92  end subroutine atmos_phy_cp_driver_resume
93 
94  !-----------------------------------------------------------------------------
96  subroutine atmos_phy_cp_driver( update_flag )
97  use scale_rm_statistics, only: &
99  stat_total
100  use scale_history, only: &
101  hist_in
102  use scale_atmos_hydrometeor, only: &
104  use scale_atmos_phy_cp, only: &
106  use scale_atmos_phy_mp, only: &
107  qs_mp, &
108  qe_mp
109  use mod_atmos_admin, only: &
111  use mod_atmos_vars, only: &
112  dens => dens_av, &
113  momz => momz_av, &
114  momx => momx_av, &
115  momy => momy_av, &
116  rhot => rhot_av, &
117  qtrc => qtrc_av, &
118  dens_t => dens_tp, &
119  momz_t => momz_tp, &
120  momx_t => momx_tp, &
121  momy_t => momy_tp, &
122  rhot_t => rhot_tp, &
123  rhoq_t => rhoq_tp
124  use mod_atmos_phy_cp_vars, only: &
125  dens_t_cp => atmos_phy_cp_dens_t, &
126  momz_t_cp => atmos_phy_cp_momz_t, &
127  momx_t_cp => atmos_phy_cp_momx_t, &
128  momy_t_cp => atmos_phy_cp_momy_t, &
129  rhot_t_cp => atmos_phy_cp_rhot_t, &
130  rhoq_t_cp => atmos_phy_cp_rhoq_t, &
131  mflx_cloudbase => atmos_phy_cp_mflx_cloudbase, &
132  sflx_rain => atmos_phy_cp_sflx_rain, & ! convective rain [kg/m2/s]
133  cloudtop => atmos_phy_cp_cloudtop, & ! cloud top height [m]
134  cloudbase => atmos_phy_cp_cloudbase, & ! cloud base height [m]
135  cldfrac_dp => atmos_phy_cp_cldfrac_dp, & ! cloud fraction (deep convection) (0-1)
136  cldfrac_sh => atmos_phy_cp_cldfrac_sh, & ! cloud fraction (shallow convection) (0-1)
137  kf_nca => atmos_phy_cp_kf_nca, & ! advection/cumulus convection timescale/dt for KF [step]
138  kf_w0avg => atmos_phy_cp_kf_w0avg ! rannning mean vertical wind velocity for KF [m/s]
139  implicit none
140 
141  logical, intent(in) :: update_flag
142 
143  real(RP) :: total ! dummy
144 
145  integer :: k, i, j, iq
146  !---------------------------------------------------------------------------
147 
148  if ( update_flag ) then
149 
150  call atmos_phy_cp( dens, & ! [IN]
151  momz, & ! [IN]
152  momx, & ! [IN]
153  momy, & ! [IN]
154  rhot, & ! [IN]
155  qtrc, & ! [IN]
156  dens_t_cp, & ! [INOUT]
157  momz_t_cp, & ! [INOUT]
158  momx_t_cp, & ! [INOUT]
159  momy_t_cp, & ! [INOUT]
160  rhot_t_cp, & ! [INOUT]
161  rhoq_t_cp, & ! [INOUT]
162  mflx_cloudbase, & ! [INOUT]
163  sflx_rain, & ! [OUT]
164  cloudtop, & ! [OUT]
165  cloudbase, & ! [OUT]
166  cldfrac_dp, & ! [OUT]
167  cldfrac_sh, & ! [OUT]
168  kf_nca, & ! [OUT]
169  kf_w0avg ) ! [OUT]
170 
171  ! tentative reset
172 !OCL XFILL
173  do j = js, je
174  do i = is, ie
175  do k = ks, ke
176  momz_t_cp(k,i,j) = 0.0_rp
177  momx_t_cp(k,i,j) = 0.0_rp
178  momy_t_cp(k,i,j) = 0.0_rp
179  enddo
180  enddo
181  enddo
182 
183 !OCL XFILL
184  do j = js, je
185  do i = is, ie
186  mflx_cloudbase(i,j) = 0.0_rp
187  enddo
188  enddo
189 
190  ! diagnose tendency of number concentration
191  call atmos_hydrometeor_diagnose_number_concentration( rhoq_t_cp(:,:,:,:) ) ! [INOUT]
192 
193  call hist_in( mflx_cloudbase(:,:), 'CBMFX', 'cloud base mass flux', 'kg/m2/s', nohalo=.true. )
194  call hist_in( sflx_rain(:,:), 'RAIN_CP', 'surface rain rate by CP', 'kg/m2/s', nohalo=.true. )
195  call hist_in( sflx_rain(:,:), 'PREC_CP', 'surface precipitation rate by CP', 'kg/m2/s', nohalo=.true. )
196  call hist_in( cloudtop(:,:), 'CUMHGT', 'CP cloud top height', 'm', nohalo=.true. )
197  call hist_in( cloudbase(:,:), 'CUBASE', 'CP cloud base height', 'm', nohalo=.true. )
198  call hist_in( cldfrac_dp(:,:,:), 'CUMFRC_DP', 'CP cloud fraction (deep)', '1', nohalo=.true. )
199  call hist_in( cldfrac_sh(:,:,:), 'CUMFRC_SH', 'CP cloud fraction (shallow)', '1', nohalo=.true. )
200 
201  call hist_in( kf_nca(:,:), 'kf_nca', 'advection or cumulus convection timescale for KF', 's', nohalo=.true. )
202  call hist_in( kf_w0avg(:,:,:), 'kf_w0avg', 'rannning mean vertical wind velocity for KF', 'kg/m2/s', nohalo=.true. )
203 
204  call hist_in( dens_t_cp(:,:,:), 'DENS_t_CP', 'tendency DENS in CP', 'kg/m3/s' , nohalo=.true. )
205  call hist_in( momz_t_cp(:,:,:), 'MOMZ_t_CP', 'tendency MOMZ in CP', 'kg/m2/s2' , nohalo=.true. )
206  call hist_in( momx_t_cp(:,:,:), 'MOMX_t_CP', 'tendency MOMX in CP', 'kg/m2/s2' , nohalo=.true. )
207  call hist_in( momy_t_cp(:,:,:), 'MOMY_t_CP', 'tendency MOMY in CP', 'kg/m2/s2' , nohalo=.true. )
208  call hist_in( rhot_t_cp(:,:,:), 'RHOT_t_CP', 'tendency RHOT in CP', 'K*kg/m3/s', nohalo=.true. )
209 
210  do iq = qs_mp, qe_mp
211  call hist_in( rhoq_t_cp(:,:,:,iq), trim(tracer_name(iq))//'_t_CP', &
212  'tendency rho*'//trim(tracer_name(iq))//'in CP', 'kg/m3/s', nohalo=.true. )
213  enddo
214 
215  endif
216 
217  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
218  do j = js, je
219  do i = is, ie
220  do k = ks, ke
221  dens_t(k,i,j) = dens_t(k,i,j) + dens_t_cp(k,i,j)
222  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_cp(k,i,j)
223  momx_t(k,i,j) = momx_t(k,i,j) + momx_t_cp(k,i,j)
224  momy_t(k,i,j) = momy_t(k,i,j) + momy_t_cp(k,i,j)
225  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_cp(k,i,j)
226  enddo
227  enddo
228  enddo
229 
230  do iq = qs_mp, qe_mp
231  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
232  do j = js, je
233  do i = is, ie
234  do k = ks, ke
235  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_cp(k,i,j,iq)
236  enddo
237  enddo
238  enddo
239  enddo
240 
241  if ( statistics_checktotal ) then
242  call stat_total( total, dens_t_cp(:,:,:), 'DENS_t_CP' )
243  call stat_total( total, momz_t_cp(:,:,:), 'MOMZ_t_CP' )
244  call stat_total( total, momx_t_cp(:,:,:), 'MOMX_t_CP' )
245  call stat_total( total, momy_t_cp(:,:,:), 'MOMY_t_CP' )
246  call stat_total( total, rhot_t_cp(:,:,:), 'RHOT_t_CP' )
247 
248  do iq = qs_mp, qe_mp
249  call stat_total( total, rhoq_t_cp(:,:,:,iq), trim(tracer_name(iq))//'_t_CP' )
250  enddo
251  endif
252 
253  return
254  end subroutine atmos_phy_cp_driver
255 
256 end module mod_atmos_phy_cp_driver
module ATMOS admin
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, public is
start point of inner domain: x, local
logical, public atmos_sw_phy_cp
subroutine, public atmos_phy_cp_driver_resume
Redume.
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:,:), allocatable, target, public momz
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_momy_t
integer, public je
end point of inner domain: y, local
module Atmosphere / Physics Cumulus
subroutine, public atmos_phy_cp_setup(CP_TYPE)
Setup Cumulus parameterization.
module ATMOSPHERE / Physics Cumulus Parameterization
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_momx_t
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public momy_tp
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module ATMOSPHERE / Physics Cumulus
module ATMOSPHERE / Physics Cloud Microphysics
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
real(rp), dimension(:,:,:), allocatable, target, public momx
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public rhot_tp
character(len=h_short), public atmos_phy_cp_type
real(rp), dimension(:,:,:), allocatable, target, public dens
subroutine, public atmos_hydrometeor_diagnose_number_concentration(QTRC)
character(len=h_short), dimension(qa_max), public tracer_name
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_dens_t
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_cp_rhoq_t
module Statistics
module grid index
real(rp), dimension(:,:,:), pointer, public momx_av
module TRACER
subroutine, public atmos_phy_cp_driver_setup
Setup.
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
integer, public js
start point of inner domain: y, local
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_momz_t
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhot_t
subroutine, public atmos_phy_cp_driver(update_flag)
Driver.
character(len=h_short), public atmos_phy_mp_type
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public momx_tp
real(rp), dimension(:,:,:), allocatable, target, public momy
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:156
procedure(cp), pointer, public atmos_phy_cp
real(rp), dimension(:,:,:), allocatable, public momz_tp
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
module HISTORY
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:,:,:), pointer, public rhot_av
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_mflx_cloudbase
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(:,:,:), pointer, public momy_av
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:204
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_rain
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc