SABER
tools_asa007.F90
Go to the documentation of this file.
1 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/tools_asa007.fypp"
2 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/../instrumentation.fypp" 1
3 # 1 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/../subr_list.fypp" 1
4 !----------------------------------------------------------------------
5 ! Header: subr_list
6 !> Subroutines/functions list
7 ! Author: Benjamin Menetrier
8 ! Licensing: this code is distributed under the CeCILL-C license
9 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
10 !----------------------------------------------------------------------
11 
12 # 726 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/../subr_list.fypp"
13 # 2 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/../instrumentation.fypp" 2
14 !----------------------------------------------------------------------
15 ! Header: instrumentation
16 !> Instrumentation functions
17 ! Author: Benjamin Menetrier
18 ! Licensing: this code is distributed under the CeCILL-C license
19 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
20 !----------------------------------------------------------------------
21 
22 # 112 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/../instrumentation.fypp"
23 # 2 "/Users/miesch/JEDI/code/working_copy/public/fv3-bundle/saber/src/saber/external/tools_asa007.fypp" 2
24 !----------------------------------------------------------------------
25 ! Module: tools_asa007
26 !> Inverse of symmetric positive definite matrix routines
27 ! Source: https://people.sc.fsu.edu/~jburkardt/f_src/asa007/asa007.html
28 ! Author: Michael Healy
29 ! Original licensing: none
30 ! Modified by Alan Miller
31 ! Fortran 90 version by John Burkardt
32 ! Modified by Benjamin Menetrier for BUMP
33 ! Licensing: this code is distributed under the CeCILL-C license
34 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
35 !----------------------------------------------------------------------
37 
38 use tools_const, only: zero,one
39 use tools_kinds, only: kind_real
40 use tools_repro, only: inf,infeq
41 use type_mpl, only: mpl_type
42 
43 
44 implicit none
45 
46 real(kind_real),parameter :: eta = 1.0e-9_kind_real !< Small parameter for the Cholesky decomposition
47 
48 interface cholesky
49  module procedure asa007_cholesky
50 end interface
51 interface syminv
52  module procedure asa007_syminv
53 end interface
54 
55 private
56 public :: cholesky,syminv
57 
58 contains
59 
60 !----------------------------------------------------------------------
61 ! Subroutine: asa007_cholesky
62 !> Compute cholesky decomposition
63 !----------------------------------------------------------------------
64 subroutine asa007_cholesky(mpl,n,nn,a,u)
65 
66 implicit none
67 
68 ! Passed variables
69 type(mpl_type),intent(inout) :: mpl !< MPI data
70 integer,intent(in) :: n !< Matrix rank
71 integer,intent(in) :: nn !< Half-matrix size (n*(n-1)/2)
72 real(kind_real),intent(in) :: a(nn) !< Matrix
73 real(kind_real),intent(out) :: u(nn) !< Matrix square-root
74 
75 ! Local variables
76 integer :: i,icol,ii,irow,j,k,kk,l,m
77 real(kind_real) :: w,x
78 
79 ! Set name
80 
81 
82 ! Probe in
83 
84 
85 ! Initialization
86 ii = 0
87 j = 1
88 k = 0
89 if (nn/=(n*(n+1))/2) then
90  call mpl%abort('asa007_cholesky','wrong size in Cholesky decomposition')
91 end if
92 w = zero
93 
94 ! Factorize column by column, ICOL = column number
95 do icol=1,n
96  ii = ii+icol
97  x = eta**2*a(ii)
98  l = 0
99  kk = 0
100 
101  ! IROW = row number within column ICOL
102  do irow=1,icol
103  kk = kk+irow
104  k = k+1
105  w = a(k)
106  m = j
107  do i=1,irow-1
108  l = l+1
109  w = w-u(l)*u(m)
110  m = m+1
111  end do
112  l = l+1
113  if (irow==icol) exit
114  if (abs(u(l))>zero) then
115  u(k) = w/u(l)
116  else
117  u(k) = zero
118  if (inf(abs(x*a(k)),w**2)) call mpl%abort('asa007_cholesky','matrix is not positive semi-definite')
119  end if
120  end do
121 
122  ! End of row, estimate relative accuracy of diagonal element
123  if (infeq(abs(w),abs(eta*a(k)))) then
124  u(k) = zero
125  else
126  if (w<zero) call mpl%abort('asa007_cholesky','matrix is not positive semi-definite')
127  u(k) = sqrt(w)
128  end if
129  j = j+icol
130 end do
131 
132 ! Probe out
133 
134 
135 end subroutine asa007_cholesky
136 
137 !----------------------------------------------------------------------
138 ! Subroutine: asa007_syminv
139 !> Compute inverse of a symmetric matrix
140 !----------------------------------------------------------------------
141 subroutine asa007_syminv(mpl,n,nn,a,c)
142 
143 implicit none
144 
145 ! Passed variables
146 type(mpl_type),intent(inout) :: mpl !< MPI data
147 integer,intent(in) :: n !< Matrix rank
148 integer,intent(in) :: nn !< Half-matrix size (n*(n-1)/2)
149 real(kind_real),intent(in) :: a(nn) !< Matrix
150 real(kind_real),intent(out) :: c(nn) !< Matrix inverse
151 
152 ! Local variables
153 integer :: i,icol,irow,j,jcol,k,l,mdiag,ndiag,nrow
154 real(kind_real) :: w(n),x
155 
156 ! Set name
157 
158 
159 ! Probe in
160 
161 
162 ! Initialization
163 nrow = n
164 if (nn/=(n*(n+1))/2) then
165  call mpl%abort('asa007_syminv','wrong size in matrix inversion')
166 end if
167 w = zero
168 
169 ! Compute the Cholesky factorization of A
170 call cholesky(mpl,n,nn,a,c)
171 
172 ! Invert C and form the product (Cinv)' * Cinv, where Cinv is the inverse of C, row by row starting with the last row
173 irow = nrow
174 ndiag = nn
175 
176 do
177  if (abs(c(ndiag))>zero) then
178  ! General case
179  l = ndiag
180  do i=irow,nrow
181  w(i) = c(l)
182  l = l+i
183  end do
184  icol = nrow
185  jcol = nn
186  mdiag = nn
187  do
188  l = jcol
189  if (icol==irow) then
190  x = one/w(irow)
191  else
192  x = zero
193  end if
194  k = nrow
195  do while (irow<k)
196  x = x-w(k)*c(l)
197  k = k-1
198  l = l-1
199  if (mdiag<l) l = l-k+1
200  end do
201  c(l) = x/w(irow)
202  if (icol<=irow) exit
203  mdiag = mdiag-icol
204  icol = icol-1
205  jcol = jcol-1
206  end do
207  else
208  ! Special case, zero diagonal element
209  l = ndiag
210  do j=irow,nrow
211  c(l) = zero
212  l = l+j
213  end do
214  end if
215  ndiag = ndiag-irow
216  irow = irow-1
217  if (irow<=0) exit
218 end do
219 
220 ! Probe out
221 
222 
223 end subroutine asa007_syminv
224 
225 end module tools_asa007
Subroutines/functions list.
real(kind_real), parameter eta
Small parameter for the Cholesky decomposition.
subroutine asa007_syminv(mpl, n, nn, a, c)
Compute inverse of a symmetric matrix.
subroutine asa007_cholesky(mpl, n, nn, a, u)
Compute cholesky decomposition.
Subroutines/functions list.
Definition: tools_const.F90:31
real(kind_real), parameter, public one
One.
Definition: tools_const.F90:42
real(kind_real), parameter, public zero
Zero.
Definition: tools_const.F90:37
Kinds definition.
Definition: tools_kinds.F90:9
integer, parameter, public kind_real
Real kind alias for the whole code.
Definition: tools_kinds.F90:25
Generic ranks, dimensions and types.
Definition: tools_repro.F90:42
Generic ranks, dimensions and types.
Definition: type_mpl.F90:42