CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
h2ophys.f
1
3
5 module h2ophys
6
7 implicit none
8
9 private
10
11 public :: h2ophys_init, h2ophys_run
12
13 contains
14
15 subroutine h2ophys_init(h2o_phys, errmsg, errflg)
16
17 implicit none
18 logical, intent(in) :: h2o_phys
19 character(len=*), intent(out) :: errmsg
20 integer, intent(out) :: errflg
21
22 ! Initialize CCPP error handling variables
23 errmsg = ''
24 errflg = 0
25
26 if (.not.h2o_phys) then
27 write (errmsg,'(*(a))') 'Logic error: h2o_phys == .false.'
28 errflg = 1
29 return
30 endif
31 end subroutine h2ophys_init
32
40 subroutine h2ophys_run(im, levs, kh2o, dt, h2o, ph2o, prsl, &
41 & h2opltc, h2o_coeff, me, &
42 & errmsg, errflg)
43!
44! May 2015 - Shrinivas Moorthi - Adaptation of NRL H2O physics for
45! stratosphere and mesosphere
46!
47! this code assumes that both prsl and ph2o are from bottom to top
48! as are all other variables
49!
50 use machine , only : kind_phys
51 implicit none
52! interface variables
53 integer, intent(in) :: im, levs, kh2o, h2o_coeff, me
54 real(kind=kind_phys), intent(in) :: dt
55 real(kind=kind_phys), intent(inout) :: h2o(:,:)
56 real(kind=kind_phys), intent(in) :: ph2o(:)
57 real(kind=kind_phys), intent(in) :: prsl(:,:)
58 real(kind=kind_phys), intent(in) :: h2opltc(:,:,:)
59 !real(kind=kind_phys), intent(inout) :: h2op(im,levs,h2o_coeff)
60 character(len=*), intent(out) :: errmsg
61 integer, intent(out) :: errflg
62! local variables
63 integer k,kmax,kmin,l,i,j
64 logical flg(im)
65 real(kind=kind_phys) pmax, pmin, tem, temp
66 real(kind=kind_phys) wk1(im), wk2(im), wk3(im), pltc(im,h2o_coeff)
67 &, h2oib(im)
68 real, parameter :: prsmax=10000.0, pmaxl=log(prsmax)
69!
70! initialize CCPP error handling variables
71 errmsg = ''
72 errflg = 0
73!
74! write(1000+me,*)' in h2ophys im=', im, levs, kh2o, dt
75 do l=1,levs
76 pmin = 1.0e10
77 pmax = -1.0e10
78!
79 do i=1,im
80 wk1(i) = log(prsl(i,l))
81 pmin = min(wk1(i), pmin)
82 pmax = max(wk1(i), pmax)
83 pltc(i,:) = 0.0
84 enddo
85 if (pmin < pmaxl) then
86 kmax = 1
87 kmin = 1
88 do k=1,kh2o-1
89 if (pmin < ph2o(k)) kmax = k
90 if (pmax < ph2o(k)) kmin = k
91 enddo
92!
93 do k=kmin,kmax
94 temp = 1.0 / (ph2o(k) - ph2o(k+1))
95 do i=1,im
96 flg(i) = .false.
97 if (wk1(i) < ph2o(k) .and. wk1(i) >= ph2o(k+1)) then
98 flg(i) = .true.
99 wk2(i) = (wk1(i) - ph2o(k+1)) * temp
100 wk3(i) = 1.0 - wk2(i)
101 endif
102 enddo
103 do j=1,h2o_coeff
104 do i=1,im
105 if (flg(i)) then
106 pltc(i,j) = wk2(i) * h2opltc(i,k,j)
107 & + wk3(i) * h2opltc(i,k+1,j)
108 endif
109 enddo
110 enddo
111 enddo
112!
113 do j=1,h2o_coeff
114 do i=1,im
115 if (wk1(i) < ph2o(kh2o)) then
116 pltc(i,j) = h2opltc(i,kh2o,j)
117 endif
118 if (wk1(i) >= ph2o(1)) then
119 pltc(i,j) = h2opltc(i,1,j)
120 endif
121 enddo
122 enddo
123 endif
124 do i=1,im
125 if (prsl(i,l) < prsmax .and. pltc(i,2) /= 0.0) then
126 h2oib(i) = h2o(i,l) ! no filling
127 tem = 1.0 / pltc(i,2) ! 1/teff
128 h2o(i,l) = (h2oib(i) + (pltc(i,1)+pltc(i,3)*tem)*dt)
129 & / (1.0 + tem*dt)
130 endif
131
132! if (i == 1) write(1000+me,*)' h2oib=',h2oib(i),' pltc1=',
133! &pltc(i,1),' pltc2=', pltc(i,2),' tem=',tem ,' dt=',dt
134! &,' l=',l
135 enddo
136!
137! if (ldiag3d) then ! h2o change diagnostics
138! do i=1,im
139! h2op(i,l,1) = h2op(i,l,1) + pltc(i,1)*dt
140! h2op(i,l,2) = h2op(i,l,2) + (h2oo(i,l) - h2oib(i))
141! enddo
142! endif
143 enddo ! vertical loop
144!
145 return
146 end subroutine h2ophys_run
148
149 end module h2ophys
subroutine, public h2ophys_run(im, levs, kh2o, dt, h2o, ph2o, prsl, h2opltc, h2o_coeff, me, errmsg, errflg)
Definition h2ophys.f:43
This module contains the CCPP-compliant H2O physics for stratosphere and mesosphere.
Definition h2ophys.f:5