-
Notifications
You must be signed in to change notification settings - Fork 1
/
pfit.f90
208 lines (177 loc) · 6.52 KB
/
pfit.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
! pfit.f90: Module for polynomial least-squares fitting
! http://infty.net/pfit/pfit.html
! v0.8.1
!
! Copyright (c) 2010-2013 Christopher N. Gilbreth
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
module mod_pfit
integer, parameter, private :: rk = selected_real_kind(p=15)
contains
subroutine pfit(x,y,sig,a,coeff,cov)
! Fit data to a polynomial a_0 + a_1 x + ... + a_d x**d
! Inputs:
! x(1:npt) - abscissas
! y(1:npt) - data values
! sig(1:npt) - data errors
! Outputs:
! a(1:d+1) - max. likelihood parameters
! coeff(1:npt,1:d+1) - coefficients giving the max. likelihood parameters
! in terms of the data:
! a(i) = \Sum_{j} coeff(j,i) * y(j)
! cov(1:n,1:d+1) - Covariance matrix, cov(i,j) = Cov(a(i),a(j))
! The estimated error in a(i) is sqrt(Cov(a(i),a(i))).
! Notes:
! This routine uses a QR decomposition method, which should be more
! numerically stable than solving the normal equations.
implicit none
real(rk), intent(in) :: x(:), y(:), sig(:)
real(rk), intent(out) :: a(:)
real(rk), intent(out) :: coeff(:,:), cov(:,:)
real(rk), allocatable :: work(:), C(:,:), Q(:,:), R(:,:), b(:)
integer :: ipiv(size(a)), lwork
integer :: i,j,k,d,npt,ifail
npt = size(x) ! Number of data points
d = size(a)-1 ! Max degree of polynomial
if (size(y) .ne. npt) stop "Error 1 in pfit"
if (size(sig) .ne. npt) stop "Error 2 in pfit"
if (d+1 .gt. npt) stop "Error 4 in pfit"
if (size(coeff,1) .ne. npt) stop "Error 6 in pfit"
if (size(coeff,2) .ne. d+1) stop "Error 5 in pfit"
if (size(cov,1) .ne. d+1) stop "Error 7 in pfit"
if (size(cov,2) .ne. d+1) stop "Error 8 in pfit"
allocate(C(npt,d+1), Q(npt,d+1), R(d+1,d+1), b(d+1))
! Vandermonde matrix
do j=1,d+1
do i=1,npt
C(i,j) = x(i)**(j-1)/sig(i)
end do
end do
! QR decomposition
call DQRF(C,Q,R,work)
! Inversion of R factor
R = dinverse(R)
! Compute max-likelihood parameters
! a = R^-1 Q^T y/σ
b = 0.d0
do j=1,d+1
do k=1,npt
b(j) = b(j) + Q(k,j) * y(k) / sig(k)
end do
end do
a = 0.d0
do i=1,d+1
do j=1,d+1
a(i) = a(i) + R(i,j) * b(j)
end do
end do
! Compute coefficient matrix coeff such that a(i) = Σ_j coeff(j,i) y(j)
! Here a(i) = R^{-1}(i,j) Q(k,j) y(k)/σ(k)
! So coeff(k,i) = R^{-1}(i,j) Q(k,j) / σ(k)
coeff = 0.d0
do i=1,d+1
do j=1,d+1
do k=1,npt
coeff(k,i) = coeff(k,i) + R(i,j) * Q(k,j) / sig(k)
end do
end do
end do
! Compute covariance matrix Cov(a(i),a(j)) = Σ_k C(k,i) C(k,j) σ(k)^2
cov = 0.d0
do j=1,d+1
do i=1,d+1
do k=1,npt
cov(i,j) = cov(i,j) + coeff(k,i) * coeff(k,j) * sig(k)**2
end do
end do
end do
end subroutine pfit
subroutine DQRF(A,Q,R,work)
! Compute the QR factorization of a general real matrix A:
! A = Q R
! where Q is unitary and R is upper triangular, using the LAPACK routine
! zgeqrf.
! Inputs:
! A: Matrix to be factorized, m x n
! Ouputs:
! Q: Unitary matrix, m x m
! R: Upper triangular, n x n
! Input/output:
! work: real(8) allocatable workspace array. If unallocated, this
! routine will allocate it to an appropriate size. If allocated,
! it is assumed to be the correct size for this problem.
implicit none
real(8), intent(in) :: A(:,:)
real(8), intent(out) :: Q(:,:), R(:,:)
real(8), allocatable :: work(:)
integer :: m, n, lwork, ierr, i, j
real(8) :: tau(size(A,2)), qwork(1)
real(8) :: A1(size(A,1),size(A,2))
external DGEQRF
external dorgqr
m = size(A,1)
n = size(A,2)
if (m .lt. n) stop "Error in DQRF: m < n"
if (size(Q,1) .ne. m) stop "Error in DQRF (2)"
if (size(Q,2) .ne. n) stop "Error in DQRF (3)"
if (size(R,1) .ne. n) stop "Error in DQRF (4)"
if (size(R,2) .ne. n) stop "Error in DQRF (5)"
A1 = A
if (.not. allocated(work)) then
! Compute size of workspace
lwork = -1
call DGEQRF(m, n, A1, m, TAU, qwork, LWORK, ierr)
if (ierr .ne. 0) stop "Error calling DGEQRF (1)"
lwork = qwork(1)
allocate(work(lwork))
end if
lwork = size(work)
call dgeqrf(m,n,A1,m,tau,work,lwork,ierr)
if (ierr .ne. 0) stop "Error calling DGEQRF (2)"
R = 0.d0
do j=1,n
do i=1,j
R(i,j) = A1(i,j)
end do
end do
Q(:,1:n) = A1
call dorgqr(m,n,n,Q,m,tau,work,lwork,ierr)
if (ierr .ne. 0) stop "Error calling DORGQR"
end subroutine DQRF
function dinverse(A)
! Invert a square matrix
implicit none
real(8), intent(in) :: A(:,:)
real(8) :: dinverse(size(A,1),size(A,1))
integer :: ipiv(size(A,1)), ierr, lwork
real(8), allocatable :: work(:)
real(8) :: work1(1)
external dgetri
dinverse = A
call dgetrf(size(A,1), size(A,1), dinverse, size(A,1), ipiv, ierr)
if (ierr .ne. 0) stop "Error computing LU decomposition for matrix inverse."
lwork = -1
call dgetri(size(A,1), dinverse, size(A,1), ipiv, work1, lwork, ierr)
if (ierr.ne.0) stop "Error allocating space for dgetri"
lwork = int(work1(1),kind(lwork))
allocate(work(max(1,lwork)))
call dgetri(size(A,1), dinverse, size(A,1), ipiv, work, lwork, ierr)
if (ierr .ne. 0) stop "Error calling zgetri."
end function dinverse
end module mod_pfit