ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
mp_radar Module Reference

Functions/Subroutines

subroutine, public radar_init
 
subroutine, public rayleigh_soak_wetgraupel (x_g, a_geo, b_geo, fmelt, meltratio_outside, m_w, m_i, lambda, c_back, mixingrule, matrix, inclusion, host, hostmatrix, hostinclusion)
 
real(kind=kind_phys) function wgamma (y)
 
real(kind=kind_phys) function gammln (xx)
 
complex(kind=r8kind) function get_m_mix_nested (m_a, m_i, m_w, volair, volice, volwater, mixingrule, host, matrix, inclusion, hostmatrix, hostinclusion, cumulerror)
 
complex(kind=r8kind) function get_m_mix (m_a, m_i, m_w, volair, volice, volwater, mixingrule, matrix, inclusion, error)
 
complex(kind=r8kind) function m_complex_maxwellgarnett (vol1, vol2, vol3, m1, m2, m3, inclusion, error)
 
complex(kind=r8kind) function m_complex_water_ray (lambda, t)
 
complex(kind=r8kind) function m_complex_ice_maetzler (lambda, t)
 

Variables

integer, parameter kind_phys = c_double
 
integer, parameter, private r4kind = selected_real_kind(6)
 
integer, parameter, private r8kind = selected_real_kind(12)
 
integer, parameter, public nrbins = 50
 
integer, parameter, public slen = 20
 
character(len=slen), public mixingrulestring_s
 
character(len=slen), public matrixstring_s
 
character(len=slen), public inclusionstring_s
 
character(len=slen), public hoststring_s
 
character(len=slen), public hostmatrixstring_s
 
character(len=slen), public hostinclusionstring_s
 
character(len=slen), public mixingrulestring_g
 
character(len=slen), public matrixstring_g
 
character(len=slen), public inclusionstring_g
 
character(len=slen), public hoststring_g
 
character(len=slen), public hostmatrixstring_g
 
character(len=slen), public hostinclusionstring_g
 
complex(kind=r8kind), public m_w_0
 
complex(kind=r8kind), public m_i_0
 
double precision, dimension(nrbins+1), public xxdx
 
double precision, dimension(nrbins), public xxds
 
double precision, dimension(nrbins), public xdts
 
double precision, dimension(nrbins), public xxdg
 
double precision, dimension(nrbins), public xdtg
 
double precision, parameter, public lamda_radar = 0.10
 
double precision, public k_w
 
double precision, public pi5
 
double precision, public lamda4
 
double precision, dimension(nrbins+1), public simpson
 
double precision, dimension(3), parameter, public basis = (/1.d0/3.d0, 4.d0/3.d0, 1.d0/3.d0/)
 
real(kind=kind_phys), dimension(4), public xcre
 
real(kind=kind_phys), dimension(4), public xcse
 
real(kind=kind_phys), dimension(4), public xcge
 
real(kind=kind_phys), dimension(4), public xcrg
 
real(kind=kind_phys), dimension(4), public xcsg
 
real(kind=kind_phys), dimension(4), public xcgg
 
real(kind=kind_phys), public xam_r
 
real(kind=kind_phys), public xbm_r
 
real(kind=kind_phys), public xmu_r
 
real(kind=kind_phys), public xobmr
 
real(kind=kind_phys), public xam_s
 
real(kind=kind_phys), public xbm_s
 
real(kind=kind_phys), public xmu_s
 
real(kind=kind_phys), public xoams
 
real(kind=kind_phys), public xobms
 
real(kind=kind_phys), public xocms
 
real(kind=kind_phys), public xam_g
 
real(kind=kind_phys), public xbm_g
 
real(kind=kind_phys), public xmu_g
 
real(kind=kind_phys), public xoamg
 
real(kind=kind_phys), public xobmg
 
real(kind=kind_phys), public xocmg
 
real(kind=kind_phys), public xorg2
 
real(kind=kind_phys), public xosg2
 
real(kind=kind_phys), public xogg2
 
character(len=256) radar_debug
 
double precision, parameter, public melt_outside_s = 0.9d0
 
double precision, parameter, public melt_outside_g = 0.9d0
 

Function/Subroutine Documentation

◆ gammln()

real(kind=kind_phys) function mp_radar::gammln ( real(kind=kind_phys), intent(in)  xx)
private
311  implicit none
312 !(C) Copr. 1986-92 Numerical Recipes Software 2.02
313 !=================================================================================================================
314 
315 !--- inout arguments:
316  real(kind=kind_phys),intent(in):: xx
317 
318 !--- local variables:
319  integer:: j
320 
321  double precision,parameter:: stp = 2.5066282746310005d0
322  double precision,dimension(6), parameter:: &
323  cof = (/76.18009172947146d0, -86.50532032941677d0, &
324  24.01409824083091d0, -1.231739572450155d0, &
325  .1208650973866179d-2, -.5395239384953d-5/)
326  double precision:: ser,tmp,x,y
327 
328 !-----------------------------------------------------------------------------------------------------------------
329 
330 !--- returns the value ln(gamma(xx)) for xx > 0.
331  x = xx
332  y = x
333  tmp = x+5.5d0
334  tmp = (x+0.5d0)*log(tmp)-tmp
335  ser = 1.000000000190015d0
336  do j = 1,6
337  y=y+1.d0
338  ser=ser+cof(j)/y
339  enddo
340 
341  gammln=tmp+log(stp*ser/x)
342 

Referenced by wgamma().

Here is the caller graph for this function:

◆ get_m_mix()

complex(kind=r8kind) function mp_radar::get_m_mix ( complex(kind=r8kind), intent(in)  m_a,
complex(kind=r8kind), intent(in)  m_i,
complex(kind=r8kind), intent(in)  m_w,
double precision, intent(in)  volair,
double precision, intent(in)  volice,
double precision, intent(in)  volwater,
character(len=*), intent(in)  mixingrule,
character(len=*), intent(in)  matrix,
character(len=*), intent(in)  inclusion,
integer, intent(out)  error 
)
private
490  implicit none
491 !=================================================================================================================
492 
493 !--- input arguments:
494  character(len=*),intent(in):: mixingrule, matrix, inclusion
495 
496  complex(kind=R8KIND), intent(in):: m_a, m_i, m_w
497 
498  double precision, intent(in):: volice, volair, volwater
499 
500 !--- output arguments:
501  integer,intent(out):: error
502 
503 !-----------------------------------------------------------------------------------------------------------------
504  error = 0
505  get_m_mix = cmplx(1.0d0,0.0d0)
506 
507  if (mixingrule .eq. 'maxwellgarnett') then
508  if (matrix .eq. 'ice') then
509  get_m_mix = m_complex_maxwellgarnett(volice, volair, volwater, &
510  m_i, m_a, m_w, inclusion, error)
511  elseif (matrix .eq. 'water') then
512  get_m_mix = m_complex_maxwellgarnett(volwater, volair, volice, &
513  m_w, m_a, m_i, inclusion, error)
514  elseif (matrix .eq. 'air') then
515  get_m_mix = m_complex_maxwellgarnett(volair, volwater, volice, &
516  m_a, m_w, m_i, inclusion, error)
517  else
518  write(radar_debug,*) 'GET_M_MIX: unknown matrix: ', matrix
519 ! call physics_message(radar_debug)
520  error = 1
521  endif
522 
523  else
524  write(radar_debug,*) 'GET_M_MIX: unknown mixingrule: ', mixingrule
525 ! call physics_message(radar_debug)
526  error = 2
527  endif
528 
529  if (error .ne. 0) then
530  write(radar_debug,*) 'GET_M_MIX: error encountered'
531 ! call physics_message(radar_debug)
532  endif
533 

Referenced by get_m_mix_nested().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_m_mix_nested()

complex(kind=r8kind) function mp_radar::get_m_mix_nested ( complex(kind=r8kind), intent(in)  m_a,
complex(kind=r8kind), intent(in)  m_i,
complex(kind=r8kind), intent(in)  m_w,
double precision, intent(in)  volair,
double precision, intent(in)  volice,
double precision, intent(in)  volwater,
character(len=*), intent(in)  mixingrule,
character(len=*), intent(in)  host,
character(len=*), intent(in)  matrix,
character(len=*), intent(in)  inclusion,
character(len=*), intent(in)  hostmatrix,
character(len=*), intent(in)  hostinclusion,
integer, intent(out)  cumulerror 
)
private
349  implicit none
350 !=================================================================================================================
351 
352 !--- input arguments:
353  character(len=*),intent(in):: mixingrule, host, matrix, &
354  inclusion, hostmatrix, hostinclusion
355 
356  complex(kind=R8KIND),intent(in):: m_a, m_i, m_w
357 
358  double precision,intent(in):: volice, volair, volwater
359 
360 !--- output arguments:
361  integer,intent(out):: cumulerror
362 
363 !--- local variables:
364  integer:: error
365 
366  complex(kind=R8KIND):: mtmp
367 
368  double precision:: vol1, vol2
369 
370 !-----------------------------------------------------------------------------------------------------------------
371 
372 !..Folded: ( (m1 + m2) + m3), where m1,m2,m3 could each be air, ice, or water
373  cumulerror = 0
374  get_m_mix_nested = cmplx(1.0d0,0.0d0)
375 
376  if (host .eq. 'air') then
377  if (matrix .eq. 'air') then
378  write(radar_debug,*) 'GET_M_MIX_NESTED: bad matrix: ', matrix
379 ! call physics_message(radar_debug)
380  cumulerror = cumulerror + 1
381  else
382  vol1 = volice / max(volice+volwater,1d-10)
383  vol2 = 1.0d0 - vol1
384  mtmp = get_m_mix(m_a, m_i, m_w, 0.0d0, vol1, vol2, &
385  mixingrule, matrix, inclusion, error)
386  cumulerror = cumulerror + error
387 
388  if (hostmatrix .eq. 'air') then
389  get_m_mix_nested = get_m_mix(m_a, mtmp, 2.0*m_a, &
390  volair, (1.0d0-volair), 0.0d0, mixingrule, &
391  hostmatrix, hostinclusion, error)
392  cumulerror = cumulerror + error
393  elseif (hostmatrix .eq. 'icewater') then
394  get_m_mix_nested = get_m_mix(m_a, mtmp, 2.0*m_a, &
395  volair, (1.0d0-volair), 0.0d0, mixingrule, &
396  'ice', hostinclusion, error)
397  cumulerror = cumulerror + error
398  else
399  write(radar_debug,*) 'GET_M_MIX_NESTED: bad hostmatrix: ', hostmatrix
400 ! call physics_message(radar_debug)
401  cumulerror = cumulerror + 1
402  endif
403  endif
404 
405  elseif (host .eq. 'ice') then
406 
407  if (matrix .eq. 'ice') then
408  write(radar_debug,*) 'GET_M_MIX_NESTED: bad matrix: ', matrix
409 ! call physics_message(radar_debug)
410  cumulerror = cumulerror + 1
411  else
412  vol1 = volair / max(volair+volwater,1d-10)
413  vol2 = 1.0d0 - vol1
414  mtmp = get_m_mix(m_a, m_i, m_w, vol1, 0.0d0, vol2, &
415  mixingrule, matrix, inclusion, error)
416  cumulerror = cumulerror + error
417 
418  if (hostmatrix .eq. 'ice') then
419  get_m_mix_nested = get_m_mix(mtmp, m_i, 2.0*m_a, &
420  (1.0d0-volice), volice, 0.0d0, mixingrule, &
421  hostmatrix, hostinclusion, error)
422  cumulerror = cumulerror + error
423  elseif (hostmatrix .eq. 'airwater') then
424  get_m_mix_nested = get_m_mix(mtmp, m_i, 2.0*m_a, &
425  (1.0d0-volice), volice, 0.0d0, mixingrule, &
426  'air', hostinclusion, error)
427  cumulerror = cumulerror + error
428  else
429  write(radar_debug,*) 'GET_M_MIX_NESTED: bad hostmatrix: ', hostmatrix
430 ! call physics_message(radar_debug)
431  cumulerror = cumulerror + 1
432  endif
433  endif
434 
435  elseif (host .eq. 'water') then
436 
437  if (matrix .eq. 'water') then
438  write(radar_debug,*) 'GET_M_MIX_NESTED: bad matrix: ', matrix
439 ! call physics_message(radar_debug)
440  cumulerror = cumulerror + 1
441  else
442  vol1 = volair / max(volice+volair,1d-10)
443  vol2 = 1.0d0 - vol1
444  mtmp = get_m_mix(m_a, m_i, m_w, vol1, vol2, 0.0d0, &
445  mixingrule, matrix, inclusion, error)
446  cumulerror = cumulerror + error
447 
448  if (hostmatrix .eq. 'water') then
449  get_m_mix_nested = get_m_mix(2*m_a, mtmp, m_w, &
450  0.0d0, (1.0d0-volwater), volwater, mixingrule, &
451  hostmatrix, hostinclusion, error)
452  cumulerror = cumulerror + error
453  elseif (hostmatrix .eq. 'airice') then
454  get_m_mix_nested = get_m_mix(2*m_a, mtmp, m_w, &
455  0.0d0, (1.0d0-volwater), volwater, mixingrule, &
456  'ice', hostinclusion, error)
457  cumulerror = cumulerror + error
458  else
459  write(radar_debug,*) 'GET_M_MIX_NESTED: bad hostmatrix: ', hostmatrix
460 ! call physics_message(radar_debug)
461  cumulerror = cumulerror + 1
462  endif
463  endif
464 
465  elseif (host .eq. 'none') then
466 
467  get_m_mix_nested = get_m_mix(m_a, m_i, m_w, &
468  volair, volice, volwater, mixingrule, &
469  matrix, inclusion, error)
470  cumulerror = cumulerror + error
471 
472  else
473  write(radar_debug,*) 'GET_M_MIX_NESTED: unknown matrix: ', host
474 ! call physics_message(radar_debug)
475  cumulerror = cumulerror + 1
476  endif
477 
478  if (cumulerror .ne. 0) then
479  write(radar_debug,*) 'get_m_mix_nested: error encountered'
480 ! call physics_message(radar_debug)
481  get_m_mix_nested = cmplx(1.0d0,0.0d0)
482  endif
483 

Referenced by rayleigh_soak_wetgraupel().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ m_complex_ice_maetzler()

complex(kind=r8kind) function mp_radar::m_complex_ice_maetzler ( double precision, intent(in)  lambda,
double precision, intent(in)  t 
)
private
634  implicit none
635 !=================================================================================================================
636 
637 !complex refractive index of ice as function of Temperature T
638 ![deg C] and radar wavelength lambda [m]; valid for
639 !lambda in [0.0001,30] m; T in [-250.0,0.0] C
640 !Original comment from the Matlab-routine of Prof. Maetzler:
641 !Function for calculating the relative permittivity of pure ice in
642 !the microwave region, according to C. Maetzler, "Microwave
643 !properties of ice and snow", in B. Schmitt et al. (eds.) Solar
644 !System Ices, Astrophys. and Space Sci. Library, Vol. 227, Kluwer
645 !Academic Publishers, Dordrecht, pp. 241-257 (1998). Input:
646 !TK = temperature (K), range 20 to 273.15
647 !f = frequency in GHz, range 0.01 to 3000
648 
649 !--- input arguments:
650  double precision,intent(in):: t,lambda
651 
652 !--- local variables:
653  double precision:: f,c,tk,b1,b2,b,deltabeta,betam,beta,theta,alfa
654 
655 !-----------------------------------------------------------------------------------------------------------------
656 
657  c = 2.99d8
658  tk = t + 273.16
659  f = c / lambda * 1d-9
660 
661  b1 = 0.0207
662  b2 = 1.16d-11
663  b = 335.0d0
664  deltabeta = exp(-10.02 + 0.0364*(tk-273.16))
665  betam = (b1/tk) * ( exp(b/tk) / ((exp(b/tk)-1)**2) ) + b2*f*f
666  beta = betam + deltabeta
667  theta = 300. / tk - 1.
668  alfa = (0.00504d0 + 0.0062d0*theta) * exp(-22.1d0*theta)
669  m_complex_ice_maetzler = 3.1884 + 9.1e-4*(tk-273.16)
670  m_complex_ice_maetzler = m_complex_ice_maetzler &
671  + cmplx(0.0d0, (alfa/f + beta*f))
672  m_complex_ice_maetzler = sqrt(conjg(m_complex_ice_maetzler))
673 

Referenced by radar_init().

Here is the caller graph for this function:

◆ m_complex_maxwellgarnett()

complex(kind=r8kind) function mp_radar::m_complex_maxwellgarnett ( double precision, intent(in)  vol1,
double precision, intent(in)  vol2,
double precision, intent(in)  vol3,
complex(kind=r8kind), intent(in)  m1,
complex(kind=r8kind), intent(in)  m2,
complex(kind=r8kind), intent(in)  m3,
character(len=*), intent(in)  inclusion,
integer, intent(out)  error 
)
private
539  implicit none
540 !=================================================================================================================
541 
542 !--- input arguments:
543  character(len=*),intent(in):: inclusion
544 
545  complex(kind=R8KIND),intent(in):: m1,m2,m3
546 
547  double precision,intent(in):: vol1,vol2,vol3
548 
549 
550 !--- output arguments:
551  integer,intent(out):: error
552 
553 !--- local variables:
554  complex(kind=R8KIND) :: beta2, beta3, m1t, m2t, m3t
555 
556 !-----------------------------------------------------------------------------------------------------------------
557 
558  error = 0
559 
560  if (dabs(vol1+vol2+vol3-1.0d0) .gt. 1d-6) then
561  write(radar_debug,*) 'M_COMPLEX_MAXWELLGARNETT: sum of the ', &
562  'partial volume fractions is not 1...ERROR'
563 ! call physics_message(radar_debug)
564  m_complex_maxwellgarnett = cmplx(-999.99d0,-999.99d0)
565  error = 1
566  return
567  endif
568 
569  m1t = m1**2
570  m2t = m2**2
571  m3t = m3**2
572 
573  if (inclusion .eq. 'spherical') then
574  beta2 = 3.0d0*m1t/(m2t+2.0d0*m1t)
575  beta3 = 3.0d0*m1t/(m3t+2.0d0*m1t)
576  elseif (inclusion .eq. 'spheroidal') then
577  beta2 = 2.0d0*m1t/(m2t-m1t) * (m2t/(m2t-m1t)*log(m2t/m1t)-1.0d0)
578  beta3 = 2.0d0*m1t/(m3t-m1t) * (m3t/(m3t-m1t)*log(m3t/m1t)-1.0d0)
579  else
580  write(radar_debug,*) 'M_COMPLEX_MAXWELLGARNETT: ', 'unknown inclusion: ', inclusion
581 ! call physics_message(radar_debug)
582  m_complex_maxwellgarnett=cmplx(-999.99d0,-999.99d0,kind=r8kind)
583  error = 1
584  return
585  endif
586 
587  m_complex_maxwellgarnett = sqrt(((1.0d0-vol2-vol3)*m1t + vol2*beta2*m2t + vol3*beta3*m3t) / &
588  (1.0d0-vol2-vol3+vol2*beta2+vol3*beta3))
589 

Referenced by get_m_mix().

Here is the caller graph for this function:

◆ m_complex_water_ray()

complex(kind=r8kind) function mp_radar::m_complex_water_ray ( double precision, intent(in)  lambda,
double precision, intent(in)  t 
)
private
594  implicit none
595 !=================================================================================================================
596 
597 !complex refractive Index of Water as function of Temperature T
598 ![deg C] and radar wavelength lambda [m]; valid for
599 !lambda in [0.001,1.0] m; T in [-10.0,30.0] deg C
600 !after Ray (1972)
601 
602 !--- input arguments:
603  double precision,intent(in):: t,lambda
604 
605 !--- local variables:
606  double precision,parameter:: pix=3.1415926535897932384626434d0
607  double precision:: epsinf,epss,epsr,epsi
608  double precision:: alpha,lambdas,sigma,nenner
609  complex(kind=R8KIND),parameter:: i = (0d0,1d0)
610 
611 !-----------------------------------------------------------------------------------------------------------------
612 
613  epsinf = 5.27137d0 + 0.02164740d0 * t - 0.00131198d0 * t*t
614  epss = 78.54d+0 * (1.0 - 4.579d-3 * (t - 25.0) &
615  + 1.190d-5 * (t - 25.0)*(t - 25.0) &
616  - 2.800d-8 * (t - 25.0)*(t - 25.0)*(t - 25.0))
617  alpha = -16.8129d0/(t+273.16) + 0.0609265d0
618  lambdas = 0.00033836d0 * exp(2513.98d0/(t+273.16)) * 1e-2
619 
620  nenner = 1.d0+2.d0*(lambdas/lambda)**(1d0-alpha)*sin(alpha*pix*0.5) &
621  + (lambdas/lambda)**(2d0-2d0*alpha)
622  epsr = epsinf + ((epss-epsinf) * ((lambdas/lambda)**(1d0-alpha) &
623  * sin(alpha*pix*0.5)+1d0)) / nenner
624  epsi = ((epss-epsinf) * ((lambdas/lambda)**(1d0-alpha) &
625  * cos(alpha*pix*0.5)+0d0)) / nenner &
626  + lambda*1.25664/1.88496
627 
628  m_complex_water_ray = sqrt(cmplx(epsr,-epsi))
629 

Referenced by radar_init().

Here is the caller graph for this function:

◆ radar_init()

subroutine, public mp_radar::radar_init
74  implicit none
75 !=================================================================================================================
76 
77  integer:: n
78 
79 !-----------------------------------------------------------------------------------------------------------------
80 
81  pi5 = 3.14159*3.14159*3.14159*3.14159*3.14159
82  lamda4 = lamda_radar*lamda_radar*lamda_radar*lamda_radar
83  m_w_0 = m_complex_water_ray(lamda_radar, 0.0d0)
84  m_i_0 = m_complex_ice_maetzler(lamda_radar, 0.0d0)
85  k_w = (abs( (m_w_0*m_w_0 - 1.0) /(m_w_0*m_w_0 + 2.0) ))**2
86 
87  do n = 1, nrbins+1
88  simpson(n) = 0.0d0
89  enddo
90  do n = 1, nrbins-1, 2
91  simpson(n) = simpson(n) + basis(1)
92  simpson(n+1) = simpson(n+1) + basis(2)
93  simpson(n+2) = simpson(n+2) + basis(3)
94  enddo
95 
96  do n = 1, slen
97  mixingrulestring_s(n:n) = char(0)
98  matrixstring_s(n:n) = char(0)
99  inclusionstring_s(n:n) = char(0)
100  hoststring_s(n:n) = char(0)
101  hostmatrixstring_s(n:n) = char(0)
102  hostinclusionstring_s(n:n) = char(0)
103  mixingrulestring_g(n:n) = char(0)
104  matrixstring_g(n:n) = char(0)
105  inclusionstring_g(n:n) = char(0)
106  hoststring_g(n:n) = char(0)
107  hostmatrixstring_g(n:n) = char(0)
108  hostinclusionstring_g(n:n) = char(0)
109  enddo
110 
111  mixingrulestring_s = 'maxwellgarnett'
112  hoststring_s = 'air'
113  matrixstring_s = 'water'
114  inclusionstring_s = 'spheroidal'
115  hostmatrixstring_s = 'icewater'
116  hostinclusionstring_s = 'spheroidal'
117 
118  mixingrulestring_g = 'maxwellgarnett'
119  hoststring_g = 'air'
120  matrixstring_g = 'water'
121  inclusionstring_g = 'spheroidal'
122  hostmatrixstring_g = 'icewater'
123  hostinclusionstring_g = 'spheroidal'
124 
125 !..Create bins of snow (from 100 microns up to 2 cm).
126  xxdx(1) = 100.d-6
127  xxdx(nrbins+1) = 0.02d0
128  do n = 2, nrbins
129  xxdx(n) = dexp(real(n-1,kind=r8kind)/real(nrbins,kind=r8kind) &
130  * dlog(xxdx(nrbins+1)/xxdx(1)) +dlog(xxdx(1)))
131  enddo
132  do n = 1, nrbins
133  xxds(n) = dsqrt(xxdx(n)*xxdx(n+1))
134  xdts(n) = xxdx(n+1) - xxdx(n)
135  enddo
136 
137 !..create bins of graupel (from 100 microns up to 5 cm).
138  xxdx(1) = 100.d-6
139  xxdx(nrbins+1) = 0.05d0
140  do n = 2, nrbins
141  xxdx(n) = dexp(real(n-1,kind=r8kind)/real(nrbins,kind=r8kind) &
142  * dlog(xxdx(nrbins+1)/xxdx(1)) +dlog(xxdx(1)))
143  enddo
144  do n = 1, nrbins
145  xxdg(n) = dsqrt(xxdx(n)*xxdx(n+1))
146  xdtg(n) = xxdx(n+1) - xxdx(n)
147  enddo
148 
149 
150 !.. The calling program must set the m(D) relations and gamma shape
151 !.. parameter mu for rain, snow, and graupel. Easily add other types
152 !.. based on the template here. For majority of schemes with simpler
153 !.. exponential number distribution, mu=0.
154 
155  xcre(1) = 1. + xbm_r
156  xcre(2) = 1. + xmu_r
157  xcre(3) = 4. + xmu_r
158  xcre(4) = 7. + xmu_r
159  do n = 1, 4
160  xcrg(n) = wgamma(xcre(n))
161  enddo
162  xorg2 = 1./xcrg(2)
163 
164  xcse(1) = 1. + xbm_s
165  xcse(2) = 1. + xmu_s
166  xcse(3) = 4. + xmu_s
167  xcse(4) = 7. + xmu_s
168  do n = 1, 4
169  xcsg(n) = wgamma(xcse(n))
170  enddo
171  xosg2 = 1./xcsg(2)
172 
173  xcge(1) = 1. + xbm_g
174  xcge(2) = 1. + xmu_g
175  xcge(3) = 4. + xmu_g
176  xcge(4) = 7. + xmu_g
177  do n = 1, 4
178  xcgg(n) = wgamma(xcge(n))
179  enddo
180  xogg2 = 1./xcgg(2)
181 
182  xobmr = 1./xbm_r
183  xoams = 1./xam_s
184  xobms = 1./xbm_s
185  xocms = xoams**xobms
186  xoamg = 1./xam_g
187  xobmg = 1./xbm_g
188  xocmg = xoamg**xobmg
189 
double real
Definition: ERF_OrbCosZenith.H:9

Referenced by mp_wsm6::mp_wsm6_init().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rayleigh_soak_wetgraupel()

subroutine, public mp_radar::rayleigh_soak_wetgraupel ( double precision, intent(in)  x_g,
double precision, intent(in)  a_geo,
double precision, intent(in)  b_geo,
double precision, intent(in)  fmelt,
double precision, intent(in)  meltratio_outside,
complex(kind=r8kind), intent(in)  m_w,
complex(kind=r8kind), intent(in)  m_i,
double precision, intent(in)  lambda,
double precision, intent(out)  c_back,
character(len=*), intent(in)  mixingrule,
character(len=*), intent(in)  matrix,
character(len=*), intent(in)  inclusion,
character(len=*), intent(in)  host,
character(len=*), intent(in)  hostmatrix,
character(len=*), intent(in)  hostinclusion 
)
195  implicit none
196 !=================================================================================================================
197 
198 !--- input arguments:
199  character(len=*), intent(in):: mixingrule, matrix, inclusion, &
200  host, hostmatrix, hostinclusion
201 
202  complex(kind=R8KIND),intent(in):: m_w, m_i
203 
204  double precision, intent(in):: x_g, a_geo, b_geo, fmelt, lambda, meltratio_outside
205 
206 !--- output arguments:
207  double precision,intent(out):: c_back
208 
209 !--- local variables:
210  integer:: error
211 
212  complex(kind=R8KIND):: m_core, m_air
213 
214  double precision, parameter:: pix=3.1415926535897932384626434d0
215  double precision:: d_large, d_g, rhog, x_w, xw_a, fm, fmgrenz, &
216  volg, vg, volair, volice, volwater, &
217  meltratio_outside_grenz, mra
218 
219 !-----------------------------------------------------------------------------------------------------------------
220 
221 !refractive index of air:
222  m_air = (1.0d0,0.0d0)
223 
224 !Limiting the degree of melting --- for safety:
225  fm = dmax1(dmin1(fmelt, 1.0d0), 0.0d0)
226 !Limiting the ratio of (melting on outside)/(melting on inside):
227  mra = dmax1(dmin1(meltratio_outside, 1.0d0), 0.0d0)
228 
229 !The relative portion of meltwater melting at outside should increase
230 !from the given input value (between 0 and 1)
231 !to 1 as the degree of melting approaches 1,
232 !so that the melting particle "converges" to a water drop.
233 !Simplest assumption is linear:
234  mra = mra + (1.0d0-mra)*fm
235 
236  x_w = x_g * fm
237 
238  d_g = a_geo * x_g**b_geo
239 
240  if(d_g .ge. 1d-12) then
241 
242  vg = pix/6. * d_g**3
243  rhog = dmax1(dmin1(x_g / vg, 900.0d0), 10.0d0)
244  vg = x_g / rhog
245 
246  meltratio_outside_grenz = 1.0d0 - rhog / 1000.
247 
248  if (mra .le. meltratio_outside_grenz) then
249  !..In this case, it cannot happen that, during melting, all the
250  !.. air inclusions within the ice particle get filled with
251  !.. meltwater. This only happens at the end of all melting.
252  volg = vg * (1.0d0 - mra * fm)
253 
254  else
255  !..In this case, at some melting degree fm, all the air
256  !.. inclusions get filled with meltwater.
257  fmgrenz=(900.0-rhog)/(mra*900.0-rhog+900.0*rhog/1000.)
258 
259  if (fm .le. fmgrenz) then
260  !.. not all air pockets are filled:
261  volg = (1.0 - mra * fm) * vg
262  else
263  !..all air pockets are filled with meltwater, now the
264  !.. entire ice sceleton melts homogeneously:
265  volg = (x_g - x_w) / 900.0 + x_w / 1000.
266  endif
267 
268  endif
269 
270  d_large = (6.0 / pix * volg) ** (1./3.)
271  volice = (x_g - x_w) / (volg * 900.0)
272  volwater = x_w / (1000. * volg)
273  volair = 1.0 - volice - volwater
274 
275  !..complex index of refraction for the ice-air-water mixture
276  !.. of the particle:
277  m_core = get_m_mix_nested(m_air, m_i, m_w, volair, volice, &
278  volwater, mixingrule, host, matrix, inclusion, &
279  hostmatrix, hostinclusion, error)
280  if (error .ne. 0) then
281  c_back = 0.0d0
282  return
283  endif
284 
285  !..rayleigh-backscattering coefficient of melting particle:
286  c_back = (abs((m_core**2-1.0d0)/(m_core**2+2.0d0)))**2 &
287  * pi5 * d_large**6 / lamda4
288 
289  else
290  c_back = 0.0d0
291  endif
292 

Referenced by mp_wsm6::refl10cm_wsm6().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wgamma()

real(kind=kind_phys) function mp_radar::wgamma ( real(kind=kind_phys), intent(in)  y)
private
297  implicit none
298 !=================================================================================================================
299 
300 !--- input arguments:
301  real(kind=kind_phys),intent(in):: y
302 
303 !-----------------------------------------------------------------------------------------------------------------
304 
305  wgamma = exp(gammln(y))
306 

Referenced by radar_init().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ basis

double precision, dimension(3), parameter, public mp_radar::basis = (/1.d0/3.d0, 4.d0/3.d0, 1.d0/3.d0/)
52  double precision, dimension(3), parameter, public:: basis = &
53  (/1.d0/3.d0, 4.d0/3.d0, 1.d0/3.d0/)

Referenced by radar_init().

◆ hostinclusionstring_g

character(len=slen), public mp_radar::hostinclusionstring_g

◆ hostinclusionstring_s

character(len=slen), public mp_radar::hostinclusionstring_s

◆ hostmatrixstring_g

character(len=slen), public mp_radar::hostmatrixstring_g

◆ hostmatrixstring_s

character(len=slen), public mp_radar::hostmatrixstring_s

◆ hoststring_g

character(len=slen), public mp_radar::hoststring_g

◆ hoststring_s

character(len=slen), public mp_radar::hoststring_s

◆ inclusionstring_g

character(len=slen), public mp_radar::inclusionstring_g

◆ inclusionstring_s

character(len=slen), public mp_radar::inclusionstring_s

◆ k_w

double precision, public mp_radar::k_w
49  double precision,public:: k_w,pi5,lamda4

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ kind_phys

integer, parameter mp_radar::kind_phys = c_double
6  integer, parameter :: kind_phys = c_double

◆ lamda4

double precision, public mp_radar::lamda4

◆ lamda_radar

double precision, parameter, public mp_radar::lamda_radar = 0.10
48  double precision,parameter,public:: lamda_radar = 0.10 ! in meters

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ m_i_0

complex(kind=r8kind), public mp_radar::m_i_0

◆ m_w_0

complex(kind=r8kind), public mp_radar::m_w_0
44  complex(kind=R8KIND),public:: m_w_0, m_i_0

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ matrixstring_g

character(len=slen), public mp_radar::matrixstring_g

◆ matrixstring_s

character(len=slen), public mp_radar::matrixstring_s

◆ melt_outside_g

double precision, parameter, public mp_radar::melt_outside_g = 0.9d0
66  double precision,parameter,public:: melt_outside_g = 0.9d0

Referenced by mp_wsm6::refl10cm_wsm6().

◆ melt_outside_s

double precision, parameter, public mp_radar::melt_outside_s = 0.9d0
65  double precision,parameter,public:: melt_outside_s = 0.9d0

Referenced by mp_wsm6::refl10cm_wsm6().

◆ mixingrulestring_g

character(len=slen), public mp_radar::mixingrulestring_g

◆ mixingrulestring_s

character(len=slen), public mp_radar::mixingrulestring_s
38  character(len=slen), public:: &
39  mixingrulestring_s, matrixstring_s, inclusionstring_s, &
40  hoststring_s, hostmatrixstring_s, hostinclusionstring_s, &
41  mixingrulestring_g, matrixstring_g, inclusionstring_g, &
42  hoststring_g, hostmatrixstring_g, hostinclusionstring_g

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ nrbins

integer, parameter, public mp_radar::nrbins = 50
36  integer,parameter,public:: nrbins = 50

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ pi5

double precision, public mp_radar::pi5

◆ r4kind

integer, parameter, private mp_radar::r4kind = selected_real_kind(6)
private
33  integer, parameter, private :: R4KIND = selected_real_kind(6)

◆ r8kind

integer, parameter, private mp_radar::r8kind = selected_real_kind(12)
private
34  integer, parameter, private :: R8KIND = selected_real_kind(12)

Referenced by m_complex_maxwellgarnett(), and radar_init().

◆ radar_debug

character(len=256) mp_radar::radar_debug
private
63  character(len=256):: radar_debug

Referenced by get_m_mix(), get_m_mix_nested(), and m_complex_maxwellgarnett().

◆ simpson

double precision, dimension(nrbins+1), public mp_radar::simpson
51  double precision, dimension(nrbins+1), public:: simpson

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ slen

integer, parameter, public mp_radar::slen = 20
37  integer,parameter,public:: slen = 20

Referenced by radar_init().

◆ xam_g

real(kind=kind_phys), public mp_radar::xam_g
58  real(kind=kind_phys),public:: xam_g,xbm_g,xmu_g,xoamg,xobmg,xocmg

Referenced by mp_wsm6::mp_wsm6_init(), radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ xam_r

real(kind=kind_phys), public mp_radar::xam_r
56  real(kind=kind_phys),public:: xam_r,xbm_r,xmu_r,xobmr

Referenced by mp_wsm6::mp_wsm6_init(), and mp_wsm6::refl10cm_wsm6().

◆ xam_s

real(kind=kind_phys), public mp_radar::xam_s
57  real(kind=kind_phys),public:: xam_s,xbm_s,xmu_s,xoams,xobms,xocms

Referenced by mp_wsm6::mp_wsm6_init(), radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ xbm_g

real(kind=kind_phys), public mp_radar::xbm_g

◆ xbm_r

real(kind=kind_phys), public mp_radar::xbm_r

◆ xbm_s

real(kind=kind_phys), public mp_radar::xbm_s

◆ xcge

real(kind=kind_phys), dimension(4), public mp_radar::xcge

◆ xcgg

real(kind=kind_phys), dimension(4), public mp_radar::xcgg

◆ xcre

real(kind=kind_phys), dimension(4), public mp_radar::xcre
55  real(kind=kind_phys),public,dimension(4):: xcre,xcse,xcge,xcrg,xcsg,xcgg

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ xcrg

real(kind=kind_phys), dimension(4), public mp_radar::xcrg

◆ xcse

real(kind=kind_phys), dimension(4), public mp_radar::xcse

◆ xcsg

real(kind=kind_phys), dimension(4), public mp_radar::xcsg

◆ xdtg

double precision, dimension(nrbins), public mp_radar::xdtg

◆ xdts

double precision, dimension(nrbins), public mp_radar::xdts

◆ xmu_g

real(kind=kind_phys), public mp_radar::xmu_g

◆ xmu_r

real(kind=kind_phys), public mp_radar::xmu_r

◆ xmu_s

real(kind=kind_phys), public mp_radar::xmu_s

◆ xoamg

real(kind=kind_phys), public mp_radar::xoamg

Referenced by radar_init().

◆ xoams

real(kind=kind_phys), public mp_radar::xoams

Referenced by radar_init().

◆ xobmg

real(kind=kind_phys), public mp_radar::xobmg

◆ xobmr

real(kind=kind_phys), public mp_radar::xobmr

Referenced by radar_init().

◆ xobms

real(kind=kind_phys), public mp_radar::xobms

◆ xocmg

real(kind=kind_phys), public mp_radar::xocmg

◆ xocms

real(kind=kind_phys), public mp_radar::xocms

◆ xogg2

real(kind=kind_phys), public mp_radar::xogg2

Referenced by radar_init().

◆ xorg2

real(kind=kind_phys), public mp_radar::xorg2
59  real(kind=kind_phys),public:: xorg2,xosg2,xogg2

Referenced by radar_init().

◆ xosg2

real(kind=kind_phys), public mp_radar::xosg2

Referenced by radar_init().

◆ xxdg

double precision, dimension(nrbins), public mp_radar::xxdg

◆ xxds

double precision, dimension(nrbins), public mp_radar::xxds
47  double precision,dimension(nrbins),public:: xxds,xdts,xxdg,xdtg

Referenced by radar_init(), and mp_wsm6::refl10cm_wsm6().

◆ xxdx

double precision, dimension(nrbins+1), public mp_radar::xxdx
46  double precision,dimension(nrbins+1),public:: xxdx

Referenced by radar_init().