FV3 Bundle
tridiagonal.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
21 
22 ! <CONTACT EMAIL="Isaac.Held@noaa.gov">
23 ! Isaac Held
24 ! </CONTACT>
25 ! <CONTACT EMAIL="Bruce.Wyman@noaa.gov">
26 ! Bruce Wyman
27 ! </CONTACT>
28 
29 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
30 
31 ! <OVERVIEW>
32 ! Solves the tridiagonal system of equations.
33 ! </OVERVIEW>
34 ! <DESCRIPTION>
35 ! The following schematic represents the system of equations solved,
36 ! where X is the solution.
37 ! <PRE>
38 ! | B(1) A(1) 0 0 ....... 0 | |X(1)| |D(1)|
39 ! | C(2) B(2) A(2) 0 ....... 0 | |X(2)| |D(2)|
40 ! | 0 C(3) B(3) A(3) 0 ....... 0 | | .. | | .. |
41 ! | .......................................... | | .. | = | .. |
42 ! | .......................................... | | .. | | .. |
43 ! | C(N-2) B(N-2) A(N-2) 0 | | .. | | .. |
44 ! | 0 C(N-1) B(N-1) A(N-1)| | .. | | .. |
45 ! | 0 0 C(N) B(N) | |X(N)| |D(N)|
46 !
47 ! </PRE>
48 ! To solve this system
49 ! <PRE>
50 ! call tri_invert(X,D,A,B,C)
51 !
52 ! real, intent(out), dimension(:,:,:) :: X
53 ! real, intent(in), dimension(:,:,:) :: D
54 ! real, optional, dimension(:,:,:) :: A,B,C
55 ! </PRE>
56 ! For simplicity (?), A and C are assumed to be dimensioned the same size
57 ! as B, D, and X, although any input values for A(N) and C(1) are ignored.
58 ! (some checks are needed here)
59 !
60 ! If A is not present, it is assumed that the matrix (A,B.C) has not been changed
61 ! since the last call to tri_invert.
62 !
63 ! To release memory,
64 ! <PRE>
65 ! call close_tridiagonal
66 ! </PRE>
67 ! The following module variables are used to retain the relevant information
68 ! if one recalls tri_invert without changing (A,B,C)
69 
70 
71 ! </DESCRIPTION>
72 
73 !--------------------------------------------------------------------------
74 real, private, allocatable, dimension(:,:,:) :: e,g,cc
75 real, private, allocatable, dimension(:,:) :: bb
76 logical, private :: init_tridiagonal = .false.
77 !--------------------------------------------------------------------------
78 
79 contains
80 
81 !--------------------------------------------------------------------------
82 
83 ! <SUBROUTINE NAME="tri_invert">
84 
85 ! <OVERVIEW>
86 ! Sets up and solves the tridiagonal system of equations.
87 ! </OVERVIEW>
88 ! <DESCRIPTION>
89 ! Sets up and solves the tridiagonal system of equations.
90 ! </DESCRIPTION>
91 ! <TEMPLATE>
92 ! call tri_invert ( x,d,a,b,c)
93 ! </TEMPLATE>
94 
95 ! <IN NAME="d" TYPE="real" DIM="(:,:,:)">
96 ! The right-hand side term, see the schematic above.
97 ! </IN>
98 ! <OUT NAME="x" TYPE="real" DIM="(:,:,:)">
99 ! The solution to the tridiagonal system of equations.
100 ! </OUT>
101 ! <INOUT NAME="a,b,c" TYPE="real" DIM="(:,:,:)">
102 ! The left-hand-side terms (matrix), see the schematic above.
103 ! If A is not present, it is assumed that the matrix (A,B.C)
104 ! has not been changed since the last call to tri_invert.
105 ! </INOUT>
106 
107 ! <NOTE>
108 ! For simplicity, A and C are assumed to be dimensioned the same size
109 ! as B, D, and X, although any input values for A(N) and C(1) are ignored.
110 ! There are no checks to make sure the sizes agree.
111 ! </NOTE>
112 ! <NOTE>
113 ! The value of A(N) is modified on output, and B and C are unchanged.
114 ! </NOTE>
115 
116 subroutine tri_invert(x,d,a,b,c)
118 implicit none
119 
120 real, intent(out), dimension(:,:,:) :: x
121 real, intent(in), dimension(:,:,:) :: d
122 real, optional, dimension(:,:,:) :: a,b,c
123 
124 real, dimension(size(x,1),size(x,2),size(x,3)) :: f
125 integer :: k
126 
127 if(present(a)) then
128  init_tridiagonal = .true.
129 
130  if(allocated(e)) deallocate(e)
131  if(allocated(g)) deallocate(g)
132  if(allocated(bb)) deallocate(bb)
133  if(allocated(cc)) deallocate(cc)
134  allocate(e(size(x,1),size(x,2),size(x,3)))
135  allocate(g(size(x,1),size(x,2),size(x,3)))
136  allocate(bb(size(x,1),size(x,2)))
137  allocate(cc(size(x,1),size(x,2),size(x,3)))
138 
139  e(:,:,1) = - a(:,:,1)/b(:,:,1)
140  a(:,:,size(x,3)) = 0.0
141 
142  do k= 2,size(x,3)
143  g(:,:,k) = 1.0/(b(:,:,k)+c(:,:,k)*e(:,:,k-1))
144  e(:,:,k) = - a(:,:,k)*g(:,:,k)
145  end do
146  cc = c
147  bb = 1.0/b(:,:,1)
148 
149 end if
150 
151 ! if(.not.init_tridiagonal) error
152 
153 f(:,:,1) = d(:,:,1)*bb
154 do k= 2, size(x,3)
155  f(:,:,k) = (d(:,:,k) - cc(:,:,k)*f(:,:,k-1))*g(:,:,k)
156 end do
157 
158 x(:,:,size(x,3)) = f(:,:,size(x,3))
159 do k = size(x,3)-1,1,-1
160  x(:,:,k) = e(:,:,k)*x(:,:,k+1)+f(:,:,k)
161 end do
162 
163 return
164 end subroutine tri_invert
165 ! </SUBROUTINE>
166 
167 !-----------------------------------------------------------------
168 
169 ! <SUBROUTINE NAME="close_tridiagonal">
170 
171 ! <OVERVIEW>
172 ! Releases memory used by the solver.
173 ! </OVERVIEW>
174 ! <DESCRIPTION>
175 ! Releases memory used by the solver.
176 ! </DESCRIPTION>
177 ! <TEMPLATE>
178 ! call close_tridiagonal
179 ! </TEMPLATE>
180 
181 ! <NOTE>
182 ! There are no arguments to this routine.
183 ! </NOTE>
184 
185 subroutine close_tridiagonal
187 implicit none
188 
189 deallocate(e)
190 deallocate(g)
191 deallocate(bb)
192 deallocate(cc)
193 
194 return
195 end subroutine close_tridiagonal
196 ! </SUBROUTINE>
197 
198 !----------------------------------------------------------------
199 
200 
201 
202 end module tridiagonal_mod
203 
204 ! <INFO>
205 
206 ! <BUG>
207 ! Optional arguments A,B,C have no intent declaration,
208 ! so the default intent is inout. The value of A(N) is modified
209 ! on output, and B and C are unchanged.
210 ! </BUG>
211 ! <NOTE>
212 ! The following private allocatable arrays save the relevant information
213 ! if one recalls tri_invert without changing (A,B,C):
214 ! <PRE>
215 ! allocate ( e (size(x,1), size(x,2), size(x,3)) )
216 ! allocate ( g (size(x,1), size(x,2), size(x,3)) )
217 ! allocate ( cc (size(x,1), size(x,2), size(x,3)) )
218 ! allocate ( bb (size(x,1), size(x,2)) )
219 ! </PRE>
220 ! This storage is deallocated when close_tridiagonal is called.
221 ! </NOTE>
222 ! <FUTURE>
223 ! Maybe a cleaner version?
224 ! </FUTURE>
225 
226 ! </INFO>
subroutine close_tridiagonal
real, dimension(:,:), allocatable, private bb
Definition: tridiagonal.F90:75
real, dimension(:,:,:), allocatable, private e
Definition: tridiagonal.F90:74
logical, private init_tridiagonal
Definition: tridiagonal.F90:76
subroutine tri_invert(x, d, a, b, c)
real, dimension(:,:,:), allocatable, private g
Definition: tridiagonal.F90:74
real, dimension(:,:,:), allocatable, private cc
Definition: tridiagonal.F90:74