forked from mjlaine/mcmcf90
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmatutils.F90
1847 lines (1652 loc) · 48.3 KB
/
matutils.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
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!!! $Id: matutils.F90,v 1.49 2017/08/10 07:18:33 mjlaine Exp $
!!! ------------------------------------------------------------------------
!!! this file is part of mcmc library
!!! File: matutils.F90
!!! Purpose: utility subroutines for the mcmc module
!!! matrix utilities
!!! error messages
!!! allocate utils
!!! file load and save
!!!
!!! Calls BLAS and LAPACK for matrix and vector operations, and
!!! auxiliary files dchud.f dchdd.f (from LINPACK) for Cholesky
!!! up/downdate
!!!
!!! Marko Laine 2001 <[email protected]>
!!! Copyrights licensed under a MIT License.
!!! See the accompanying LICENSE.txt file for terms.
!!!
module matutils
!#define EPPES 1 ! uncomment this if compiling the EPPES version
#define FILEUNIT 1234
use mcmcprec ! precision and machine constants
implicit none
character(len=*), parameter :: utils_version="$Revision: 1.49 $ Compiled: "//__DATE__
include 'lapack_inc.h'
interface writedata
module procedure writedata_mat, writedata_vec, writedata_scal
end interface
interface readdata
module procedure readdata_vec, readdata_mat, readdata_n, readdata_x
end interface
interface reallocate
module procedure reallocate_rv, reallocate_rm
end interface
interface loaddata
module procedure loaddata_mat, loaddata_vec, loaddata_x, &
loaddata_int, loaddata_ints
end interface
interface loaddata2
module procedure loaddata_mata, loaddata_veca
end interface
interface printmat
module procedure printmat_mat, printmat_vec, printmat_ivec
end interface
!!! A problem with output format here, the old G20.11 produces numbers like 0.12323200000-100
!!! with missing "E" is 3 numbers in the exponent, on th eother hand some older gfortran (4.3.4)
!!! does not allow for G0
#ifdef __GFORTRAN__
character(len=*), parameter :: gfrmstr = 'G0'
#else
character(len=*), parameter :: gfrmstr = 'G0'
#endif
character(len=526), save :: msgbuff
integer, save :: UtilInfoLevel = 0
integer, save :: UtilErrorLevel = 0
integer, parameter :: InfoLevelMinimal = 0
integer, parameter :: InfoLevelDebug = 2
integer, parameter :: InfoLevelInformational = 3
private :: gfrmstr, UtilInfoLevel, UtilErrorLevel, msgbuff
contains
!!! we need these matrix multiplication utilities to use the Cholesky
!!! factor provided by LAPACK routine, where C=R'R, R upper diagonal
!!! and the lower diagonal entries are arbitrary.
!!! TODO: general matmulu(A,B, trans) for AB A'B and AX, also need A'BA
!!! with B matrix or vector
!!!
!!! vector matrix multiplication
!!! p = x'*v (default) or p = x*v , x upper diagonal
!!!
function matmulu(x,v, trans) result(p)
implicit none
real(kind=dbl), intent(in) :: v(:)
real(kind=dbl), intent(in) :: x(:,:)
character(len=1), intent(in), optional :: trans
real(kind=dbl), dimension(size(v)) :: p
integer :: n,m
character(len=1) t
t = 'T'
if (present(trans)) then
if (trans == 'T' .or. trans == 't') then
t=trans
else
t='n'
end if
end if
n = size(v,1)
m = size(x,1)
if (n .ne. m) then
stop 'matmulu: incorrect dimensions'
end if
call dcopy(n,v,1,p,1) ! copy v to p p = v
call dtrmv('u',t,'n',n,x,n,p,1)
end function matmulu
!!!
!!! same but with lower diagonal x
!!!
function matmull(x,v) result(p)
implicit none
real(kind=dbl), intent(in) :: v(:)
real(kind=dbl), intent(in) :: x(:,:)
real(kind=dbl), dimension(size(v)) :: p
integer :: n,m
n = size(v,1)
m = size(x,1)
if (n .ne. m) then
stop 'matmull: incorrect dimensions'
end if
! call dsymv('l',n,1.0d0,x,n,v,1,0.0d0,p,1)
call dcopy(n,v,1,p,1) ! copy v to p p = v
call dtrmv('l','n','n',n,x,n,p,1)
end function matmull
!!!
!!! same but with general x (n * n)
!!! p = x*v (default) or x = x'*v
!!!
function matmulx(x,v, trans) result(p)
implicit none
real(kind=dbl), intent(in) :: v(:)
real(kind=dbl), intent(in) :: x(:,:)
character(len=1), intent(in), optional :: trans
real(kind=dbl), dimension(size(v)) :: p
integer :: n,m
character(len=1) t
t = 'N'
if (present(trans)) then
if (trans == 'T' .or. trans == 't') then
t=trans
else
t='n'
end if
end if
n = size(v,1)
m = size(x,1)
if (n .ne. m) then
stop 'matmulx: incorrect dimensions'
end if
call dgemv(t,n,n,1.0d0,x,n,v,1,0.0d0,p,1)
end function matmulx
!!!
!!! vector matrix multiplication
!!! p = x*v, x symmetric, using upper part
!!!
function matmuls(x,v) result(p)
implicit none
real(kind=dbl), intent(in) :: v(:)
real(kind=dbl), intent(in) :: x(:,:)
real(kind=dbl), dimension(size(v)) :: p
integer :: n,m
n = size(v,1)
m = size(x,1)
if (n .ne. m) then
stop 'matmuls: incorrect dimensions'
end if
call dsymv('u',n,1.0d0,x,n,v,1,0.0d0,p,1)
end function matmuls
!!!
!!! matrix matrix multiplication
!!! p = x*y, x upper diagonal, using upper part
!!!
!!!!! DO NOT USE YET !!!!
function matmulxy(x,y,side0,trans0) result(p)
implicit none
real(kind=dbl), intent(in) :: y(:,:)
real(kind=dbl), intent(in) :: x(:,:)
character(len=1), intent(in), optional :: side0, trans0
real(kind=dbl), dimension(size(x,1),size(y,2)) :: p
integer :: mx,nx,my,ny
character(len=1) side, trans
side = 'L'
if (present(side0)) then
if (side0 == 'R' .or. side0 == 'r') then
side=side0
end if
end if
trans = 'n'
if (present(trans0)) then
if (trans0 == 'T' .or. trans0 == 't') then
trans=trans0
end if
end if
mx = size(x,1)
nx = size(x,2)
my = size(y,1)
ny = size(y,2)
if (nx .ne. my) then
stop 'matmulxy: incorrect dimensions'
end if
! call dsymm(side,'u',mx,ny,1.0d0,x,mx,y,my,0.0d0,p,mx)
! ( SIDE, UPLO, M, N, ALPHA, A, LDA, B, LDB, BETA, C, LDC )
p = y
call dtrmm(side,'u',trans,'n',my,ny,1.0d0,x,mx,p,my)
!( SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB )
end function matmulxy
!!!
!!! calculates covariance matrix of x,
!!! with optional (frequency) weights in w,
!!! returns cmat and optionally mean and sum of weights in xmean and wsum
!!! updates cmat, cmean and wsum, if update=.true. and wsum>0
!!!
subroutine covmat(x, cmat, w, xmean, wsum, update)
implicit none
real(kind=dbl), intent(in) :: x(:,:)
real(kind=dbl), intent(inout) :: cmat(:,:)
real(kind=dbl), intent(in), optional :: w(:)
real(kind=dbl), intent(inout), optional :: xmean(size(cmat,2))
real(kind=dbl), intent(inout), optional :: wsum
logical, intent(in), optional :: update
! local variables
real(kind=dbl) :: xmean2(size(cmat,2))
integer :: n,p,i,j
real(kind=dbl) :: wsum2, w2, w3
logical :: doupdate
n=size(x,1)
p=size(x,2)
if (size(cmat,1) /= p .or. size(cmat,2) /= p) then
stop 'covmat: invalid cmat'
end if
!! optional weight in w
if (present(w) .and. size(w) == n) then
w2 = -1.0e0_dbl
wsum2 = sum(w)
else if (present(w) .and. size(w) == 1) then
w2 = w(1)
wsum2 = dble(n)*w2
else if (.not. present(w)) then
w2 = 1.0e0_dbl
wsum2 = dble(n)
else
stop 'covmat: invalid weight in w'
end if
!! update
doupdate=.false.
if (present(update) .and. present(wsum)) then
if (update .and. (wsum > 0.0_dbl)) then
doupdate=.true.
if (.not. present(xmean)) then
stop 'invalid input in covmat'
end if
else
doupdate=.false.
end if
end if
if (doupdate) then
!! covariance and mean update
!! one line update at a time
do i=1,n
xmean2 = x(i,:)-xmean ! center x and store it in xmean2
if (w2 == -1.0e0_dbl) then
w3 = w(i)
else
w3 = w2
end if
! cmat = cmat + w3/(wsum+w3-1.0d0) &
! * (wsum/(wsum+w3) &
! * matmul(reshape(xmean2,(/p,1/)), &
! reshape(xmean2,(/1,p/))) &
! - cmat)
!! new version, as matmul(reshape(x,(/p,1/)),reshape(x,(/1,p/)))
!! seems to generate run time error in gfortran 8.1.0
cmat = cmat + w3/(wsum+w3-1.0d0) &
* (wsum/(wsum+w3) &
* outer_product(xmean2,xmean2) &
- cmat)
!! or simplify to this (needs test)
! cmat = (wsum-1.0d0)/(w3+wsum-1.0d0)*cmat + &
! w3*wsum/(w3+wsum-1.0d0)/(w3+wsum) &
! * matmul(reshape(xmean2,(/p,1/)), &
! reshape(xmean2,(/1,p/)))
xmean = xmean + w3/(wsum+w3)*xmean2
wsum = w3 + wsum
end do
else !!! no update
do i=1,p
if (w2 == -1.0e0_dbl) then
xmean2(i) = sum(x(:,i)*w)/wsum2
else
xmean2(i) = sum(x(:,i)*w2)/wsum2
end if
end do
do i=1,p
do j=1,i
if ( w2 == -1.0e0_dbl ) then
cmat(i,j) = dot_product( (x(:,i)-xmean2(i)), &
(x(:,j)-xmean2(j))*w) / (wsum2 - 1.0_dbl)
else
cmat(i,j) = dot_product( (x(:,i)-xmean2(i)), &
(x(:,j)-xmean2(j))*w2) / (wsum2 - 1.0_dbl)
end if
if (i/=j) cmat(j,i) = cmat(i,j)
end do
end do
if (present(xmean)) then
xmean=xmean2
end if
if (present(wsum)) then
wsum=wsum2
end if
end if ! update
end subroutine covmat
!!!
!!! Calculates upper diagonal Cholesky transformation
!!! of symmetric matrix 'cmat'
subroutine covtor(cmat,R,info)
implicit none
real(kind=dbl) :: cmat(:,:)
real(kind=dbl), optional :: R(:,:)
integer, intent(out), optional :: info
integer :: info2, n
n = size(cmat,1)
if (size(cmat,2) /= n) then
stop 'size error in covtor'
end if
if (present(R)) then
if (size(R,1) /= n .or. size(R,1) /= size(R,2)) then
stop 'size error in R in covtor'
end if
R = cmat
call dpotrf('u',n,R,n,info2)
else
call dpotrf('u',n,cmat,n,info2)
end if
if (present(info)) then
info=info2
else if (info2 .ne. 0) then
stop 'error in Chol'
end if
end subroutine covtor
!!!
!!! test for svd sqrt
!!!
subroutine covtor_svd(cmat,R,condmax,info)
implicit none
real(kind=dbl) :: cmat(:,:)
real(kind=dbl) :: R(:,:)
real(kind=dbl), intent(in) :: condmax
integer, intent(out), optional :: info
integer :: info2, n, i
real(kind=dbl), allocatable :: work(:)
real(kind=dbl), allocatable :: u(:,:), s(:)
integer :: lwork
real(kind=dbl) :: tol, cond
n = size(cmat,1)
if (size(cmat,2) /= n) then
stop 'size error in covtor'
end if
lwork = 5*n
allocate(work(lwork),u(n,n),s(n), stat=info2)
if (info2 /= 0) then
stop 'memory allocation in svd'
end if
if (size(R,1) /= n .or. size(R,1) /= size(R,2)) then
stop 'size error in R in covtor'
end if
R = cmat
!!! call dgesvd(jobu,jobvt,m,n,a,lda,s,u,ldu,vt,ldvt,work,lwork,info)
call dgesvd('A','N',n,n,R,n,s,u,n,u,n,work,lwork,info2)
! write(*,*) 'S : ', real(s(max(1,n-3):n))
if (s(1) == 0.0d0) then
write(*,*) 'SVD: s(1)=0'
if (present(info)) then
info=n
end if
return
end if
if (info2 > 0) then
write(*,*) 'svd info = ', info2
end if
tol = dble(n)*s(1)*2.2204d-16 !!!! check
tol = dble(n)*s(1)*1.0d-6
cond = condmax
tol = s(1)/cond
if (s(n) <= tol) then
write(*,*) 'adjusting svd, tol:',real(tol)
do i=1,n
if (s(i) < tol) then
s(i) = tol
end if
end do
info2 = -1
end if
! R = u*sqrt(s)
R = u
!! scale R so that R'*R = cmat (or R*R' check!)
!! for scam we would need to have the unscaled version!
do i=1,n
!!! R(1:n,i) = sqrt(s(i))*R(1:n,i)
call dscal(n,sqrt(s(i)),R(1:n,i),1)
end do
if (present(info)) then
info=info2
else if (info2 < 0 .or. info2 > n-2) then
write(*,*) 'svd info = ', info2, ' Stopping'
stop
end if
deallocate(work,u,s)
end subroutine covtor_svd
!!!
!!! invert symmetric matrix using Cholesky decomposition
!!!
function invertmat_sym(x,info0) result(y)
implicit none
real(kind=dbl), intent(in) :: x(:,:)
integer, optional, intent(inout) :: info0
real(kind=dbl) :: y(size(x,1),size(x,2))
integer(kind=ik4) :: npar, info
npar = size(x,1)
if (npar .ne. size(x,2)) then
write(*,*) 'ERROR: x must be square'
stop
end if
y = x
call dpotrf('u',npar,y,npar,info)
if (info .ne. 0) then
if (present(info0)) then
info0 = info
return
else
write(*,*) 'ERROR: cannot factor cmat, info = ', info
stop
end if
end if
call dpotri('u',npar,y,npar,info)
if (info .ne. 0) then
if (present(info0)) then
info0 = info
return
else
write(*,*) 'ERROR: cannot invert cmat, info = ', info
stop
end if
end if
!! fill the lower diagonal with the upper
call copylower(y)
if (present(info0)) info0=0
end function invertmat_sym
!!!
!!! invert symmetric matrix using dsyrti
!!!
function invertmat_sym2(x,info0) result(y)
implicit none
real(kind=dbl), intent(in) :: x(:,:)
integer, optional, intent(inout) :: info0
real(kind=dbl) :: y(size(x,1),size(x,2))
real(kind=dbl) :: work0(1)
real(kind=dbl),allocatable :: work(:)
integer :: ipiv(size(x,1))
integer(kind=ik4) :: npar, lwork, info
npar = size(x,1)
if (npar .ne. size(x,2)) then
write(*,*) 'ERROR: x must be square'
stop
end if
y = x
! call dpotrf('u',npar,y,npar,info)
call dsytrf('u', npar, y, npar, ipiv, work0, -1, info)
lwork = max(npar,int(work0(1)))
allocate(work(lwork))
call dsytrf('u', npar, y, npar, ipiv, work, lwork, info)
if (info .ne. 0) then
if (present(info0)) then
info0 = info
deallocate(work)
return
else
write(*,*) 'ERROR: cannot factor cmat, info = ', info
stop
end if
end if
! call dpotri('u',npar,y,npar,info)
call dsytri('u', npar, y, npar, ipiv, work, info)
deallocate(work)
if (info .ne. 0) then
if (present(info0)) then
info0 = info
return
else
write(*,*) 'ERROR: cannot invert cmat, info = ', info
stop
end if
end if
!! fill the lower diagonal with the upper
call copylower(y)
if (present(info0)) info0=0
end function invertmat_sym2
!!!
!!! outer product of two vectors
!!!
function outer_product(x1,x2) result(y)
implicit none
real(kind=dbl), intent(in) :: x1(:), x2(:)
real(kind=dbl) :: y(size(x1),size(x2))
y = spread(x1,2,size(x2))*spread(x2,1,size(x1))
end function outer_product
!!!
!!! Mahalanobis distance of rows of x from mu given cmat
function mahalanobis(x,mu,cmat) result(d)
implicit none
real(kind=dbl) :: x(:,:), mu(:), cmat(:,:), d(size(x,1))
real(kind=dbl) :: cinv(size(cmat,1),size(cmat,2)), x0(size(x,1),size(x,2))
integer :: i, n
n = size(x,1)
! remove mean
do i=1,n
x0(i,:) = x(i,:)-mu
end do
! invert cmat
cinv=invertmat_sym(cmat)
! Mahalanobis distance for each row in x
d = sum(matmul(x0,cinv)*x0,2)
end function mahalanobis
!!!
!!! test for svd sqrt
!!!
subroutine scam_svd(cmat,R,std,condmax,info)
implicit none
real(kind=dbl) :: cmat(:,:)
real(kind=dbl) :: R(:,:)
real(kind=dbl) :: std(size(R,1))
real(kind=dbl), intent(in) :: condmax
integer, intent(out), optional :: info
integer :: info2, n, i
real(kind=dbl), allocatable :: work(:)
real(kind=dbl), allocatable :: u(:,:), s(:)
integer :: lwork
real(kind=dbl) :: tol, cond
n = size(cmat,1)
if (size(cmat,2) /= n) then
stop 'size error in covtor'
end if
lwork = 5*n
allocate(work(lwork),u(n,n),s(n), stat=info2)
if (info2 /= 0) then
stop 'memory allocation in svd'
end if
if (size(R,1) /= n .or. size(R,1) /= size(R,2)) then
stop 'size error in R in covtor'
end if
R = cmat
!!! call dgesvd(jobu,jobvt,m,n,a,lda,s,u,ldu,vt,ldvt,work,lwork,info)
call dgesvd('A','N',n,n,R,n,s,u,n,u,n,work,lwork,info2)
! write(*,*) 'S : ', real(s(max(1,n-3):n))
if (s(1) == 0.0d0) then
write(*,*) 'SVD: s(1)=0'
if (present(info)) then
info=n
end if
return
end if
if (info2 > 0) then
write(*,*) 'svd info = ', info2
end if
tol = dble(n)*s(1)*2.2204d-16 !!!! check
! tol = dble(n)*s(1)*1.0d-6
cond = condmax
tol = s(1)/cond
if (s(n) <= tol) then
write(*,*) 'adjusting svd, tol:',real(tol)
do i=1,n
if (s(i) < tol) then
s(i) = tol
end if
end do
info2 = -1
end if
! R = u*sqrt(s)
R = u
std = sqrt(s)
if (present(info)) then
info=info2
else if (info2 < 0 .or. info2 > n-2) then
write(*,*) 'svd info = ', info2, ' Stopping'
stop
end if
deallocate(work,u,s)
end subroutine scam_svd
!!!
!!! Cholesky update using dchud from Linpack
!!!
!!! for the EPPES version we do not want to include dchud.f
#ifndef EPPES
subroutine cholupdate(r,x,info)
implicit none
real(kind=dbl), intent(inout) :: r(:,:)
real(kind=dbl), intent(in) :: x(:)
integer, intent(out), optional :: info
integer ldr,p,ldz,nz
real(kind=dbl) rho(1), c(size(r,2))
real(kind=dbl) z(1,1),y(1),s(size(r,2))
interface
subroutine dchud(r,ldr,p,x,z,ldz,nz,y,rho,c,s)
implicit none
integer ldr,p,ldz,nz
double precision rho(nz),c(p)
double precision r(ldr,p),x(p),z(ldz,nz),y(nz),s(p)
end subroutine dchud
end interface
ldr = size(r,1)
p = size(r,2)
ldz = 0
nz = 0
call dchud(r,ldr,p,x,z,ldz,nz,y,rho,c,s)
if (present(info)) info = 0
end subroutine cholupdate
#endif
!!!
!!! Cholesky downdate using dchdd from Linpack
!!!
#ifndef EPPES
subroutine choldowndate(r,x,info)
implicit none
real(kind=dbl), intent(inout) :: r(:,:)
real(kind=dbl), intent(in) :: x(:)
integer, intent(out), optional :: info
integer ldr,p,ldz,nz
real(kind=dbl) rho(1), c(size(r,2))
real(kind=dbl) z(1,1),y(1),s(size(r,2))
integer info0
interface
subroutine dchdd(r,ldr,p,x,z,ldz,nz,y,rho,c,s,info)
implicit none
integer ldr,p,ldz,nz,info
double precision r(ldr,p),x(p),z(ldz,nz),y(nz),s(p)
double precision rho(nz),c(p)
end subroutine dchdd
end interface
ldr = size(r,1)
p = size(r,2)
ldz = 0
nz = 0
call dchdd(r,ldr,p,x,z,ldz,nz,y,rho,c,s,info0)
if (present(info)) then
info = info0
elseif (info0.ne.0) then
write(*,*) 'error in coldowndate, info:', info
stop
end if
end subroutine choldowndate
#endif
!!!
!!!
!!! error routines NOT IN USE YET
!!!
!!! SetMessageLevel - set doerror and message levels
!!!
subroutine SetMessageLevel(infolevel,errorlevel)
implicit none
integer, intent(in), optional :: infolevel
integer, intent(in), optional :: errorlevel
if (present(infolevel)) then
UtilInfoLevel = max(0,infolevel)
end if
if (present(errorlevel)) then
UtilErrorLevel = max(0,errorlevel)
end if
end subroutine SetMessageLevel
!!!
!!!GetInfoLevel - return current InfoLevel
!!!
integer function getinfolevel()
implicit none
getinfolevel = UtilInfoLevel
end function getinfolevel
!!!
!!!GetErrorLevel - return current ErrorLevel
!!!
integer function geterrorlevel()
implicit none
geterrorlevel = UtilErrorLevel
end function geterrorlevel
!!!
!!! DoError - error message and action
!!! action = 0, stop
!!! action = 1,
subroutine doerror(creason,elevel,action)
implicit none
character(len=*) , intent(in) :: creason ! character mesassage
integer, intent(in), optional :: elevel ! error level
integer, intent(in), optional :: action ! error action
integer :: act, el
if (present(elevel)) then
el = elevel
else
el = 0
end if
if (present(action)) then
act = action
else
act = 0
end if
if (el >= UtilErrorLevel) then
write(*,*) 'ERROR: ',trim(creason)
if (act .eq. 0) then
stop 'error stop'
end if
end if
return
end subroutine doerror
!!!
!!! Message - send text to console
!!!
subroutine message(text,ilevel)
implicit none
character(len=*), intent(in) :: text
integer, intent(in), optional :: ilevel
integer :: il
if(.not.present(ilevel)) then
il = 0
else
il = ilevel
end if
if (il >= UtilInfoLevel) then
write(*,*) text
end if
end subroutine message
!!! note routines
subroutine note_txt(verb,txt)
implicit none
integer, intent(in) :: verb
character(len=*), intent(in) :: txt
write(msgbuff,*) 'note: ',trim(txt)
call message(msgbuff,verb)
end subroutine note_txt
subroutine note_x(verb,txt,x)
implicit none
integer, intent(in) :: verb
character(len=*), intent(in) :: txt
real(kind=dbl), intent(in) :: x
write(msgbuff,*) 'note: ',trim(txt),' ',real(x)
call message(msgbuff,verb)
end subroutine note_x
subroutine note_c(verb,txt,c)
implicit none
integer, intent(in) :: verb
character(len=*), intent(in) :: txt
character(len=*), intent(in) :: c
write(msgbuff,*) 'note: ',trim(txt),' ',trim(c)
call message(msgbuff,verb)
end subroutine note_c
!!!
!!! write data matrix to file, with optional locking
!!!
subroutine writedata_mat(file,xmat,stat,uselock)
implicit none
character(len=*), intent(in) :: file
real(kind=dbl), intent(in) :: xmat(:,:)
integer, intent(out), optional :: stat
logical, intent(in), optional :: uselock
integer :: m,n, fstat, i, j
logical :: uselock2
if (len_trim(file) == 0) then
if (present(stat)) stat=-1
return
end if
uselock2 = .false.
if (present(uselock)) uselock2 = uselock
m = size(xmat,1)
n = size(xmat,2)
if (n < 1 .or. m < 1) then
write(*,*) 'Error in writedata, empty matrix, file:',trim(file)
return
end if
if (uselock2) then
call open_with_lock(file=file,unit=FILEUNIT,status='replace',iostat=fstat)
else
#ifdef _WIN32
open(unit=FILEUNIT, file=file, status='replace', &
recordtype='stream_cr', iostat=fstat)
#else
open(unit=FILEUNIT, file=file, status='replace', iostat=fstat)
#endif
end if
if (fstat /= 0) then
if (present(stat)) then
stat = fstat
return
else
write(*,*) 'Error opening file ', trim(file)
stop
end if
end if
do i=1,m
do j=1,n
! write(FILEUNIT,'('//gfrmstr//',x)', iostat=fstat, advance='NO') real(xmat(i,j))
write(FILEUNIT,'('//gfrmstr//',x)', iostat=fstat, advance='NO') xmat(i,j)
if (fstat /= 0) then
write(*,*) 'Write error on file ', trim(file)
stop
end if
end do
#ifdef __PGI
!! PGI fortran adds extra newline at eof
if (i<m) write(FILEUNIT,*)
#else
write(FILEUNIT,*)
#endif
end do
call close_with_lock(unit=FILEUNIT,file=file,iostat=fstat,uselock=uselock2)
if (present(stat)) stat=fstat
return
end subroutine writedata_mat
!!!
!!! write vector to file
!!!
subroutine writedata_vec(file,xvec,stat)
implicit none
character(len=*), intent(in) :: file
real(kind=dbl), intent(in) :: xvec(:)
integer, intent(out), optional :: stat
! real(kind=dbl) :: xmat(size(xvec),1)
integer :: i, n, fstat
if (len_trim(file) == 0) then
if (present(stat)) stat=-1
return
end if
n = size(xvec)
! xmat = reshape(xvec,(/n,1/)) !! stack overflow to happen here
! call writedata(file,xmat,stat)
#ifdef _WIN32
open(unit=FILEUNIT, file=file, status='replace', &
recordtype='stream_cr', iostat=fstat)
#else
open(unit=FILEUNIT, file=file, status='replace', iostat=fstat)
#endif
if (fstat /= 0) then
if (present(stat)) then
stat = fstat
return
else
write(*,*) 'Error opening file ', trim(file)
stop
end if
end if
do i=1,n
write(FILEUNIT,'('//gfrmstr//')', iostat=fstat, advance='YES') xvec(i)
if (fstat /= 0) then
write(*,*) 'Write error on file ', trim(file)
stop
end if
end do
close(unit=FILEUNIT)
if (present(stat)) stat=fstat
end subroutine writedata_vec
!!!
!!! write dbl scalar to ascii file
!!!
subroutine writedata_scal(file,x,stat)
implicit none
character(len=*), intent(in) :: file
real(kind=dbl), intent(in) :: x
integer, intent(out), optional :: stat
real(kind=dbl) :: xvec(1)
xvec(1) = x
call writedata(file,xvec,stat)
end subroutine writedata_scal
!!!
!!! write vector x to unit and try to avoid a newline
!!! this is not used anywhere?
subroutine writevec(unit,x,stat)
implicit none
integer, intent(in) :: unit
real(kind=dbl), intent(in) :: x(:)
integer, intent(out), optional :: stat
integer :: n, fstat, i
n = size(x)
do i=1,n
write(unit,'('//gfrmstr//',x)', iostat=fstat, advance='NO') x(i)
if (fstat /= 0) then
if (present(stat)) then
stat=fstat
return
else
write(*,*) 'Write error on unit ', unit
stop
end if
end if
end do
write(unit,*) ! newline
if (present(stat)) stat = 0
return
end subroutine writevec