CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_diag_post.F90
1
3
6
7 contains
8
11#if 0
12
15#endif
16 subroutine sfc_diag_post_run (im, lsm, lsm_noahmp, opt_diag, dry, lssav, dtf, con_eps, con_epsm1, pgr,&
17 vegtype,t2mmp,q2mp, t2m, q2m, u10m, v10m, tmpmin, tmpmax, spfhmin, spfhmax, &
18 wind10mmax, u10mmax, v10mmax, dpt2m, errmsg, errflg)
19
20 use machine, only: kind_phys, kind_dbl_prec
21
22 implicit none
23
24 integer, intent(in) :: im, lsm, lsm_noahmp,opt_diag
25 integer, dimension(:), intent(in) :: vegtype ! vegetation type (integer index)
26 logical, intent(in) :: lssav
27 real(kind=kind_phys), intent(in) :: dtf, con_eps, con_epsm1
28 logical , dimension(:), intent(in) :: dry
29 real(kind=kind_phys), dimension(:), intent(in) :: pgr, u10m, v10m
30 real(kind=kind_phys), dimension(:), intent(inout) :: t2m, q2m, tmpmin, tmpmax, spfhmin, spfhmax
31 real(kind=kind_phys), dimension(:), intent(in), optional :: t2mmp, q2mp
32 real(kind=kind_phys), dimension(:), intent(inout) :: wind10mmax, u10mmax, v10mmax, dpt2m
33
34 character(len=*), intent(out) :: errmsg
35 integer, intent(out) :: errflg
36
37 integer :: i
38 real(kind=kind_dbl_prec) :: tem ! made dbl prec always, JM 20211104
39
40 ! Initialize CCPP error handling variables
41 errmsg = ''
42 errflg = 0
43
44 if (lsm == lsm_noahmp) then
45! over shrublands use opt_diag=2
46 do i=1, im
47 if(dry(i)) then
48 if (vegtype(i) == 6 .or. vegtype(i) == 7 &
49 .or. vegtype(i) == 16) then
50 t2m(i) = t2mmp(i)
51 q2m(i) = q2mp(i)
52 endif
53 endif
54 enddo
55
56 if (opt_diag == 2 .or. opt_diag == 3) then
57 do i=1,im
58 if(dry(i)) then
59 t2m(i) = t2mmp(i)
60 q2m(i) = q2mp(i)
61 endif
62 enddo
63 endif
64 endif
65
66 if (lssav) then
67 do i=1,im
68 tmpmax(i) = max(tmpmax(i),t2m(i))
69 tmpmin(i) = min(tmpmin(i),t2m(i))
70 spfhmax(i) = max(spfhmax(i),q2m(i))
71 spfhmin(i) = min(spfhmin(i),q2m(i))
72 enddo
73 ! Find max wind speed then decompose
74 do i=1, im
75 tem = sqrt(u10m(i)*u10m(i) + v10m(i)*v10m(i))
76 if (tem > wind10mmax(i)) then
77 wind10mmax(i) = tem
78 u10mmax(i) = u10m(i)
79 v10mmax(i) = v10m(i)
80 endif
81 ! Compute dew point, first using vapor pressure
82 tem = max(pgr(i) * q2m(i) / ( con_eps - con_epsm1 *q2m(i)), 1.d-8)
83 dpt2m(i) = 243.5_kind_dbl_prec / &
84 ( ( 17.67_kind_dbl_prec / log(tem/611.2_kind_dbl_prec) ) - 1.) + 273.14
85 enddo
86 endif
87
88 end subroutine sfc_diag_post_run
90 end module sfc_diag_post
subroutine sfc_diag_post_run(im, lsm, lsm_noahmp, opt_diag, dry, lssav, dtf, con_eps, con_epsm1, pgr, vegtype, t2mmp, q2mp, t2m, q2m, u10m, v10m, tmpmin, tmpmax, spfhmin, spfhmax, wind10mmax, u10mmax, v10mmax, dpt2m, errmsg, errflg)
This module contains code related to the surface diagnostic scheme.