46 public :: atmos_hydrostatic_buildrho
47 public :: atmos_hydrostatic_buildrho_real
48 public :: atmos_hydrostatic_buildrho_atmos
49 public :: atmos_hydrostatic_buildrho_bytemp
50 public :: atmos_hydrostatic_buildrho_bytemp_atmos
56 public :: atmos_hydrostatic_barometric_law_mslp
57 public :: atmos_hydrostatic_barometric_law_pres
59 interface atmos_hydrostatic_buildrho
62 end interface atmos_hydrostatic_buildrho
64 interface atmos_hydrostatic_buildrho_real
67 end interface atmos_hydrostatic_buildrho_real
69 interface atmos_hydrostatic_buildrho_atmos
72 end interface atmos_hydrostatic_buildrho_atmos
74 interface atmos_hydrostatic_buildrho_bytemp
77 end interface atmos_hydrostatic_buildrho_bytemp
79 interface atmos_hydrostatic_buildrho_bytemp_atmos
82 end interface atmos_hydrostatic_buildrho_bytemp_atmos
84 interface atmos_hydrostatic_barometric_law_mslp
86 module procedure atmos_hydrostatic_barometric_law_mslp_2d
87 end interface atmos_hydrostatic_barometric_law_mslp
89 interface atmos_hydrostatic_barometric_law_pres
90 module procedure atmos_hydrostatic_barometric_law_pres_0d
91 module procedure atmos_hydrostatic_barometric_law_pres_2d
92 end interface atmos_hydrostatic_barometric_law_pres
102 private :: atmos_hydrostatic_buildrho_atmos_2d
108 integer,
private,
parameter :: itelim = 100
109 real(RP),
private :: criteria
110 logical,
private :: hydrostatic_uselapserate = .false.
111 integer,
private :: hydrostatic_buildrho_real_kref = 1
114 real(RP),
private :: cv_qv
115 real(RP),
private :: cv_qc
128 namelist / param_atmos_hydrostatic / &
129 hydrostatic_uselapserate, &
130 hydrostatic_buildrho_real_kref
136 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[HYDROSTATIC] / Categ[ATMOS SHARE] / Origin[SCALElib]' 140 read(
io_fid_conf,nml=param_atmos_hydrostatic,iostat=ierr)
142 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 143 elseif( ierr > 0 )
then 144 write(*,*)
'xxx Not appropriate names in namelist PARAM_ATMOS_HYDROSTATIC. Check!' 152 if(
io_l )
write(
io_fid_log,*)
'*** use lapse rate for estimation of surface temperature? : ', hydrostatic_uselapserate
153 if(
io_l )
write(
io_fid_log,*)
'*** buildrho conversion criteria : ', criteria
155 if ( thermodyn_type ==
'EXACT' )
then 158 elseif( thermodyn_type ==
'SIMPLE' )
then 184 real(RP),
intent(out) :: dens(
ka)
185 real(RP),
intent(out) :: temp(
ka)
186 real(RP),
intent(out) :: pres(
ka)
187 real(RP),
intent(in) :: pott(
ka)
188 real(RP),
intent(in) :: qv (
ka)
189 real(RP),
intent(in) :: qc (
ka)
191 real(RP),
intent(out) :: temp_sfc
192 real(RP),
intent(in) :: pres_sfc
193 real(RP),
intent(in) :: pott_sfc
194 real(RP),
intent(in) :: qv_sfc
195 real(RP),
intent(in) :: qc_sfc
200 real(RP) :: CVtot_sfc
201 real(RP) :: CPovCV_sfc
206 real(RP) :: CVovCP_sfc, CPovR, CVovCP, RovCV
207 real(RP) :: dens_s, dhyd, dgrd
214 rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
216 cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
219 cpovcv_sfc = ( cvtot_sfc + rtot_sfc ) / cvtot_sfc
221 rtot = rdry * ( 1.0_rp - qv(
ks) - qc(
ks) ) &
223 cvtot = cvdry * ( 1.0_rp - qv(
ks) - qc(
ks) ) &
226 cpovcv = ( cvtot + rtot ) / cvtot
229 cvovcp_sfc = 1.0_rp / cpovcv_sfc
230 dens_sfc = p00 / rtot_sfc / pott_sfc * ( pres_sfc/p00 )**cvovcp_sfc
231 temp_sfc = pres_sfc / ( dens_sfc * rtot_sfc )
234 if ( hydrostatic_uselapserate )
then 236 cpovr = ( cvtot + rtot ) / rtot
237 cvovcp = 1.0_rp / cpovcv
240 pres(
ks) = p00 * ( temp(
ks)/pott(
ks) )**cpovr
241 dens(
ks) = p00 / rtot / pott(
ks) * ( pres(
ks)/p00 )**cvovcp
252 if ( abs(dens(
ks)-dens_s) <= criteria )
then 259 dhyd = + ( p00 * ( dens_sfc * rtot_sfc * pott_sfc / p00 )**cpovcv_sfc &
260 - p00 * ( dens_s * rtot * pott(
ks) / p00 )**cpovcv ) /
grid_cz(
ks) &
261 - grav * 0.5_rp * ( dens_sfc + dens_s )
263 dgrd = - p00 * ( rtot * pott(
ks) / p00 )**cpovcv /
grid_cz(
ks) &
264 * cpovcv * dens_s**rovcv &
267 dens(
ks) = dens_s - dhyd/dgrd
269 if( dens(
ks)*0.0_rp /= 0.0_rp)
exit 272 if ( .NOT. converged )
then 273 write(*,*)
'xxx [buildrho 1D sfc] iteration not converged!', &
274 dens(
ks),ite,dens_s,dhyd,dgrd
311 real(RP),
intent(out) :: dens(
ka,
ia,
ja)
312 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
313 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
314 real(RP),
intent(in) :: pott(
ka,
ia,
ja)
315 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
316 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
318 real(RP),
intent(out) :: temp_sfc(1,
ia,
ja)
319 real(RP),
intent(in) :: pres_sfc(1,
ia,
ja)
320 real(RP),
intent(in) :: pott_sfc(1,
ia,
ja)
321 real(RP),
intent(in) :: qv_sfc (1,
ia,
ja)
322 real(RP),
intent(in) :: qc_sfc (1,
ia,
ja)
326 real(RP) :: dens_sfc (1,
ia,
ja)
327 real(RP) :: pott_toa (
ia,
ja)
328 real(RP) :: qv_toa (
ia,
ja)
329 real(RP) :: qc_toa (
ia,
ja)
330 real(RP) :: dens_1D (
ka)
332 real(RP) :: Rtot_sfc (
ia,
ja)
333 real(RP) :: CVtot_sfc (
ia,
ja)
334 real(RP) :: CPovCV_sfc(
ia,
ja)
335 real(RP) :: Rtot (
ia,
ja)
336 real(RP) :: CVtot (
ia,
ja)
337 real(RP) :: CPovCV (
ia,
ja)
339 real(RP) :: CVovCP_sfc, CPovR, CVovCP
358 pott_toa(i,j) = pott(
ke,i,j)
359 qv_toa(i,j) = qv(
ke,i,j)
360 qc_toa(i,j) = qc(
ke,i,j)
367 rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
368 + rvap * qv_sfc(1,i,j)
369 cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
370 + cv_qv * qv_sfc(1,i,j) &
371 + cv_qc * qc_sfc(1,i,j)
372 cpovcv_sfc(i,j) = ( cvtot_sfc(i,j) + rtot_sfc(i,j) ) / cvtot_sfc(i,j)
378 rtot(i,j) = rdry * ( 1.0_rp - qv(
ks,i,j) - qc(
ks,i,j) ) &
380 cvtot(i,j) = cvdry * ( 1.0_rp - qv(
ks,i,j) - qc(
ks,i,j) ) &
381 + cv_qv * qv(
ks,i,j) &
383 cpovcv(i,j) = ( cvtot(i,j) + rtot(i,j) ) / cvtot(i,j)
390 cvovcp_sfc = 1.0_rp / cpovcv_sfc(i,j)
391 dens_sfc(1,i,j) = p00 / rtot_sfc(i,j) / pott_sfc(1,i,j) * ( pres_sfc(1,i,j)/p00 )**cvovcp_sfc
392 temp_sfc(1,i,j) = pres_sfc(1,i,j) / ( dens_sfc(1,i,j) * rtot_sfc(i,j) )
397 if ( hydrostatic_uselapserate )
then 401 cpovr = ( cvtot(i,j) + rtot(i,j) ) / rtot(i,j)
402 cvovcp = 1.0_rp / cpovcv(i,j)
405 pres(
ks,i,j) = p00 * ( temp(
ks,i,j)/pott(
ks,i,j) )**cpovr
406 dens(
ks,i,j) = p00 / rtot(i,j) / pott(
ks,i,j) * ( pres(
ks,i,j)/p00 )**cvovcp
412 call atmos_hydrostatic_buildrho_atmos_2d( dens(
ks,:,:), &
436 call atmos_hydrostatic_buildrho_atmos_2d( dens(
ke+1,:,:), &
450 dens( 1:
ks-1,:,:) = 0.d0
451 dens(
ke+2:
ka ,:,:) = 0.d0
456 dens(:,i,j) = dens_1d(:)
498 real(RP),
intent(out) :: dens(
ka,
ia,
ja)
499 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
500 real(RP),
intent(inout) :: pres(
ka,
ia,
ja)
501 real(RP),
intent(in) :: pott(
ka,
ia,
ja)
502 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
503 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
507 real(RP) :: pott_toa(
ia,
ja)
508 real(RP) :: qv_toa (
ia,
ja)
509 real(RP) :: qc_toa (
ia,
ja)
533 pott_toa(i,j) = pott(
ke,i,j)
534 qv_toa(i,j) = qv(
ke,i,j)
535 qc_toa(i,j) = qc(
ke,i,j)
539 kref = hydrostatic_buildrho_real_kref +
ks - 1
544 rtot = rdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
545 + rvap * qv(kref,i,j)
546 cvtot = cvdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
547 + cv_qv * qv(kref,i,j) &
548 + cv_qc * qc(kref,i,j)
549 cvovcp = cvtot / ( cvtot + rtot )
550 dens(kref,i,j) = p00 / ( rtot * pott(kref,i,j) ) * ( pres(kref,i,j)/p00 )**cvovcp
572 call atmos_hydrostatic_buildrho_atmos_2d( dens(
ke+1,:,:), &
586 dens( 1:
ks-1,:,:) = 0.d0
587 dens(
ke+2:
ka ,:,:) = 0.d0
611 real(RP),
intent(out) :: dens_L2
612 real(RP),
intent(out) :: temp_L2
613 real(RP),
intent(out) :: pres_L2
614 real(RP),
intent(in) :: pott_L2
615 real(RP),
intent(in) :: qv_L2
616 real(RP),
intent(in) :: qc_L2
617 real(RP),
intent(in) :: dens_L1
618 real(RP),
intent(in) :: pott_L1
619 real(RP),
intent(in) :: qv_L1
620 real(RP),
intent(in) :: qc_L1
621 real(RP),
intent(in) :: dz
622 integer,
intent(in) :: k
624 real(RP) :: Rtot_L1 , Rtot_L2
625 real(RP) :: CVtot_L1 , CVtot_L2
626 real(RP) :: CPovCV_L1, CPovCV_L2
629 real(RP) :: dens_s, dhyd, dgrd
634 rtot_l1 = rdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
636 cvtot_l1 = cvdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
639 cpovcv_l1 = ( cvtot_l1 + rtot_l1 ) / cvtot_l1
641 rtot_l2 = rdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
643 cvtot_l2 = cvdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
646 cpovcv_l2 = ( cvtot_l2 + rtot_l2 ) / cvtot_l2
648 rovcv = rtot_l2 / cvtot_l2
655 if ( abs(dens_l2-dens_s) <= criteria )
then 662 dhyd = + ( p00 * ( dens_l1 * rtot_l1 * pott_l1 / p00 )**cpovcv_l1 &
663 - p00 * ( dens_s * rtot_l2 * pott_l2 / p00 )**cpovcv_l2 ) / dz &
664 - grav * 0.5_rp * ( dens_l1 + dens_s )
666 dgrd = - p00 * ( rtot_l2 * pott_l2 / p00 )**cpovcv_l2 / dz &
667 * cpovcv_l2 * dens_s**rovcv &
670 dens_l2 = dens_s - dhyd/dgrd
672 if( dens_l2*0.0_rp /= 0.0_rp)
exit 675 if ( .NOT. converged )
then 676 write(*,*)
'xxx [buildrho 0D atmos] iteration not converged!', &
677 k,dens_l2,ite,dens_s,dhyd,dgrd
681 pres_l2 = p00 * ( dens_l2 * rtot_l2 * pott_l2 / p00 )**cpovcv_l2
682 temp_l2 = pres_l2 / ( dens_l2 * rtot_l2 )
700 real(RP),
intent(inout) :: dens(
ka)
701 real(RP),
intent(out) :: temp(
ka)
702 real(RP),
intent(out) :: pres(
ka)
703 real(RP),
intent(in) :: pott(
ka)
704 real(RP),
intent(in) :: qv (
ka)
705 real(RP),
intent(in) :: qc (
ka)
707 real(RP) :: Rtot (
ka)
708 real(RP) :: CVtot (
ka)
709 real(RP) :: CPovCV(
ka)
712 real(RP) :: dens_s, dhyd, dgrd
720 rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
722 cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
725 cpovcv(k) = ( cvtot(k) + rtot(k) ) / cvtot(k)
729 rovcv = rtot(k) / cvtot(k)
736 if ( abs(dens(k)-dens_s) <= criteria )
then 743 dhyd = + ( p00 * ( dens(k-1) * rtot(k-1) * pott(k-1) / p00 )**cpovcv(k-1) &
744 - p00 * ( dens_s * rtot(k ) * pott(k ) / p00 )**cpovcv(k ) ) /
grid_fdz(k-1) &
745 - grav * 0.5_rp * ( dens(k-1) + dens_s )
747 dgrd = - p00 * ( rtot(k) * pott(k) / p00 )**cpovcv(k) /
grid_fdz(k-1) &
748 * cpovcv(k) * dens_s**rovcv &
751 dens(k) = dens_s - dhyd/dgrd
753 if( dens(k)*0.0_rp /= 0.0_rp)
exit 756 if ( .NOT. converged )
then 757 write(*,*)
'xxx [buildrho 1D atmos] iteration not converged!', &
758 k,dens(k),ite,dens_s,dhyd,dgrd
764 pres(k) = p00 * ( dens(k) * rtot(k) * pott(k) / p00 )**cpovcv(k)
765 temp(k) = pres(k) / ( dens(k) * rtot(k) )
768 dens( 1:
ks-1) = dens(
ks)
770 pres( 1:
ks-1) = pres(
ks)
772 temp( 1:
ks-1) = temp(
ks)
780 subroutine atmos_hydrostatic_buildrho_atmos_2d( &
797 real(RP),
intent(out) :: dens_L2(
ia,
ja)
798 real(RP),
intent(out) :: temp_L2(
ia,
ja)
799 real(RP),
intent(out) :: pres_L2(
ia,
ja)
800 real(RP),
intent(in) :: pott_L2(
ia,
ja)
801 real(RP),
intent(in) :: qv_L2 (
ia,
ja)
802 real(RP),
intent(in) :: qc_L2 (
ia,
ja)
803 real(RP),
intent(in) :: dens_L1(
ia,
ja)
804 real(RP),
intent(in) :: pott_L1(
ia,
ja)
805 real(RP),
intent(in) :: qv_L1 (
ia,
ja)
806 real(RP),
intent(in) :: qc_L1 (
ia,
ja)
807 real(RP),
intent(in) :: dz (
ia,
ja)
808 integer,
intent(in) :: k
810 real(RP) :: Rtot_L1 (
ia,
ja), Rtot_L2 (
ia,
ja)
811 real(RP) :: CVtot_L1 (
ia,
ja), CVtot_L2 (
ia,
ja)
812 real(RP) :: CPovCV_L1(
ia,
ja), CPovCV_L2(
ia,
ja)
815 real(RP) :: dens_s, dhyd, dgrd
824 rtot_l1(i,j) = rdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
826 cvtot_l1(i,j) = cvdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
827 + cv_qv * qv_l1(i,j) &
829 cpovcv_l1(i,j) = ( cvtot_l1(i,j) + rtot_l1(i,j) ) / cvtot_l1(i,j)
831 rtot_l2(i,j) = rdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
833 cvtot_l2(i,j) = cvdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
834 + cv_qv * qv_l2(i,j) &
836 cpovcv_l2(i,j) = ( cvtot_l2(i,j) + rtot_l2(i,j) ) / cvtot_l2(i,j)
842 rovcv = rtot_l2(i,j) / cvtot_l2(i,j)
845 dens_l2(i,j) = dens_l1(i,j)
849 if ( abs(dens_l2(i,j)-dens_s) <= criteria )
then 854 dens_s = dens_l2(i,j)
856 dhyd = + ( p00 * ( dens_l1(i,j) * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) &
857 - p00 * ( dens_s * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) ) / dz(i,j) &
858 - grav * 0.5_rp * ( dens_l1(i,j) + dens_s )
860 dgrd = - p00 * ( rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) / dz(i,j) &
861 * cpovcv_l2(i,j) * dens_s**rovcv &
864 dens_l2(i,j) = dens_s - dhyd/dgrd
866 if( dens_l2(i,j)*0.0_rp /= 0.0_rp)
exit 869 if ( .NOT. converged )
then 870 write(*,*)
'xxx [buildrho 2D atmos] iteration not converged!', &
871 i,j,k,dens_l2(i,j),ite,dens_s,dhyd,dgrd
879 pres_l2(i,j) = p00 * ( dens_l2(i,j) * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j)
880 temp_l2(i,j) = pres_l2(i,j) / ( dens_l2(i,j) * rtot_l2(i,j) )
885 end subroutine atmos_hydrostatic_buildrho_atmos_2d
902 real(RP),
intent(inout) :: dens(
ka,
ia,
ja)
903 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
904 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
905 real(RP),
intent(in) :: pott(
ka,
ia,
ja)
906 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
907 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
908 real(RP),
intent(in) :: dz (
ka,
ia,
ja)
909 integer,
intent(in),
optional :: kref_in
911 real(RP) :: Rtot (
ka,
ia,
ja)
912 real(RP) :: CVtot (
ka,
ia,
ja)
913 real(RP) :: CPovCV(
ka,
ia,
ja)
916 real(RP) :: dens_s, dhyd, dgrd
924 if (
present(kref_in) )
then 933 rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
935 cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
936 + cv_qv * qv(k,i,j) &
938 cpovcv(k,i,j) = ( cvtot(k,i,j) + rtot(k,i,j) ) / cvtot(k,i,j)
946 rovcv = rtot(k,i,j) / cvtot(k,i,j)
949 dens(k,i,j) = dens(k-1,i,j)
953 if ( abs(dens(k,i,j)-dens_s) <= criteria )
then 960 dhyd = + ( p00 * ( dens(k-1,i,j) * rtot(k-1,i,j) * pott(k-1,i,j) / p00 )**cpovcv(k-1,i,j) &
961 - p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) ) / dz(k,i,j) &
962 - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s )
964 dgrd = - p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k,i,j) &
965 * cpovcv(k,i,j) * dens_s**rovcv &
968 dens(k,i,j) = dens_s - dhyd/dgrd
970 if( dens(k,i,j)*0.0_rp /= 0.0_rp)
exit 973 if ( .NOT. converged )
then 974 write(*,*)
'xxx [buildrho 3D atmos] iteration not converged!', &
975 k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
985 pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
986 temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
1013 real(RP),
intent(out) :: dens_L1(
ia,
ja)
1014 real(RP),
intent(out) :: temp_L1(
ia,
ja)
1015 real(RP),
intent(out) :: pres_L1(
ia,
ja)
1016 real(RP),
intent(in) :: pott_L1(
ia,
ja)
1017 real(RP),
intent(in) :: qv_L1 (
ia,
ja)
1018 real(RP),
intent(in) :: qc_L1 (
ia,
ja)
1019 real(RP),
intent(in) :: dens_L2(
ia,
ja)
1020 real(RP),
intent(in) :: pott_L2(
ia,
ja)
1021 real(RP),
intent(in) :: qv_L2 (
ia,
ja)
1022 real(RP),
intent(in) :: qc_L2 (
ia,
ja)
1023 real(RP),
intent(in) :: dz (
ia,
ja)
1024 integer,
intent(in) :: k
1026 real(RP) :: Rtot_L1 (
ia,
ja), Rtot_L2 (
ia,
ja)
1027 real(RP) :: CVtot_L1 (
ia,
ja), CVtot_L2 (
ia,
ja)
1028 real(RP) :: CPovCV_L1(
ia,
ja), CPovCV_L2(
ia,
ja)
1031 real(RP) :: dens_s, dhyd, dgrd
1033 logical :: converged
1040 rtot_l1(i,j) = rdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1042 cvtot_l1(i,j) = cvdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1043 + cv_qv * qv_l1(i,j) &
1044 + cv_qc * qc_l1(i,j)
1045 cpovcv_l1(i,j) = ( cvtot_l1(i,j) + rtot_l1(i,j) ) / cvtot_l1(i,j)
1047 rtot_l2(i,j) = rdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1049 cvtot_l2(i,j) = cvdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1050 + cv_qv * qv_l2(i,j) &
1051 + cv_qc * qc_l2(i,j)
1052 cpovcv_l2(i,j) = ( cvtot_l2(i,j) + rtot_l2(i,j) ) / cvtot_l2(i,j)
1058 rovcv = rtot_l1(i,j) / cvtot_l1(i,j)
1061 dens_l1(i,j) = dens_l2(i,j)
1065 if ( abs(dens_l1(i,j)-dens_s) <= criteria )
then 1070 dens_s = dens_l1(i,j)
1072 dhyd = + ( p00 * ( dens_s * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) &
1073 - p00 * ( dens_l2(i,j) * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) ) / dz(i,j) &
1074 - grav * 0.5_rp * ( dens_s + dens_l2(i,j) )
1076 dgrd = + p00 * ( rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) / dz(i,j) &
1077 * cpovcv_l1(i,j) * dens_s**rovcv &
1080 dens_l1(i,j) = dens_s - dhyd/dgrd
1082 if( dens_l1(i,j)*0.0_rp /= 0.0_rp)
exit 1085 if ( .NOT. converged )
then 1086 write(*,*)
'xxx [buildrho 2D rev atmos] iteration not converged!', &
1087 i,j,k,dens_l1(i,j),ite,dens_s,dhyd,dgrd
1095 pres_l1(i,j) = p00 * ( dens_l1(i,j) * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j)
1096 temp_l1(i,j) = pres_l1(i,j) / ( dens_l1(i,j) * rtot_l1(i,j) )
1118 real(RP),
intent(inout) :: dens(
ka,
ia,
ja)
1119 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
1120 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
1121 real(RP),
intent(in) :: pott(
ka,
ia,
ja)
1122 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
1123 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
1124 real(RP),
intent(in) :: dz (
ka,
ia,
ja)
1125 integer,
intent(in),
optional :: kref_in
1127 real(RP) :: Rtot (
ka,
ia,
ja)
1128 real(RP) :: CVtot (
ka,
ia,
ja)
1129 real(RP) :: CPovCV(
ka,
ia,
ja)
1132 real(RP) :: dens_s, dhyd, dgrd
1134 logical :: converged
1140 if (
present(kref_in) )
then 1149 rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1151 cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1152 + cv_qv * qv(k,i,j) &
1154 cpovcv(k,i,j) = ( cvtot(k,i,j) + rtot(k,i,j) ) / cvtot(k,i,j)
1161 do k = kref-1,
ks, -1
1162 rovcv = rtot(k,i,j) / cvtot(k,i,j)
1165 dens(k,i,j) = dens(k+1,i,j)
1169 if ( abs(dens(k,i,j)-dens_s) <= criteria )
then 1174 dens_s = dens(k,i,j)
1176 dhyd = + ( p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) &
1177 - p00 * ( dens(k+1,i,j) * rtot(k+1,i,j) * pott(k+1,i,j) / p00 )**cpovcv(k+1,i,j) ) / dz(k+1,i,j) &
1178 - grav * 0.5_rp * ( dens_s + dens(k+1,i,j) )
1180 dgrd = + p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k+1,i,j) &
1181 * cpovcv(k,i,j) * dens_s**rovcv &
1184 dens(k,i,j) = dens_s - dhyd/dgrd
1186 if( dens(k,i,j)*0.0_rp /= 0.0_rp)
exit 1189 if ( .NOT. converged )
then 1190 write(*,*)
'xxx [buildrho 3D rev atmos] iteration not converged!', &
1191 k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1201 pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
1202 temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
1228 real(RP),
intent(out) :: dens(
ka)
1229 real(RP),
intent(out) :: pott(
ka)
1230 real(RP),
intent(out) :: pres(
ka)
1231 real(RP),
intent(in) :: temp(
ka)
1232 real(RP),
intent(in) :: qv (
ka)
1233 real(RP),
intent(in) :: qc (
ka)
1235 real(RP),
intent(out) :: pott_sfc
1236 real(RP),
intent(in) :: pres_sfc
1237 real(RP),
intent(in) :: temp_sfc
1238 real(RP),
intent(in) :: qv_sfc
1239 real(RP),
intent(in) :: qc_sfc
1241 real(RP) :: dens_sfc
1243 real(RP) :: Rtot_sfc
1244 real(RP) :: CVtot_sfc
1248 real(RP) :: RovCP_sfc
1249 real(RP) :: dens_s, dhyd, dgrd
1251 logical :: converged
1256 rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1258 cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1262 rtot = rdry * ( 1.0_rp - qv(
ks) - qc(
ks) ) &
1264 cvtot = cvdry * ( 1.0_rp - qv(
ks) - qc(
ks) ) &
1269 rovcp_sfc = rtot_sfc / ( cvtot_sfc + rtot_sfc )
1270 dens_sfc = pres_sfc / ( rtot_sfc * temp_sfc )
1271 pott_sfc = temp_sfc * ( p00/pres_sfc )**rovcp_sfc
1279 if ( abs(dens(
ks)-dens_s) <= criteria )
then 1286 dhyd = + ( dens_sfc * rtot_sfc * temp_sfc &
1288 - grav * 0.5_rp * ( dens_sfc + dens_s )
1293 dens(
ks) = dens_s - dhyd/dgrd
1295 if( dens(
ks)*0.0_rp /= 0.0_rp)
exit 1298 if ( .NOT. converged )
then 1299 write(*,*)
'xxx [buildrho bytemp 1D sfc] iteration not converged!', &
1300 dens(
ks),ite,dens_s,dhyd,dgrd
1333 real(RP),
intent(out) :: dens(
ka,
ia,
ja)
1334 real(RP),
intent(out) :: pott(
ka,
ia,
ja)
1335 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
1336 real(RP),
intent(in) :: temp(
ka,
ia,
ja)
1337 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
1338 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
1340 real(RP),
intent(out) :: pott_sfc(1,
ia,
ja)
1341 real(RP),
intent(in) :: pres_sfc(1,
ia,
ja)
1342 real(RP),
intent(in) :: temp_sfc(1,
ia,
ja)
1343 real(RP),
intent(in) :: qv_sfc (1,
ia,
ja)
1344 real(RP),
intent(in) :: qc_sfc (1,
ia,
ja)
1346 real(RP) :: dens_sfc (1,
ia,
ja)
1348 real(RP) :: Rtot_sfc (
ia,
ja)
1349 real(RP) :: CVtot_sfc (
ia,
ja)
1350 real(RP) :: Rtot (
ia,
ja)
1351 real(RP) :: CVtot (
ia,
ja)
1353 real(RP) :: RovCP_sfc
1355 real(RP) :: dens_s, dhyd, dgrd
1357 logical :: converged
1366 rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1367 + rvap * qv_sfc(1,i,j)
1368 cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1369 + cv_qv * qv_sfc(1,i,j) &
1370 + cv_qc * qc_sfc(1,i,j)
1376 rtot(i,j) = rdry * ( 1.0_rp - qv(
ks,i,j) - qc(
ks,i,j) ) &
1378 cvtot(i,j) = cvdry * ( 1.0_rp - qv(
ks,i,j) - qc(
ks,i,j) ) &
1379 + cv_qv * qv(
ks,i,j) &
1380 + cv_qc * qc(
ks,i,j)
1387 rovcp_sfc = rtot_sfc(i,j) / ( cvtot_sfc(i,j) + rtot_sfc(i,j) )
1388 dens_sfc(1,i,j) = pres_sfc(1,i,j) / ( rtot_sfc(i,j) * temp_sfc(1,i,j) )
1389 pott_sfc(1,i,j) = temp_sfc(1,i,j) / ( p00/pres_sfc(1,i,j) )**rovcp_sfc
1399 dens(
ks,i,j) = dens_sfc(1,i,j)
1403 if ( abs(dens(
ks,i,j)-dens_s) <= criteria )
then 1408 dens_s = dens(
ks,i,j)
1410 dhyd = + ( dens_sfc(1,i,j) * rtot_sfc(i,j) * temp_sfc(1,i,j) &
1411 - dens_s * rtot(i,j) * temp(
ks,i,j) ) / dz &
1412 - grav * 0.5_rp * ( dens_sfc(1,i,j) + dens_s )
1414 dgrd = - rtot(i,j) * temp(
ks,i,j) / dz &
1417 dens(
ks,i,j) = dens_s - dhyd/dgrd
1419 if( dens(
ks,i,j)*0.0_rp /= 0.0_rp)
exit 1422 if ( .NOT. converged )
then 1423 write(*,*)
'xxx [buildrho bytemp 3D sfc] iteration not converged!', &
1424 i,j,dens(
ks,i,j),ite,dens_s,dhyd,dgrd
1454 real(RP),
intent(inout) :: dens(
ka)
1455 real(RP),
intent(out) :: pott(
ka)
1456 real(RP),
intent(out) :: pres(
ka)
1457 real(RP),
intent(in) :: temp(
ka)
1458 real(RP),
intent(in) :: qv (
ka)
1459 real(RP),
intent(in) :: qc (
ka)
1461 real(RP) :: Rtot (
ka)
1462 real(RP) :: CVtot (
ka)
1465 real(RP) :: dens_s, dhyd, dgrd
1467 logical :: converged
1473 rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
1475 cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
1487 if ( abs(dens(k)-dens_s) <= criteria )
then 1494 dhyd = + ( dens(k-1) * rtot(k-1) * temp(k-1) &
1495 - dens_s * rtot(k ) * temp(k ) ) /
grid_fdz(k-1) &
1496 - grav * 0.5_rp * ( dens(k-1) + dens_s )
1498 dgrd = - rtot(k) * temp(k) /
grid_fdz(k-1) &
1501 dens(k) = dens_s - dhyd/dgrd
1503 if( dens(k)*0.0_rp /= 0.0_rp)
exit 1506 if ( .NOT. converged )
then 1507 write(*,*)
'xxx [buildrho bytemp 1D atmos] iteration not converged!', &
1508 k,dens(k),ite,dens_s,dhyd,dgrd
1514 rovcp = rtot(k) / ( cvtot(k) + rtot(k) )
1515 pres(k) = dens(k) * rtot(k) * temp(k)
1516 pott(k) = temp(k) * ( p00 / pres(k) )**rovcp
1535 real(RP),
intent(inout) :: dens(
ka,
ia,
ja)
1536 real(RP),
intent(out) :: pott(
ka,
ia,
ja)
1537 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
1538 real(RP),
intent(in) :: temp(
ka,
ia,
ja)
1539 real(RP),
intent(in) :: qv (
ka,
ia,
ja)
1540 real(RP),
intent(in) :: qc (
ka,
ia,
ja)
1542 real(RP) :: Rtot (
ka,
ia,
ja)
1543 real(RP) :: CVtot (
ka,
ia,
ja)
1547 real(RP) :: dens_s, dhyd, dgrd
1549 logical :: converged
1557 rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1559 cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1560 + cv_qv * qv(k,i,j) &
1572 dens(k,i,j) = dens(k-1,i,j)
1576 if ( abs(dens(k,i,j)-dens_s) <= criteria )
then 1581 dens_s = dens(k,i,j)
1583 dhyd = + ( dens(k-1,i,j) * rtot(k-1,i,j) * temp(k-1,i,j) &
1584 - dens_s * rtot(k ,i,j) * temp(k ,i,j) ) / dz &
1585 - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s )
1587 dgrd = - rtot(k,i,j) * temp(k,i,j) / dz &
1590 dens(k,i,j) = dens_s - dhyd/dgrd
1592 if( dens(k,i,j)*0.0_rp /= 0.0_rp)
exit 1595 if ( .NOT. converged )
then 1596 write(*,*)
'xxx [buildrho bytemp 3D atmos] iteration not converged!', &
1597 k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1607 rovcp = rtot(k,i,j) / ( cvtot(k,i,j) + rtot(k,i,j) )
1608 pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
1609 pott(k,i,j) = temp(k,i,j) * ( p00 / pres(k,i,j) )**rovcp
1626 real(RP),
intent(out) :: mslp
1627 real(RP),
intent(in) :: pres
1628 real(RP),
intent(in) :: temp
1629 real(RP),
intent(in) :: dz
1635 tm = temp + laps * dz * 0.5_rp
1638 mslp = pres * exp( grav * dz / ( rdry * tm ) )
1645 subroutine atmos_hydrostatic_barometric_law_mslp_2d( &
1652 real(RP),
intent(out) :: mslp(
ia,
ja)
1653 real(RP),
intent(in) :: pres(
ia,
ja)
1654 real(RP),
intent(in) :: temp(
ia,
ja)
1655 real(RP),
intent(in) :: dz (
ia,
ja)
1665 tm = temp(i,j) + laps * dz(i,j) * 0.5_rp
1668 mslp(i,j) = pres(i,j) * exp( grav * dz(i,j) / ( rdry * tm ) )
1673 end subroutine atmos_hydrostatic_barometric_law_mslp_2d
1677 subroutine atmos_hydrostatic_barometric_law_pres_0d( &
1684 real(RP),
intent(out) :: pres
1685 real(RP),
intent(in) :: mslp
1686 real(RP),
intent(in) :: temp
1687 real(RP),
intent(in) :: dz
1693 tm = temp + laps * dz * 0.5_rp
1696 pres = mslp / exp( grav * dz / ( rdry * tm ) )
1699 end subroutine atmos_hydrostatic_barometric_law_pres_0d
1703 subroutine atmos_hydrostatic_barometric_law_pres_2d( &
1710 real(RP),
intent(out) :: pres(
ia,
ja)
1711 real(RP),
intent(in) :: mslp(
ia,
ja)
1712 real(RP),
intent(in) :: temp(
ia,
ja)
1713 real(RP),
intent(in) :: dz (
ia,
ja)
1723 tm = temp(i,j) + laps * dz(i,j) * 0.5_rp
1726 pres(i,j) = mslp(i,j) / exp( grav * dz(i,j) / ( rdry * tm ) )
1731 end subroutine atmos_hydrostatic_barometric_law_pres_2d
subroutine atmos_hydrostatic_buildrho_3d(dens, temp, pres, pott, qv, qc, temp_sfc, pres_sfc, pott_sfc, qv_sfc, qc_sfc)
Build up density from surface (3D)
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_1d(dens, pott, pres, temp, qv, qc)
Build up density from lowermost atmosphere (1D)
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_hydrostatic_buildrho_atmos_rev_3d(dens, temp, pres, pott, qv, qc, dz, kref_in)
Build up density from lowermost atmosphere (3D)
logical, public io_l
output log or not? (this process)
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
real(rp), public const_laps
lapse rate of ISA [K/m]
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
subroutine atmos_hydrostatic_barometric_law_mslp_0d(mslp, pres, temp, dz)
Calculate mean sea-level pressure from barometric law (0D)
integer, public ka
of z whole cells (local, with HALO)
real(rp), public const_lapsdry
dry adiabatic lapse rate [K/m]
real(rp), public const_pre00
pressure reference [Pa]
module ATMOSPHERE / Hydrostatic barance
real(rp), public const_grav
standard acceleration of gravity [m/s2]
subroutine atmos_hydrostatic_buildrho_atmos_1d(dens, temp, pres, pott, qv, qc)
Build up density from lowermost atmosphere (1D)
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
subroutine atmos_hydrostatic_buildrho_1d(dens, temp, pres, pott, qv, qc, temp_sfc, pres_sfc, pott_sfc, qv_sfc, qc_sfc)
Build up density from surface (1D)
integer, public ks
start point of inner domain: z, local
subroutine atmos_hydrostatic_buildrho_atmos_3d(dens, temp, pres, pott, qv, qc, dz, kref_in)
Build up density from lowermost atmosphere (3D)
subroutine atmos_hydrostatic_buildrho_bytemp_3d(dens, pott, pres, temp, qv, qc, pott_sfc, pres_sfc, temp_sfc, qv_sfc, qc_sfc)
Build up density from surface (3D)
real(rp), public const_eps
small number
subroutine atmos_hydrostatic_buildrho_bytemp_1d(dens, pott, pres, temp, qv, qc, pott_sfc, pres_sfc, temp_sfc, qv_sfc, qc_sfc)
Build up density from surface (1D)
logical, public io_lnml
output log or not? (for namelist, this process)
subroutine atmos_hydrostatic_buildrho_real_3d(dens, temp, pres, pott, qv, qc)
Build up density from surface (3D), not to reverse from TOA.
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_3d(dens, pott, pres, temp, qv, qc)
Build up density from lowermost atmosphere (3D)
integer, public io_fid_conf
Config file ID.
integer, public io_fid_log
Log file ID.
character(len=h_short), public const_thermodyn_type
internal energy type
subroutine, public atmos_hydrostatic_buildrho_atmos_rev_2d(dens_L1, temp_L1, pres_L1, pott_L1, qv_L1, qc_L1, dens_L2, pott_L2, qv_L2, qc_L2, dz, k)
Build up density (2D)
subroutine, public atmos_hydrostatic_setup
Setup.
subroutine, public atmos_hydrostatic_buildrho_atmos_0d(dens_L2, temp_L2, pres_L2, pott_L2, qv_L2, qc_L2, dens_L1, pott_L1, qv_L1, qc_L1, dz, k)
Build up density (0D)
integer, public ja
of y whole cells (local, with HALO)