CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_nst_pre.f90
1
3
6
7 use machine , only : kind_phys
9 use module_nst_parameters , only : zero, one
10
11 implicit none
12
13contains
14
25 subroutine sfc_nst_pre_run &
26 (im, wet, tgice, tsfco, tsurf_wat, &
27 tseal, xt, xz, dt_cool, z_c, tref, cplflx, &
28 oceanfrac, nthreads, errmsg, errflg)
29
30 ! --- inputs:
31 integer, intent(in) :: im, nthreads
32 logical, dimension(:), intent(in) :: wet
33 real (kind=kind_phys), intent(in) :: tgice
34 real (kind=kind_phys), dimension(:), intent(in) :: tsfco, oceanfrac
35 real (kind=kind_phys), dimension(:), intent(in), optional :: xt, xz, dt_cool, z_c
36 logical, intent(in) :: cplflx
37
38 ! --- input/outputs:
39 real (kind=kind_phys), dimension(:), intent(inout) :: tsurf_wat, tseal
40 real (kind=kind_phys), dimension(:), intent(inout), optional :: tref
41
42 ! --- outputs:
43 character(len=*), intent(out) :: errmsg
44 integer, intent(out) :: errflg
45
46 ! --- locals
47 integer :: i
48 real(kind=kind_phys), parameter :: omz1 = 2.0_kind_phys
49 real(kind=kind_phys) :: tem2, dnsst
50 real(kind=kind_phys), dimension(im) :: dtzm, z_c_0
51
52 ! Initialize CCPP error handling variables
53 errmsg = ''
54 errflg = 0
55
56 do i=1,im
57 if (wet(i) .and. oceanfrac(i) > 0.0) then
58 ! tem = (oro(i)-oro_uf(i)) * rlapse
59 ! DH* 20190927 simplyfing this code because tem is zero
60 !tem = zero
61 !tseal(i) = tsfco(i) + tem
62 tseal(i) = tsfco(i)
63 !tsurf_wat(i) = tsurf_wat(i) + tem
64 ! *DH
65 endif
66 enddo
67 !
68 ! update tsfc & tref with T1 from OGCM & NSST Profile if coupled
69 !
70 if (cplflx) then
71 z_c_0 = zero
72 call get_dtzm_2d (xt, xz, dt_cool, z_c_0, wet, zero, omz1, im, 1, nthreads, dtzm)
73 do i=1,im
74 if (wet(i) .and. oceanfrac(i) > zero ) then
75 ! dnsst = tsfc_wat(i) - tref(i) ! retrive/get difference of Ts and Tf
76 tref(i) = max(tgice, tsfco(i) - dtzm(i)) ! update Tf with T1 and NSST T-Profile
77 ! tsfc_wat(i) = max(271.2,tref(i) + dnsst) ! get Ts updated due to Tf update
78 ! tseal(i) = tsfc_wat(i)
79 if (abs(xz(i)) > zero) then
80 tem2 = one / xz(i)
81 else
82 tem2 = zero
83 endif
84 tseal(i) = tref(i) + (xt(i)+xt(i)) * tem2 - dt_cool(i)
85 tsurf_wat(i) = tseal(i)
86 endif
87 enddo
88 endif
89
90 return
91 end subroutine sfc_nst_pre_run
92end module sfc_nst_pre
subroutine, public get_dtzm_2d(xt, xz, dt_cool, zc, wet, z1, z2, nx, ny, nth, dtm)
This module contains constants and parameters used in GFS near surface sea temperature scheme.
This module contains GFS NSST water property subroutines.
This module contain preparation for the near-surface sea temperature scheme.