-
Notifications
You must be signed in to change notification settings - Fork 27
/
scaling.f90
1721 lines (1565 loc) · 58.3 KB
/
scaling.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
! COPYRIGHT (c) 2014 The Science and Technology Facilities Council (STFC)
! Original date 18 December 2014, Version 1.0.0
!
! Written by: Jonathan Hogg
!
! Hungarian code derives from HSL MC64 code, but has been substantially
! altered for readability and to support rectangular matrices.
! All other code is fresh for SPRAL.
module spral_scaling
use spral_matrix_util, only : half_to_full
implicit none
private
! Top level routines
public :: auction_scale_sym, & ! Symmetric scaling by Auction algorithm
auction_scale_unsym, & ! Unsymmetric scaling by Auction algorithm
equilib_scale_sym, & ! Sym scaling by Equilibriation (MC77-like)
equilib_scale_unsym, & ! Unsym scaling by Equilibriation (MC77-lik)
hungarian_scale_sym, & ! Sym scaling by Hungarian alg (MC64-like)
hungarian_scale_unsym ! Unsym scaling by Hungarian alg (MC64-like)
! Inner routines that allow calling internals
public :: hungarian_match ! Find a matching (no pre/post-processing)
! Data types
public :: auction_options, auction_inform, &
equilib_options, equilib_inform, &
hungarian_options, hungarian_inform
integer, parameter :: wp = kind(0d0)
integer, parameter :: long = selected_int_kind(18)
real(wp), parameter :: rinf = huge(rinf)
type auction_options
integer :: max_iterations = 30000
integer :: max_unchanged(3) = (/ 10, 100, 100 /)
real :: min_proportion(3) = (/ 0.90, 0.0, 0.0 /)
real :: eps_initial = 0.01
end type auction_options
type auction_inform
integer :: flag = 0 ! success or failure
integer :: stat = 0 ! Fortran stat value on memory allocation failure
integer :: matched = 0 ! #matched rows/cols
integer :: iterations = 0 ! #iterations
integer :: unmatchable = 0 ! #classified as unmatchable
end type auction_inform
type equilib_options
integer :: max_iterations = 10
real :: tol = 1e-8
end type equilib_options
type equilib_inform
integer :: flag
integer :: stat
integer :: iterations
end type equilib_inform
type hungarian_options
logical :: scale_if_singular = .false.
end type hungarian_options
type hungarian_inform
integer :: flag
integer :: stat
integer :: matched
end type hungarian_inform
integer, parameter :: ERROR_ALLOCATION = -1
integer, parameter :: ERROR_SINGULAR = -2
integer, parameter :: WARNING_SINGULAR = 1
interface auction_scale_sym
module procedure auction_scale_sym_int32
module procedure auction_scale_sym_int64
end interface auction_scale_sym
interface auction_scale_unsym
module procedure auction_scale_unsym_int32
module procedure auction_scale_unsym_int64
end interface auction_scale_unsym
interface equilib_scale_sym
module procedure equilib_scale_sym_int32
module procedure equilib_scale_sym_int64
end interface equilib_scale_sym
interface equilib_scale_unsym
module procedure equilib_scale_unsym_int32
module procedure equilib_scale_unsym_int64
end interface equilib_scale_unsym
interface hungarian_scale_sym
module procedure hungarian_scale_sym_int32
module procedure hungarian_scale_sym_int64
end interface hungarian_scale_sym
interface hungarian_scale_unsym
module procedure hungarian_scale_unsym_int32
module procedure hungarian_scale_unsym_int64
end interface hungarian_scale_unsym
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Wrappers around scaling algorithms
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!**************************************************************
!
! Use matching-based scaling obtained using Hungarian algorithm (sym)
!
subroutine hungarian_scale_sym_int32(n, ptr, row, val, scaling, options, &
inform, match)
implicit none
integer, intent(in) :: n ! order of system
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(hungarian_options), intent(in) :: options
type(hungarian_inform), intent(out) :: inform
integer, dimension(n), optional, intent(out) :: match
integer(long), dimension(:), allocatable :: ptr64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call hungarian_scale_sym_int64(n, ptr64, row, val, scaling, options, &
inform, match=match)
end subroutine hungarian_scale_sym_int32
subroutine hungarian_scale_sym_int64(n, ptr, row, val, scaling, options, &
inform, match)
implicit none
integer, intent(in) :: n ! order of system
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(hungarian_options), intent(in) :: options
type(hungarian_inform), intent(out) :: inform
integer, dimension(n), optional, intent(out) :: match
integer, dimension(:), allocatable :: perm
real(wp), dimension(:), allocatable :: rscaling, cscaling
inform%flag = 0 ! Initialize to success
allocate(rscaling(n), cscaling(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
if (present(match)) then
call hungarian_wrapper(.true., n, n, ptr, row, val, match, rscaling, &
cscaling, options, inform)
else
allocate(perm(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
call hungarian_wrapper(.true., n, n, ptr, row, val, perm, rscaling, &
cscaling, options, inform)
end if
scaling(1:n) = exp( (rscaling(1:n) + cscaling(1:n)) / 2 )
end subroutine hungarian_scale_sym_int64
!**************************************************************
!
! Use matching-based scaling obtained using Hungarian algorithm (unsym)
!
subroutine hungarian_scale_unsym_int32(m, n, ptr, row, val, rscaling, cscaling,&
options, inform, match)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of cols
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(hungarian_options), intent(in) :: options
type(hungarian_inform), intent(out) :: inform
integer, dimension(m), optional, intent(out) :: match
integer(long), dimension(:), allocatable :: ptr64
! Copy from int32 to int64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call hungarian_scale_unsym_int64(m, n, ptr64, row, val, rscaling, cscaling, &
options, inform, match=match)
end subroutine hungarian_scale_unsym_int32
subroutine hungarian_scale_unsym_int64(m, n, ptr, row, val, rscaling, cscaling,&
options, inform, match)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of cols
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(hungarian_options), intent(in) :: options
type(hungarian_inform), intent(out) :: inform
integer, dimension(m), optional, intent(out) :: match
integer, dimension(:), allocatable :: perm
inform%flag = 0 ! Initialize to success
! Call main routine
if (present(match)) then
call hungarian_wrapper(.false., m, n, ptr, row, val, match, rscaling, &
cscaling, options, inform)
else
allocate(perm(m), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
call hungarian_wrapper(.false., m, n, ptr, row, val, perm, rscaling, &
cscaling, options, inform)
end if
! Apply post processing
rscaling(1:m) = exp( rscaling(1:m) )
cscaling(1:n) = exp( cscaling(1:n) )
end subroutine hungarian_scale_unsym_int64
!**************************************************************
!
! Call auction algorithm to get a scaling, then symmetrize it
!
subroutine auction_scale_sym_int32(n, ptr, row, val, scaling, options, inform, match)
implicit none
integer, intent(in) :: n ! order of system
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(auction_options), intent(in) :: options
type(auction_inform), intent(out) :: inform
integer, dimension(n), optional, intent(out) :: match
integer(long), dimension(:), allocatable :: ptr64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call auction_scale_sym_int64(n, ptr64, row, val, scaling, options, inform, &
match=match)
end subroutine auction_scale_sym_int32
subroutine auction_scale_sym_int64(n, ptr, row, val, scaling, options, inform, &
match)
implicit none
integer, intent(in) :: n ! order of system
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(auction_options), intent(in) :: options
type(auction_inform), intent(out) :: inform
integer, dimension(n), optional, intent(out) :: match
integer, dimension(:), allocatable :: perm
real(wp), dimension(:), allocatable :: rscaling, cscaling
inform%flag = 0 ! Initialize to sucess
! Allocate memory
allocate(rscaling(n), cscaling(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
! Call unsymmetric implementation with flag to expand half matrix
if (present(match)) then
call auction_match(.true., n, n, ptr, row, val, match, rscaling, &
cscaling, options, inform)
else
allocate(perm(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
call auction_match(.true., n, n, ptr, row, val, perm, rscaling, &
cscaling, options, inform)
end if
! Average rscaling and cscaling to get symmetric scaling
scaling(1:n) = exp( (rscaling(1:n)+cscaling(1:n))/2 )
end subroutine auction_scale_sym_int64
!**************************************************************
!
! Call auction algorithm to get a scaling (unsymmetric version)
!
subroutine auction_scale_unsym_int32(m, n, ptr, row, val, rscaling, cscaling, &
options, inform, match)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of columns
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(auction_options), intent(in) :: options
type(auction_inform), intent(out) :: inform
integer, dimension(m), optional, intent(out) :: match
integer(long), dimension(:), allocatable :: ptr64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call auction_scale_unsym_int64(m, n, ptr64, row, val, rscaling, cscaling, &
options, inform, match=match)
end subroutine auction_scale_unsym_int32
subroutine auction_scale_unsym_int64(m, n, ptr, row, val, rscaling, cscaling, &
options, inform, match)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of columns
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(auction_options), intent(in) :: options
type(auction_inform), intent(out) :: inform
integer, dimension(m), optional, intent(out) :: match
integer, dimension(:), allocatable :: perm
inform%flag = 0 ! Initialize to sucess
if (present(match)) then
call auction_match(.false., m, n, ptr, row, val, match, rscaling, &
cscaling, options, inform)
else
allocate(perm(m), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
call auction_match(.false., m, n, ptr, row, val, perm, rscaling, &
cscaling, options, inform)
end if
rscaling(1:m) = exp(rscaling(1:m))
cscaling(1:n) = exp(cscaling(1:n))
end subroutine auction_scale_unsym_int64
!*******************************
!
! Call the infinity-norm equilibriation algorithm (symmetric version)
!
subroutine equilib_scale_sym_int32(n, ptr, row, val, scaling, options, inform)
implicit none
integer, intent(in) :: n ! order of system
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(out) :: inform
integer(long), dimension(:), allocatable :: ptr64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call equilib_scale_sym_int64(n, ptr64, row, val, scaling, options, inform)
end subroutine equilib_scale_sym_int32
subroutine equilib_scale_sym_int64(n, ptr, row, val, scaling, options, inform)
implicit none
integer, intent(in) :: n ! order of system
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(n), intent(out) :: scaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(out) :: inform
inform%flag = 0 ! Initialize to sucess
call inf_norm_equilib_sym(n, ptr, row, val, scaling, options, inform)
end subroutine equilib_scale_sym_int64
!*******************************
!
! Call the infinity-norm equilibriation algorithm (unsymmetric version)
!
subroutine equilib_scale_unsym_int32(m, n, ptr, row, val, rscaling, cscaling, &
options, inform)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of cols
integer, intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(out) :: inform
integer(long), dimension(:), allocatable :: ptr64
allocate(ptr64(n+1), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
ptr64(1:n+1) = ptr(1:n+1)
call equilib_scale_unsym_int64(m, n, ptr64, row, val, rscaling, cscaling, &
options, inform)
end subroutine equilib_scale_unsym_int32
subroutine equilib_scale_unsym_int64(m, n, ptr, row, val, rscaling, cscaling, &
options, inform)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of cols
integer(long), intent(in) :: ptr(n+1) ! column pointers of A
integer, intent(in) :: row(*) ! row indices of A (lower triangle)
real(wp), intent(in) :: val(*) ! entries of A (in same order as in row).
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(out) :: inform
inform%flag = 0 ! Initialize to sucess
call inf_norm_equilib_unsym(m, n, ptr, row, val, rscaling, cscaling, &
options, inform)
end subroutine equilib_scale_unsym_int64
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Inf-norm Equilibriation Algorithm Implementation
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! We implement Algorithm 1 of:
! A Symmetry Preserving Algorithm for Matrix Scaling
! Philip Knight, Daniel Ruiz and Bora Ucar
! INRIA Research Report 7552 (Novemeber 2012)
! (This is similar to the algorithm used in MC77, but is a complete
! reimplementation from the above paper to ensure it is 100% STFC
! copyright and can be released as open source)
!
subroutine inf_norm_equilib_sym(n, ptr, row, val, scaling, options, inform)
implicit none
integer, intent(in) :: n
integer(long), dimension(n+1), intent(in) :: ptr
integer, dimension(ptr(n+1)-1), intent(in) :: row
real(wp), dimension(ptr(n+1)-1), intent(in) :: val
real(wp), dimension(n), intent(out) :: scaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(inout) :: inform
integer :: itr, r, c
integer(long) :: j
real(wp) :: v
real(wp), dimension(:), allocatable :: maxentry
allocate(maxentry(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
scaling(1:n) = 1.0
do itr = 1, options%max_iterations
! Find maximum entry in each row and col
! Recall: matrix is symmetric, but we only have half
maxentry(1:n) = 0.0
do c = 1, n
do j = ptr(c), ptr(c+1)-1
r = row(j)
v = abs(scaling(r) * val(j) * scaling(c))
maxentry(r) = max(maxentry(r), v)
maxentry(c) = max(maxentry(c), v)
end do
end do
! Update scaling (but beware empty cols)
where (maxentry(1:n) .gt. 0) &
scaling(1:n) = scaling(1:n) / sqrt(maxentry(1:n))
! Test convergence
if (maxval(abs(1-maxentry(1:n))) .lt. options%tol) exit
end do
inform%iterations = itr-1
end subroutine inf_norm_equilib_sym
!
! We implement Algorithm 1 of:
! A Symmetry Preserving Algorithm for Matrix Scaling
! Philip Knight, Daniel Ruiz and Bora Ucar
! INRIA Research Report 7552 (Novemeber 2012)
! (This is similar to the algorithm used in MC77, but is a complete
! reimplementation from the above paper to ensure it is 100% STFC
! copyright and can be released as open source)
!
subroutine inf_norm_equilib_unsym(m, n, ptr, row, val, rscaling, cscaling, &
options, inform)
implicit none
integer, intent(in) :: m
integer, intent(in) :: n
integer(long), dimension(n+1), intent(in) :: ptr
integer, dimension(ptr(n+1)-1), intent(in) :: row
real(wp), dimension(ptr(n+1)-1), intent(in) :: val
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(equilib_options), intent(in) :: options
type(equilib_inform), intent(inout) :: inform
integer :: itr, r, c
integer(long) :: j
real(wp) :: v
real(wp), dimension(:), allocatable :: rmaxentry, cmaxentry
allocate(rmaxentry(m), cmaxentry(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
rscaling(1:m) = 1.0
cscaling(1:n) = 1.0
do itr = 1, options%max_iterations
! Find maximum entry in each row and col
! Recall: matrix is symmetric, but we only have half
rmaxentry(1:m) = 0.0
cmaxentry(1:n) = 0.0
do c = 1, n
do j = ptr(c), ptr(c+1)-1
r = row(j)
v = abs(rscaling(r) * val(j) * cscaling(c))
rmaxentry(r) = max(rmaxentry(r), v)
cmaxentry(c) = max(cmaxentry(c), v)
end do
end do
! Update scaling (but beware empty cols)
where(rmaxentry(1:m).gt.0) &
rscaling(1:m) = rscaling(1:m) / sqrt(rmaxentry(1:m))
where(cmaxentry(1:n).gt.0) &
cscaling(1:n) = cscaling(1:n) / sqrt(cmaxentry(1:n))
! Test convergence
if ((maxval(abs(1-rmaxentry(1:m))) .lt. options%tol) .and. &
(maxval(abs(1-cmaxentry(1:n))) .lt. options%tol)) exit
end do
inform%iterations = itr-1
end subroutine inf_norm_equilib_unsym
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Hungarian Algorithm implementation (MC64)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!**********************************************************************
!
! This subroutine wraps the core algorithm of hungarian_match(). It provides
! pre- and post-processing to transform a maximum product assignment to a
! minimum sum assignment problem (and back again). It also has post-processing
! to handle the case of a structurally singular matrix as per Duff and Pralet
! (though the efficacy of such an approach is disputed!)
!
! This code is adapted from HSL_MC64 v2.3.1
!
subroutine hungarian_wrapper(sym, m, n, ptr, row, val, match, rscaling, &
cscaling, options, inform)
implicit none
logical, intent(in) :: sym
integer, intent(in) :: m
integer, intent(in) :: n
integer(long), dimension(n+1), intent(in) :: ptr
integer, dimension(*), intent(in) :: row
real(wp), dimension(*), intent(in) :: val
integer, dimension(m), intent(out) :: match
real(wp), dimension(m), intent(out) :: rscaling
real(wp), dimension(n), intent(out) :: cscaling
type(hungarian_options), intent(in) :: options
type(hungarian_inform), intent(out) :: inform
integer(long), allocatable :: ptr2(:)
integer, allocatable :: row2(:), iw(:), new_to_old(:), &
old_to_new(:), cperm(:)
real(wp), allocatable :: val2(:), dualu(:), dualv(:), cmax(:), cscale(:)
real(wp) :: colmax
integer :: i, j, nn, jj, k
integer(long) :: j1, j2, jlong, klong, ne
real(wp), parameter :: zero = 0.0
inform%flag = 0
inform%stat = 0
ne = ptr(n+1)-1
! Reset ne for the expanded symmetric matrix
ne = 2*ne
! Expand matrix, drop explicit zeroes and take log absolute values
allocate(ptr2(n+1), row2(ne), val2(ne), &
iw(5*n), dualu(m), dualv(n), cmax(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
klong = 1
do i = 1, n
ptr2(i) = klong
do jlong = ptr(i), ptr(i+1)-1
if (val(jlong) .eq. zero) cycle
row2(klong) = row(jlong)
val2(klong) = abs(val(jlong))
klong = klong + 1
end do
! Following log is seperated from above loop to expose expensive
! log operation to vectorization.
val2(ptr2(i):klong-1) = log(val2(ptr2(i):klong-1))
end do
ptr2(n+1) = klong
if (sym) call half_to_full(n, row2, ptr2, iw, a=val2)
! Compute column maximums
do i = 1, n
colmax = maxval(val2(ptr2(i):ptr2(i+1)-1))
cmax(i) = colmax
val2(ptr2(i):ptr2(i+1)-1) = colmax - val2(ptr2(i):ptr2(i+1)-1)
end do
call hungarian_match(m, n, ptr2, row2, val2, match, inform%matched, dualu, &
dualv, inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
if (inform%matched .ne. min(m,n)) then
! Singular matrix
if (options%scale_if_singular) then
! Just issue warning then continue
inform%flag = WARNING_SINGULAR
else
! Issue error and return identity scaling
inform%flag = ERROR_SINGULAR
rscaling(1:m) = 0
cscaling(1:n) = 0
end if
end if
if ((.not. sym) .or. (inform%matched .eq. n)) then ! Unsymmetric or symmetric and full rank
! Note that in this case m=n
rscaling(1:m) = dualu(1:m)
cscaling(1:n) = dualv(1:n) - cmax(1:n)
call match_postproc(m, n, ptr, row, val, rscaling, cscaling, &
inform%matched, match, inform%flag, inform%stat)
return
end if
! If we reach this point then structually rank deficient.
! As matching may not involve full set of rows and columns, but we need
! a symmetric matching/scaling, we can't just return the current matching.
! Instead, we build a full rank submatrix and call matching on it.
allocate(old_to_new(n),new_to_old(n),cperm(n),stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
j = inform%matched + 1
k = 0
do i = 1,m
if (match(i) < 0) then
! row i is not part of the matching
old_to_new(i) = -j
j = j + 1
else
k = k + 1
! old_to_new(i) holds the new index for variable i after
! removal of singular part and new_to_old(k) is the
! original index for k
old_to_new(i) = k
new_to_old(k) = i
end if
end do
! Overwrite ptr2, row2 and val2
ne = 0
k = 0
ptr2(1) = 1
j2 = 1
do i = 1,n
j1 = j2
j2 = ptr2(i+1)
! skip over unmatched entries
if (match(i) < 0) cycle
k = k + 1
do jlong = j1,j2-1
jj = row2(jlong)
if (match(jj) < 0) cycle
ne = ne + 1
row2(ne) = old_to_new(jj)
val2(ne) = val2(jlong)
end do
ptr2(k+1) = ne + 1
end do
! nn is order of non-singular part.
nn = k
call hungarian_match(nn, nn, ptr2, row2, val2, cperm, inform%matched, &
dualu, dualv, inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
do i = 1,n
j = old_to_new(i)
if (j < 0) then
rscaling(i) = -huge(rscaling)
else
! Note: we need to subtract col max using old matrix numbering
rscaling(i) = (dualu(j)+dualv(j)-cmax(i))/2
end if
end do
match(1:n) = -1
do i = 1,nn
j = cperm(i)
match(new_to_old(i)) = j
end do
do i = 1, n
if (match(i) .eq. -1) then
match(i) = old_to_new(i)
end if
end do
! Apply Duff and Pralet correction to unmatched row scalings
allocate(cscale(n), stat=inform%stat)
if (inform%stat .ne. 0) then
inform%flag = ERROR_ALLOCATION
return
end if
! For columns i not in the matched set I, set
! s_i = 1 / (max_{k in I} | a_ik s_k |)
! with convention that 1/0 = 1
cscale(1:n) = rscaling(1:n)
do i = 1,n
do jlong = ptr(i), ptr(i+1)-1
k = row(jlong)
if ((cscale(i) .eq. -huge(rscaling)) .and. (cscale(k) .ne. -huge(rscaling))) then
! i not in I, k in I
rscaling(i) = max(rscaling(i), log(abs(val(jlong)))+rscaling(k))
end if
if ((cscale(k) .eq. -huge(rscaling)) .and. (cscale(i) .ne. -huge(rscaling))) then
! k not in I, i in I
rscaling(k) = max(rscaling(k), log(abs(val(jlong)))+rscaling(i))
end if
end do
end do
do i = 1,n
if (cscale(i) .ne. -huge(rscaling)) cycle ! matched part
if (rscaling(i) .eq. -huge(rscaling)) then
rscaling(i) = 0.0
else
rscaling(i) = -rscaling(i)
end if
end do
! As symmetric, scaling is averaged on return, but rscaling(:) is correct,
! so just copy to cscaling to fix this
cscaling(1:n) = rscaling(1:n)
end subroutine hungarian_wrapper
!**********************************************************************
!
! Subroutine that initialize matching and (row) dual variable into a suitbale
! state for main Hungarian algorithm.
!
! The heuristic guaruntees that the generated partial matching is optimal
! on the restriction of the graph to the matched rows and columns.
subroutine hungarian_init_heurisitic(m, n, ptr, row, val, num, iperm, jperm, &
dualu, d, l, search_from)
implicit none
integer, intent(in) :: m
integer, intent(in) :: n
integer(long), dimension(n+1), intent(in) :: ptr
integer, dimension(ptr(n+1)-1), intent(in) :: row
real(wp), dimension(ptr(n+1)-1), intent(in) :: val
integer, intent(inout) :: num
integer, dimension(*), intent(inout) :: iperm
integer(long), dimension(*), intent(inout) :: jperm
real(wp), dimension(m), intent(out) :: dualu
real(wp), dimension(n), intent(out) :: d ! d(j) current improvement from
! matching in col j
integer(long), dimension(m), intent(out) :: l ! position of smallest entry
! of row
integer(long), dimension(n), intent(inout) :: search_from ! position we have
! reached in current search
integer :: i, i0, ii, j, jj
integer(long) :: k, k0, kk
real(wp) :: di, vj
!
! Set up initial matching on smallest entry in each row (as far as possible)
!
! Find smallest entry in each col, and record it
dualu(1:m) = RINF
l(1:m) = 0
do j = 1, n
do k = ptr(j),ptr(j+1)-1
i = row(k)
if (val(k) .gt. dualu(i)) cycle
dualu(i) = val(k) ! Initialize dual variables
iperm(i) = j ! Record col
l(i) = k ! Record posn in row(:)
end do
end do
! Loop over rows in turn. If we can match on smallest entry in row (i.e.
! column not already matched) then do so. Avoid matching on dense columns
! as this makes Hungarian algorithm take longer.
do i = 1, m
j = iperm(i) ! Smallest entry in row i is (i,j)
if (j .eq. 0) cycle ! skip empty rows
iperm(i) = 0
if (jperm(j) .ne. 0) cycle ! If we've already matched column j, skip this row
! Don't choose cheap assignment from dense columns
if ((ptr(j+1)-ptr(j) .gt. m/10) .and. (m .gt. 50)) cycle
! Assignment of column j to row i
num = num + 1
iperm(i) = j
jperm(j) = l(i)
end do
! If we already have a complete matching, we're already done!
if (num .eq. min(m,n)) return
!
! Scan unassigned columns; improve assignment
!
d(1:n) = 0.0
search_from(1:n) = ptr(1:n)
improve_assign: do j = 1, n
if (jperm(j) .ne. 0) cycle ! column j already matched
if (ptr(j) .gt. (ptr(j+1)-1)) cycle ! column j is empty
! Find smallest value of di = a_ij - u_i in column j
! In case of a tie, prefer first unmatched, then first matched row.
i0 = row(ptr(j))
vj = val(ptr(j)) - dualu(i0)
k0 = ptr(j)
do k = ptr(j)+1, ptr(j+1)-1
i = row(k)
di = val(k) - dualu(i)
if (di .gt. vj) cycle
if ((di .eq. vj) .and. (di .ne. RINF)) then
if ((iperm(i) .ne. 0) .or. (iperm(i0) .eq. 0)) cycle
end if
vj = di
i0 = i
k0 = k
end do
! Record value of matching on (i0,j)
d(j) = vj
! If row i is unmatched, then match on (i0,j) immediately
if (iperm(i0) .eq. 0) then
num = num + 1
jperm(j) = k0
iperm(i0) = j
search_from(j) = k0 + 1
cycle
end if
! Otherwise, row i is matched. Consider all rows i in column j that tie
! for this vj value. Such a row currently matches on (i,jj). Scan column
! jj looking for an unmatched row ii that improves value of matching. If
! one exists, then augment along length 2 path (i,j)->(ii,jj)
do k = k0, ptr(j+1)-1
i = row(k)
if ((val(k)-dualu(i)) .gt. vj) cycle ! Not a tie for vj value
jj = iperm(i)
! Scan remaining part of assigned column jj
do kk = search_from(jj), ptr(jj+1)-1
ii = row(kk)
if (iperm(ii) .gt. 0) cycle ! row ii already matched
if ((val(kk)-dualu(ii)) .le. d(jj)) then
! By matching on (i,j) and (ii,jj) we do better than existing
! matching on (i,jj) alone.
jperm(jj) = kk
iperm(ii) = jj
search_from(jj) = kk + 1
num = num + 1
jperm(j) = k
iperm(i) = j
search_from(j) = k + 1
cycle improve_assign
end if
end do
search_from(jj) = ptr(jj+1)
end do
cycle
end do improve_assign
end subroutine hungarian_init_heurisitic
!**********************************************************************
!
! Provides the core Hungarian Algorithm implementation for solving the
! minimum sum assignment problem as per Duff and Koster.
!
! This code is adapted from MC64 v 1.6.0
!
subroutine hungarian_match(m,n,ptr,row,val,iperm,num,dualu,dualv,st)
implicit none
integer, intent(in) :: m ! number of rows
integer, intent(in) :: n ! number of cols
integer, intent(out) :: num ! cardinality of the matching
integer(long), intent(in) :: ptr(n+1) ! column pointers
integer, intent(in) :: row(ptr(n+1)-1) ! row pointers
integer, intent(out) :: iperm(m) ! matching itself: row i is matched to
! column iperm(i).
real(wp), intent(in) :: val(ptr(n+1)-1) ! value of the entry that corresponds
! to row(k). All values val(k) must be non-negative.
real(wp), intent(out) :: dualu(m) ! dualu(i) is the reduced weight for row(i)
real(wp), intent(out) :: dualv(n) ! dualv(j) is the reduced weight for col(j)
integer, intent(out) :: st
integer(long), allocatable, dimension(:) :: jperm ! a(jperm(j)) is entry of
! A for matching in column j.
integer(long), allocatable, dimension(:) :: out ! a(out(i)) is the new entry
! in a on which we match going along the short path back to original col.
integer, allocatable, dimension(:) :: pr ! pr(i) is a pointer to the next
! column along the shortest path back to the original column
integer, allocatable, dimension(:) :: q ! q(1:qlen) forms a binary heap
! data structure sorted by d(q(i)) value. q(low:up) is a list of rows
! with equal d(i) which is lower or equal to smallest in the heap.
! q(up:n) is a list of already visited rows.
integer(long), allocatable, dimension(:) :: longwork
integer, allocatable, dimension(:) :: l ! l(:) is an inverse of q(:)
real(wp), allocatable, dimension(:) :: d ! d(i) is current shortest distance
! to row i from current column (d_i from Fig 4.1 of Duff and Koster paper)
integer :: i,j,jj,jord,q0,qlen,jdum,jsp
integer :: kk,up,low,lpos,k
integer(long) :: klong, isp
real(wp) :: csp,di,dmin,dnew,dq0,vj
!
! Initialization
!
allocate(jperm(n), out(n), pr(n), q(m), longwork(m), l(m), d(max(m,n)), &
stat=st)
if (st .ne. 0) return
num = 0
iperm(1:m) = 0
jperm(1:n) = 0
call hungarian_init_heurisitic(m, n, ptr, row, val, num, iperm, jperm, &
dualu, d, longwork, out)
if (num .eq. min(m,n)) go to 1000 ! If we got a complete matching, we're done
!
! Repeatedly find augmenting paths until all columns are included in the
! matching. At every step the current matching is optimal on the restriction
! of the graph to currently matched rows and columns.
!
! Main loop ... each pass round this loop is similar to Dijkstra's
! algorithm for solving the single source shortest path problem
d(1:m) = RINF
l(1:m) = 0
isp=-1; jsp=-1 ! initalize to avoid "may be used unitialized" warning
do jord = 1, n
if (jperm(jord) .ne. 0) cycle