SCALE-RM
Functions/Subroutines
scale_atmos_phy_mp_convert Module Reference

module ATMOSPHERE / Physics Cloud Microphysics - Convert More...

Functions/Subroutines

subroutine, public atmos_phy_mp_bulk2bin (xs, xe, ys, ye, dims, it, rank, handle, basename_org, dens_org, qtrc_org)
 Bulk to Bin. More...
 
subroutine, public atmos_phy_mp_bin2bulk
 Bin to Bulk. More...
 

Detailed Description

module ATMOSPHERE / Physics Cloud Microphysics - Convert

Description
Convert module for Cloud Microphysics Bulk to Bin, and Bin to Bulk
Author
Team SCALE
History
  • 2015-01-20 (Y.Sato) [new]
NAMELIST
  • PARAM_ATMOS_PHY_MP_BIN2BULK
    nametypedefault valuecomment
    SIGMA_SDF real(RP)
    R0_SDF real(RP)
    N0_SDF real(RP)
    RHO_SDF real(RP)
    MP_TYPE_OUTER character(LEN=16) "NONE"
    MP_TYPE_INNER character(LEN=16) "NONE"

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_mp_bulk2bin()

subroutine, public scale_atmos_phy_mp_convert::atmos_phy_mp_bulk2bin ( integer, intent(in)  xs,
integer, intent(in)  xe,
integer, intent(in)  ys,
integer, intent(in)  ye,
integer, dimension(:), intent(in)  dims,
integer, intent(in)  it,
integer, intent(in)  rank,
integer, intent(in)  handle,
character(len=*), intent(in)  basename_org,
real(rp), dimension(:,:,:), intent(in)  dens_org,
real(rp), dimension(:,:,:,:), intent(inout)  qtrc_org 
)

Bulk to Bin.

Definition at line 88 of file scale_atmos_phy_mp_convert.F90.

References scale_tracer::aq_name, scale_tracer::i_qc, scale_tracer::i_qg, scale_tracer::i_qi, scale_tracer::i_qr, scale_tracer::i_qs, scale_tracer::i_qv, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, dc_log::log(), scale_grid_nest::parent_imax, scale_grid_nest::parent_jmax, scale_process::prc_mpistop(), scale_tracer::qa, and scale_specfunc::sf_gamma().

Referenced by mod_realinput_scale::parentatominputscale().

88  use scale_grid_nest, only: &
89  parent_imax, &
91  use scale_specfunc, only: &
92  sf_gamma
93  use gtool_file, only: &
94  fileread
95  implicit none
96 
97  integer, intent(in) :: xs, xe
98  integer, intent(in) :: ys, ye
99  integer, intent(in) :: dims(:)
100  integer, intent(in) :: it
101  integer, intent(in) :: rank
102  integer, intent(in) :: handle
103  character(LEN=*), intent(in) :: basename_org
104  real(RP), intent(in) :: dens_org(:,:,:)
105  real(RP), intent(inout) :: qtrc_org(:,:,:,:)
106 ! real(RP), intent(in) :: dens_org(dims(1)+2,dims(2),dims(3))
107 ! real(RP), intent(inout) :: qtrc_org(dims(1)+2,dims(2),dims(3),11)
108 
109  real(RP) :: sigma_sdf(4)
110  real(RP) :: r0_sdf(4)
111  real(RP) :: n0_sdf(4)
112  real(RP) :: rho_sdf(4)
113  character(LEN=11),parameter :: fname_micpara="micpara.dat"
114  character(LEN=16) :: mp_type_outer = "NONE"
115  character(LEN=16) :: mp_type_inner = "NONE"
116  integer :: comm_world
117 
118  namelist / param_atmos_phy_mp_bin2bulk / &
119  sigma_sdf, &
120  r0_sdf, &
121  n0_sdf, &
122  rho_sdf, &
123  mp_type_outer, &
124  mp_type_inner
125 
126  real(RP), allocatable :: read3di(:,:,:)
127  real(RP), allocatable :: qtrc_tmp(:,:,:,:)
128 
129  integer :: fid_micpara
130  real(RP) :: coef0, coef1, coef2
131  real(RP) :: tmp_hyd, n_hyd, lambda_hyd
132  real(RP) :: dummy( nbin ), radc( nbin )
133  integer :: k, i, j, iq, ierr
134  integer :: nnbin, nnspc, nn
135  !---------------------------------------------------------------------------
136 
137  !--- define coefficients
138  coef0 = 4.0_rp/3.0_rp*pi
139  coef1 = 4.0_rp/3.0_rp*sqrt(pi/2.0_rp)
140 
141  !--- determine the parameters for interpolating SDF from qxx, Nxx of parent domain
142  sigma_sdf(1) = 0.2_rp
143  sigma_sdf(2) = 0.35_rp
144  sigma_sdf(3) = 0.35_rp
145  sigma_sdf(4) = 0.35_rp
146  r0_sdf(1) = 5.e-6_rp
147  r0_sdf(2) = 2.61e-6_rp
148  r0_sdf(3) = 5.e-6_rp
149  r0_sdf(4) = 2.61e-6_rp
150  n0_sdf(1) = 8.0e+6_rp
151  n0_sdf(2) = 0.0_rp
152  n0_sdf(3) = 3.0e+6_rp
153  n0_sdf(4) = 4.0e+6_rp
154  rho_sdf(1) = dwatr
155  rho_sdf(2) = dice
156  rho_sdf(3) = 100.0_rp
157  rho_sdf(4) = 400.0_rp
158 
159  !---- initiate
160  allocate( read3di( parent_imax(handle), parent_jmax(handle), dims(1) ) )
161  allocate( qtrc_tmp(dims(1)+2,dims(2),dims(3),11) )
162 
163  rewind(io_fid_conf)
164  read(io_fid_conf,nml=param_atmos_phy_mp_bin2bulk,iostat=ierr)
165  if( ierr < 0 ) then !--- missing
166  if( io_l ) write(io_fid_log,*) 'xxx Not found namelist. Check!'
167  call prc_mpistop
168  elseif( ierr > 0 ) then !--- fatal error
169  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_MP_BIN2BULK. Check!'
170  call prc_mpistop
171  endif
172  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_mp_bin2bulk)
173 
174  if( trim(mp_type_inner) /= "SUZUKI10" ) then
175  write(*,*) 'xxx MP_TYPE_INNER should be SUZUKI10 Check!'
176  write(*,*) 'Now MP_TYPE_INNER set as ', mp_type_inner
177  call prc_mpistop
178  endif
179 
180  !--- read micpara.dat (microphysical parameter) and broad cast
181  fid_micpara = io_get_available_fid()
182  !--- open parameter of cloud microphysics
183  open ( fid_micpara, file = fname_micpara, form = 'formatted', status = 'old', iostat=ierr )
184 
185  !--- micpara.dat does not exist
186  if( ierr == 0 ) then
187 
188  read( fid_micpara,* ) nnspc, nnbin
189 
190  if( nnbin /= nbin ) then
191  write(*,*) 'xxx nbin in inc_tracer and nbin in micpara.dat is different check!'
192  call prc_mpistop
193  end if
194 
195  ! grid parameter
196  if( io_l ) write(io_fid_log,*) '*** Radius of cloud ****'
197  do iq = 1, nbin
198  read( fid_micpara,* ) nn, dummy( iq ), radc( iq )
199  if( io_l ) write(io_fid_log,'(a,1x,i3,1x,a,1x,e15.7,1x,a)') &
200  "Radius of ", iq, "th cloud bin (bin center)= ", radc( iq ) , "[m]"
201  end do
202 
203  close ( fid_micpara )
204 
205  else
206 
207  write(*,*) 'xxx micpara.dat does not exist. check!'
208  call prc_mpistop
209 
210  endif
211 
212  call fileread( read3di(:,:,:), basename_org, "QV", it, rank )
213  do k = 1, dims(1)
214  qtrc_org(k+2,xs:xe,ys:ye,i_qv) = read3di(:,:,k)
215  end do
216 
217  if( trim(mp_type_outer) == "KESSLER" ) then
218 
219  if( io_l ) write(io_fid_log,*)
220  if( io_l ) write(io_fid_log,*) '+++ SDF of Bin model is created from'
221  if( io_l ) write(io_fid_log,*) '+++ Kessler type Bulk microphysical model'
222 
223  call fileread( read3di(:,:,:), basename_org, "QC", it, rank )
224  do k = 1, dims(1)
225  qtrc_tmp(k+2,xs:xe,ys:ye,i_qc) = read3di(:,:,k)
226  end do
227 ! do k = 1, dims(1)
228 ! do i = xs, xe
229 ! do j = ys, ye
230 ! qtrc_tmp(k+2,i,j,I_QC) = qtrc_tmp(k+2,i,j,I_QC)*dens_org(k+2,i,j)
231 ! enddo
232 ! enddo
233 ! enddo
234 
235  call fileread( read3di(:,:,:), basename_org, "QR", it, rank )
236  do k = 1, dims(1)
237  qtrc_tmp(k+2,xs:xe,ys:ye,i_qr) = read3di(:,:,k)
238  end do
239 ! do k = 1, dims(1)
240 ! do i = xs, xe
241 ! do j = ys, ye
242 ! qtrc_tmp(k+2,i,j,I_QR) = qtrc_tmp(k+2,i,j,I_QR)*dens_org(k+2,i,j)
243 ! enddo
244 ! enddo
245 ! enddo
246 
247  do k = 1, dims(1)
248  do i = xs, xe
249  do j = ys, ye
250 
251  do iq = 1, nbin
252  dummy(iq) = coef1 / sigma_sdf(1) * dwatr * radc( iq )**3 &
253  * exp( &
254  - ( log( radc(iq) )-log( r0_sdf(1) ) )**2*0.5_rp &
255  / sigma_sdf(1) / sigma_sdf(1) &
256  )
257  enddo
258 
259  tmp_hyd = 0.0_rp
260  do iq = 1, nbin
261  tmp_hyd = tmp_hyd + dummy(iq)
262  enddo
263 
264  coef2 = ( qtrc_tmp(k+2,i,j,i_qc)+qtrc_tmp(k+2,i,j,i_qr) ) &
265  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
266  do iq = 1, nbin
267  qtrc_org(k+2,i,j,i_qv+iq) = coef2 * dummy(iq)
268  enddo
269 
270  enddo
271  enddo
272  enddo
273 
274  elseif( trim(mp_type_outer) == "TOMITA08" ) then
275 
276  if( io_l ) write(io_fid_log,*)
277  if( io_l ) write(io_fid_log,*) '+++ SDF of Bin model is created from'
278  if( io_l ) write(io_fid_log,*) '+++ TOMITA08 Bulk microphysical model'
279 
280  call fileread( read3di(:,:,:), basename_org, "QC", it, rank )
281  do k = 1, dims(1)
282  qtrc_tmp(k+2,xs:xe,ys:ye,i_qc) = read3di(:,:,k)
283  end do
284 ! do k = 1, dims(1)
285 ! do i = xs, xe
286 ! do j = ys, ye
287 ! qtrc_tmp(k+2,i,j,I_QC) = qtrc_tmp(k+2,i,j,I_QC)*dens_org(k+2,i,j)
288 ! enddo
289 ! enddo
290 ! enddo
291 
292  call fileread( read3di(:,:,:), basename_org, "QR", it, rank )
293  do k = 1, dims(1)
294  qtrc_tmp(k+2,xs:xe,ys:ye,i_qr) = read3di(:,:,k)
295  end do
296 ! do k = 1, dims(1)
297 ! do i = xs, xe
298 ! do j = ys, ye
299 ! qtrc_tmp(k+2,i,j,I_QR) = qtrc_tmp(k+2,i,j,I_QR)*dens_org(k+2,i,j)
300 ! enddo
301 ! enddo
302 ! enddo
303 
304  call fileread( read3di(:,:,:), basename_org, "QI", it, rank )
305  do k = 1, dims(1)
306  qtrc_tmp(k+2,xs:xe,ys:ye,i_qi) = read3di(:,:,k)
307  end do
308 ! do k = 1, dims(1)
309 ! do i = xs, xe
310 ! do j = ys, ye
311 ! qtrc_tmp(k+2,i,j,I_QI) = qtrc_tmp(k+2,i,j,I_QI)*dens_org(k+2,i,j)
312 ! enddo
313 ! enddo
314 ! enddo
315 
316  call fileread( read3di(:,:,:), basename_org, "QS", it, rank )
317  do k = 1, dims(1)
318  qtrc_tmp(k+2,xs:xe,ys:ye,i_qs) = read3di(:,:,k)
319  end do
320 ! do k = 1, dims(1)
321 ! do i = xs, xe
322 ! do j = ys, ye
323 ! qtrc_tmp(k+2,i,j,I_QS) = qtrc_tmp(k+2,i,j,I_QS)*dens_org(k+2,i,j)
324 ! enddo
325 ! enddo
326 ! enddo
327 
328  call fileread( read3di(:,:,:), basename_org, "QG", it, rank )
329  do k = 1, dims(1)
330  qtrc_tmp(k+2,xs:xe,ys:ye,i_qg) = read3di(:,:,k)
331  end do
332 ! do k = 1, dims(1)
333 ! do i = xs, xe
334 ! do j = ys, ye
335 ! qtrc_tmp(k+2,i,j,I_QG) = qtrc_tmp(k+2,i,j,I_QG)*dens_org(k+2,i,j)
336 ! enddo
337 ! enddo
338 ! enddo
339 
340  if( nspc == 1 ) then !--- put all hydrometeors to liquid (warm bin)
341 
342  do k = 1, dims(1)
343  do i = xs, xe
344  do j = ys, ye
345 
346  do iq = 1, nbin
347  dummy(iq) = coef1 / sigma_sdf(1) * rho_sdf(1) * radc( iq )**3 &
348  * exp( &
349  - ( log( radc(iq) )-log( r0_sdf(1) ) )**2*0.5_rp &
350  / sigma_sdf(1) / sigma_sdf(1) &
351  )
352  enddo
353 
354  tmp_hyd = 0.0_rp
355  do iq = 1, nbin
356  tmp_hyd = tmp_hyd + dummy(iq)
357  enddo
358 
359  coef2 = ( &
360  qtrc_tmp(k+2,i,j,i_qc)+qtrc_tmp(k+2,i,j,i_qr) &
361  + qtrc_tmp(k+2,i,j,i_qs)+qtrc_tmp(k+2,i,j,i_qi) &
362  + qtrc_tmp(k+2,i,j,i_qg) &
363  ) &
364  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
365  do iq = 1, nbin
366  qtrc_org(k+2,i,j,i_qv+(il-1)*nbin+iq) = coef2 * dummy(iq)
367  enddo
368 
369  enddo
370  enddo
371  enddo
372 
373  elseif( nspc == 7 ) then !--- put each hydrometer to each category (ice bin)
374 
375  do k = 1, dims(1)
376  do i = xs, xe
377  do j = ys, ye
378 
379  !--- Rain and Cloud put into liquid bin (log-normal)
380  do iq = 1, nbin
381  dummy(iq) = coef1 / sigma_sdf(1) * rho_sdf(1) * radc( iq )**3 &
382  * exp( &
383  - ( log( radc(iq) )-log( r0_sdf(1) ) )**2*0.5_rp &
384  / sigma_sdf(1) / sigma_sdf(1) &
385  )
386  enddo
387 
388  tmp_hyd = 0.0_rp
389  do iq = 1, nbin
390  tmp_hyd = tmp_hyd + dummy(iq)
391  enddo
392 
393  coef2 = ( qtrc_tmp(k+2,i,j,i_qc)+qtrc_tmp(k+2,i,j,i_qr) ) &
394  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
395  do iq = 1, nbin
396  qtrc_org(k+2,i,j,i_qv+(il-1)*nbin+iq) = coef2 * dummy(iq)
397  enddo
398 
399  !--- Ice put into plate bin (log-normal)
400  do iq = 1, nbin
401  dummy(iq) = coef1 / sigma_sdf(2) * rho_sdf(2) * radc( iq )**3 &
402  * exp( &
403  - ( log( radc(iq) )-log( r0_sdf(2) ) )**2*0.5_rp &
404  / sigma_sdf(2) / sigma_sdf(2) &
405  )
406  enddo
407 
408  tmp_hyd = 0.0_rp
409  do iq = 1, nbin
410  tmp_hyd = tmp_hyd + dummy(iq)
411  enddo
412 
413  coef2 = qtrc_tmp(k+2,i,j,i_qi) &
414  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
415  do iq = 1, nbin
416  qtrc_org(k+2,i,j,i_qv+(ip-1)*nbin+iq) = coef2 * dummy(iq)
417  enddo
418 
419  !--- Snow put into snow bin (gamma)
420  n_hyd = coef0 * n0_sdf(3) * rho_sdf(3)
421  lambda_hyd = ( pi * rho_sdf(3) / 6.0_rp *n0_sdf(3) * sf_gamma(4.0_rp) &
422  / ( qtrc_tmp(k+2,i,j,i_qs) &
423  + (0.50_rp-sign(0.50_rp,qtrc_tmp(k+2,i,j,i_qs)-eps)) &
424  ) )**(0.25_rp)
425  do iq = 1, nbin
426  dummy(iq) = n_hyd * radc( iq )**3 &
427  * exp( -lambda_hyd * 0.5_rp * radc( iq ) )
428  enddo
429  tmp_hyd = 0.0_rp
430  do iq = 1, nbin
431  tmp_hyd = tmp_hyd + dummy(iq)
432  enddo
433 
434  coef2 = qtrc_tmp(k+2,i,j,i_qs) &
435  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
436  do iq = 1, nbin
437  qtrc_org(k+2,i,j,i_qv+(iss-1)*nbin+iq) = coef2 * dummy(iq)
438  enddo
439 
440  !--- Graupel put into Graupel bin (gamma)
441  n_hyd = coef0 * n0_sdf(4) * rho_sdf(4)
442  lambda_hyd = ( pi * rho_sdf(4) / 6.0_rp *n0_sdf(4) * sf_gamma(4.0_rp) &
443  / ( qtrc_tmp(k+2,i,j,i_qg) &
444  + (0.50_rp-sign(0.50_rp,qtrc_tmp(k+2,i,j,i_qg)-eps)) &
445  ) )**(0.25_rp)
446  do iq = 1, nbin
447  dummy(iq) = n_hyd * radc( iq )**3 &
448  * exp( -lambda_hyd * 0.5_rp * radc( iq ) )
449  enddo
450  tmp_hyd = 0.0_rp
451  do iq = 1, nbin
452  tmp_hyd = tmp_hyd + dummy(iq)
453  enddo
454 
455  coef2 = qtrc_tmp(k+2,i,j,i_qg) &
456  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
457  do iq = 1, nbin
458  qtrc_org(k+2,i,j,i_qv+(ig-1)*nbin+iq) = coef2 * dummy(iq)
459  enddo
460 
461  enddo
462  enddo
463  enddo
464 
465  endif
466 
467  elseif( trim(mp_type_outer) == "SN14" ) then
468 
469  write(*,*) 'SN14 is not supported for MP_TYPE_OUTER now'
470  write(*,*) 'Please wait'
471  call prc_mpistop
472 
473  elseif( trim(mp_type_outer) == "SUZUKI10" ) then
474 
475  if( io_l ) write(io_fid_log,*)
476  if( io_l ) write(io_fid_log,*) '+++ SDF of Bin model is created directory'
477  if( io_l ) write(io_fid_log,*) '+++ from Bin microphysical model'
478  do iq = 1, qa
479  call fileread( read3di(:,:,:), basename_org, aq_name(iq), it, rank )
480  do k = 1, dims(1)
481  qtrc_org(k+2,xs:xe,ys:ye,iq) = read3di(:,:,k)
482  end do
483  end do
484 
485  else
486 
487  write(*,*) 'MP_TYPE_OUTER should be KESSLER, TOMITA08, or SUZUKI10'
488  write(*,*) 'Please check! Now MP_TYPE_OUTER set as ', mp_type_outer
489  call prc_mpistop
490 
491  endif
492 
493  deallocate( read3di )
494  deallocate( qtrc_tmp )
495  return
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public prc_mpistop
Abort MPI.
integer, dimension(2), public parent_jmax
parent max number in y-direction
module GRID (nesting system)
real(rp) function, public sf_gamma(x)
Gamma function.
integer, public qa
integer, public i_qv
character(len=h_short), dimension(:), allocatable, public aq_name
module SPECFUNC
subroutine, public log(type, message)
Definition: dc_log.f90:133
integer, public i_qs
integer, public i_qi
integer, dimension(2), public parent_imax
parent max number in x-direction
integer, public i_qg
integer, public i_qr
integer, public i_qc
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_bin2bulk()

subroutine, public scale_atmos_phy_mp_convert::atmos_phy_mp_bin2bulk ( )

Bin to Bulk.

Definition at line 501 of file scale_atmos_phy_mp_convert.F90.

501  !---------------------------------------------------------------------------
502 
503  return