13 subroutine tridi1(l,n,cl,cm,cu,r1,au,a1)
15 use machine ,
only : kind_phys
17 integer,
parameter :: one = 1.0_kind_phys
19 real(kind=kind_phys) fk
21 real(kind=kind_phys) cl(l,2:n),cm(l,n),cu(l,n-1),r1(l,n), &
31 fk = one / (cm(i,k)-cl(i,k)*au(i,k-1))
33 a1(i,k) = fk*(r1(i,k)-cl(i,k)*a1(i,k-1))
37 fk = one / (cm(i,n)-cl(i,n)*au(i,n-1))
38 a1(i,n) = fk*(r1(i,n)-cl(i,n)*a1(i,n-1))
42 a1(i,k) = a1(i,k)-au(i,k)*a1(i,k+1)
52 subroutine tridi2(l,n,cl,cm,cu,r1,r2,au,a1,a2)
54 use machine ,
only : kind_phys
56 integer,
parameter :: one = 1.0_kind_phys
58 real(kind=kind_phys) fk
60 real(kind=kind_phys) cl(l,2:n),cm(l,n),cu(l,n-1),r1(l,n),r2(l,n), &
61 & au(l,n-1),a1(l,n),a2(l,n)
71 fk = one / (cm(i,k)-cl(i,k)*au(i,k-1))
73 a1(i,k) = fk*(r1(i,k)-cl(i,k)*a1(i,k-1))
74 a2(i,k) = fk*(r2(i,k)-cl(i,k)*a2(i,k-1))
78 fk = one / (cm(i,n)-cl(i,n)*au(i,n-1))
79 a1(i,n) = fk*(r1(i,n)-cl(i,n)*a1(i,n-1))
80 a2(i,n) = fk*(r2(i,n)-cl(i,n)*a2(i,n-1))
84 a1(i,k) = a1(i,k)-au(i,k)*a1(i,k+1)
85 a2(i,k) = a2(i,k)-au(i,k)*a2(i,k+1)
97 subroutine tridin(l,n,nt,cl,cm,cu,r1,r2,au,a1,a2)
99 use machine ,
only : kind_phys
101 integer,
parameter :: one = 1.0_kind_phys
102 integer is,k,kk,n,nt,l,i
103 real(kind=kind_phys) fk(l)
105 real(kind=kind_phys) cl(l,2:n), cm(l,n), cu(l,n-1), &
106 & r1(l,n), r2(l,n*nt), &
107 & au(l,n-1), a1(l,n), a2(l,n*nt), &
111 fk(i) = one / cm(i,1)
112 au(i,1) = fk(i)*cu(i,1)
113 a1(i,1) = fk(i)*r1(i,1)
118 a2(i,1+is) = fk(i) * r2(i,1+is)
123 fkk(i,k) = one / (cm(i,k)-cl(i,k)*au(i,k-1))
124 au(i,k) = fkk(i,k)*cu(i,k)
125 a1(i,k) = fkk(i,k)*(r1(i,k)-cl(i,k)*a1(i,k-1))
132 a2(i,k+is) = fkk(i,k)*(r2(i,k+is)-cl(i,k)*a2(i,k+is-1))
137 fk(i) = one / (cm(i,n)-cl(i,n)*au(i,n-1))
138 a1(i,n) = fk(i)*(r1(i,n)-cl(i,n)*a1(i,n-1))
143 a2(i,n+is) = fk(i)*(r2(i,n+is)-cl(i,n)*a2(i,n+is-1))
148 a1(i,k) = a1(i,k) - au(i,k)*a1(i,k+1)
155 a2(i,k+is) = a2(i,k+is) - au(i,k)*a2(i,k+is+1)
166 subroutine tridit(l,n,nt,cl,cm,cu,rt,au,at)
169 use machine ,
only : kind_phys
171 integer,
parameter :: one = 1.0_kind_phys
172 integer is,k,kk,n,nt,l,i
173 real(kind=kind_phys) fk(l)
175 real(kind=kind_phys) cl(l,2:n), cm(l,n), cu(l,n-1), &
177 & au(l,n-1), at(l,n*nt), &
181 fk(i) = one / cm(i,1)
182 au(i,1) = fk(i)*cu(i,1)
187 at(i,1+is) = fk(i) * rt(i,1+is)
192 fkk(i,k) = one / (cm(i,k)-cl(i,k)*au(i,k-1))
193 au(i,k) = fkk(i,k)*cu(i,k)
200 at(i,k+is) = fkk(i,k)*(rt(i,k+is)-cl(i,k)*at(i,k+is-1))
205 fk(i) = one / (cm(i,n)-cl(i,n)*au(i,n-1))
210 at(i,n+is) = fk(i)*(rt(i,n+is)-cl(i,n)*at(i,n+is-1))
217 at(i,k+is) = at(i,k+is) - au(i,k)*at(i,k+is+1)
subroutine tridi1(l, n, cl, cm, cu, r1, au, a1)
Routine to solve the tridiagonal system to calculate temperature and moisture at ; part of two-part p...
subroutine tridit(l, n, nt, cl, cm, cu, rt, au, at)
This subroutine solves tridiagonal problem for TKE.
subroutine tridin(l, n, nt, cl, cm, cu, r1, r2, au, a1, a2)
Routine to solve the tridiagonal system to calculate u- and v-momentum at ; part of two-part process ...
subroutine tridi2(l, n, cl, cm, cu, r1, r2, au, a1, a2)
This subroutine ..