CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
tridi.f
1
3
6 module tridi_mod
7 contains
8
13 subroutine tridi1(l,n,cl,cm,cu,r1,au,a1)
14 !
15 use machine , only : kind_phys
16 implicit none
17 integer, parameter :: one = 1.0_kind_phys
18 integer k,n,l,i
19 real(kind=kind_phys) fk
20 !
21 real(kind=kind_phys) cl(l,2:n),cm(l,n),cu(l,n-1),r1(l,n), &
22 & au(l,n-1),a1(l,n)
23 !
24 do i=1,l
25 fk = one / cm(i,1)
26 au(i,1) = fk*cu(i,1)
27 a1(i,1) = fk*r1(i,1)
28 enddo
29 do k=2,n-1
30 do i=1,l
31 fk = one / (cm(i,k)-cl(i,k)*au(i,k-1))
32 au(i,k) = fk*cu(i,k)
33 a1(i,k) = fk*(r1(i,k)-cl(i,k)*a1(i,k-1))
34 enddo
35 enddo
36 do i=1,l
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))
39 enddo
40 do k=n-1,1,-1
41 do i=1,l
42 a1(i,k) = a1(i,k)-au(i,k)*a1(i,k+1)
43 enddo
44 enddo
45 !
46 return
47 end subroutine tridi1
48
49!-----------------------------------------------------------------------
52 subroutine tridi2(l,n,cl,cm,cu,r1,r2,au,a1,a2)
53!
54 use machine , only : kind_phys
55 implicit none
56 integer, parameter :: one = 1.0_kind_phys
57 integer k,n,l,i
58 real(kind=kind_phys) fk
59!
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)
62!----------------------------------------------------------------------
63 do i=1,l
64 fk = one / cm(i,1)
65 au(i,1) = fk*cu(i,1)
66 a1(i,1) = fk*r1(i,1)
67 a2(i,1) = fk*r2(i,1)
68 enddo
69 do k=2,n-1
70 do i=1,l
71 fk = one / (cm(i,k)-cl(i,k)*au(i,k-1))
72 au(i,k) = fk*cu(i,k)
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))
75 enddo
76 enddo
77 do i=1,l
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))
81 enddo
82 do k=n-1,1,-1
83 do i=1,l
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)
86 enddo
87 enddo
88!-----------------------------------------------------------------------
89 return
90 end subroutine tridi2
91
92!-----------------------------------------------------------------------
97 subroutine tridin(l,n,nt,cl,cm,cu,r1,r2,au,a1,a2)
98!
99 use machine , only : kind_phys
100 implicit none
101 integer, parameter :: one = 1.0_kind_phys
102 integer is,k,kk,n,nt,l,i
103 real(kind=kind_phys) fk(l)
104!
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), &
108 & fkk(l,2:n-1)
109!-----------------------------------------------------------------------
110 do i=1,l
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)
114 enddo
115 do k = 1, nt
116 is = (k-1) * n
117 do i = 1, l
118 a2(i,1+is) = fk(i) * r2(i,1+is)
119 enddo
120 enddo
121 do k=2,n-1
122 do i=1,l
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))
126 enddo
127 enddo
128 do kk = 1, nt
129 is = (kk-1) * n
130 do k=2,n-1
131 do i=1,l
132 a2(i,k+is) = fkk(i,k)*(r2(i,k+is)-cl(i,k)*a2(i,k+is-1))
133 enddo
134 enddo
135 enddo
136 do i=1,l
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))
139 enddo
140 do k = 1, nt
141 is = (k-1) * n
142 do i = 1, l
143 a2(i,n+is) = fk(i)*(r2(i,n+is)-cl(i,n)*a2(i,n+is-1))
144 enddo
145 enddo
146 do k=n-1,1,-1
147 do i=1,l
148 a1(i,k) = a1(i,k) - au(i,k)*a1(i,k+1)
149 enddo
150 enddo
151 do kk = 1, nt
152 is = (kk-1) * n
153 do k=n-1,1,-1
154 do i=1,l
155 a2(i,k+is) = a2(i,k+is) - au(i,k)*a2(i,k+is+1)
156 enddo
157 enddo
158 enddo
159!-----------------------------------------------------------------------
160 return
161 end subroutine tridin
162
163!-----------------------------------------------------------------------
166 subroutine tridit(l,n,nt,cl,cm,cu,rt,au,at)
167!-----------------------------------------------------------------------
168!!
169 use machine , only : kind_phys
170 implicit none
171 integer, parameter :: one = 1.0_kind_phys
172 integer is,k,kk,n,nt,l,i
173 real(kind=kind_phys) fk(l)
174!!
175 real(kind=kind_phys) cl(l,2:n), cm(l,n), cu(l,n-1), &
176 & rt(l,n*nt), &
177 & au(l,n-1), at(l,n*nt), &
178 & fkk(l,2:n-1)
179!-----------------------------------------------------------------------
180 do i=1,l
181 fk(i) = one / cm(i,1)
182 au(i,1) = fk(i)*cu(i,1)
183 enddo
184 do k = 1, nt
185 is = (k-1) * n
186 do i = 1, l
187 at(i,1+is) = fk(i) * rt(i,1+is)
188 enddo
189 enddo
190 do k=2,n-1
191 do i=1,l
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)
194 enddo
195 enddo
196 do kk = 1, nt
197 is = (kk-1) * n
198 do k=2,n-1
199 do i=1,l
200 at(i,k+is) = fkk(i,k)*(rt(i,k+is)-cl(i,k)*at(i,k+is-1))
201 enddo
202 enddo
203 enddo
204 do i=1,l
205 fk(i) = one / (cm(i,n)-cl(i,n)*au(i,n-1))
206 enddo
207 do k = 1, nt
208 is = (k-1) * n
209 do i = 1, l
210 at(i,n+is) = fk(i)*(rt(i,n+is)-cl(i,n)*at(i,n+is-1))
211 enddo
212 enddo
213 do kk = 1, nt
214 is = (kk-1) * n
215 do k=n-1,1,-1
216 do i=1,l
217 at(i,k+is) = at(i,k+is) - au(i,k)*at(i,k+is+1)
218 enddo
219 enddo
220 enddo
221!-----------------------------------------------------------------------
222 return
223 end subroutine tridit
224 end module tridi_mod
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...
Definition tridi.f:14
subroutine tridit(l, n, nt, cl, cm, cu, rt, au, at)
This subroutine solves tridiagonal problem for TKE.
Definition tridi.f:167
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 ...
Definition tridi.f:98
subroutine tridi2(l, n, cl, cm, cu, r1, r2, au, a1, a2)
This subroutine ..
Definition tridi.f:53
This module contains routine to compute tridiagonal matrix elements for TKE, heat,...
Definition tridi.f:6